backgroundMeshDecomposition.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "meshSearch.H"
28 #include "conformationSurfaces.H"
30 #include "Time.H"
31 #include "Random.H"
32 #include "pointConversion.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(backgroundMeshDecomposition, 0);
39 }
40 
41 
42 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
43 
45 (
46  const List<label>& toProc
47 )
48 {
49  // Determine send map
50  // ~~~~~~~~~~~~~~~~~~
51 
52  // 1. Count
53  labelList nSend(Pstream::nProcs(), 0);
54 
55  forAll(toProc, i)
56  {
57  label proci = toProc[i];
58 
59  nSend[proci]++;
60  }
61 
62 
63  // 2. Size sendMap
64  labelListList sendMap(Pstream::nProcs());
65 
66  forAll(nSend, proci)
67  {
68  sendMap[proci].setSize(nSend[proci]);
69 
70  nSend[proci] = 0;
71  }
72 
73  // 3. Fill sendMap
74  forAll(toProc, i)
75  {
76  label proci = toProc[i];
77 
78  sendMap[proci][nSend[proci]++] = i;
79  }
80 
81  // 4. Send over how many I need to receive
82  labelList recvSizes;
83  Pstream::exchangeSizes(sendMap, recvSizes);
84 
85 
86  // Determine receive map
87  // ~~~~~~~~~~~~~~~~~~~~~
88 
89  labelListList constructMap(Pstream::nProcs());
90 
91  // Local transfers first
92  constructMap[Pstream::myProcNo()] = identity
93  (
94  sendMap[Pstream::myProcNo()].size()
95  );
96 
97  label constructSize = constructMap[Pstream::myProcNo()].size();
98 
99  forAll(constructMap, proci)
100  {
101  if (proci != Pstream::myProcNo())
102  {
103  label nRecv = recvSizes[proci];
104 
105  constructMap[proci].setSize(nRecv);
106 
107  for (label i = 0; i < nRecv; i++)
108  {
109  constructMap[proci][i] = constructSize++;
110  }
111  }
112  }
113 
114  return autoPtr<mapDistribute>
115  (
116  new mapDistribute
117  (
118  constructSize,
119  move(sendMap),
120  move(constructMap)
121  )
122  );
123 }
124 
125 
126 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
127 
128 void Foam::backgroundMeshDecomposition::initialRefinement()
129 {
130  volScalarField cellWeights
131  (
132  IOobject
133  (
134  "cellWeights",
135  mesh_.time().timeName(),
136  mesh_,
139  ),
140  mesh_,
142  zeroGradientFvPatchScalarField::typeName
143  );
144 
145  const conformationSurfaces& geometry = geometryToConformTo_;
146 
147  decompositionMethod& decomposer = decomposerPtr_();
148 
149  volScalarField::Internal& icellWeights = cellWeights;
150 
151  // For each cell in the mesh has it been determined if it is fully
152  // inside, outside, or overlaps the surface
153  List<volumeType> volumeStatus
154  (
155  mesh_.nCells(),
157  );
158 
159  // Surface refinement
160  {
161  while (true)
162  {
163  // Determine/update the status of each cell
164  forAll(volumeStatus, celli)
165  {
166  if (volumeStatus[celli] == volumeType::unknown)
167  {
168  treeBoundBox cellBb
169  (
170  mesh_.cells()[celli].points
171  (
172  mesh_.faces(),
173  mesh_.points()
174  )
175  );
176 
177  if (geometry.overlaps(cellBb))
178  {
179  volumeStatus[celli] = volumeType::mixed;
180  }
181  else if (geometry.inside(cellBb.midpoint()))
182  {
183  volumeStatus[celli] = volumeType::inside;
184  }
185  else
186  {
187  volumeStatus[celli] = volumeType::outside;
188  }
189  }
190  }
191 
192  {
193  labelList refCells = selectRefinementCells
194  (
195  volumeStatus,
196  cellWeights
197  );
198 
199  // Maintain 2:1 ratio
200  labelList newCellsToRefine
201  (
202  meshCutter_.consistentRefinement
203  (
204  refCells,
205  true // extend set
206  )
207  );
208 
209  forAll(newCellsToRefine, nCTRI)
210  {
211  label celli = newCellsToRefine[nCTRI];
212 
213  if (volumeStatus[celli] == volumeType::mixed)
214  {
215  volumeStatus[celli] = volumeType::unknown;
216  }
217 
218  icellWeights[celli] = max
219  (
220  1.0,
221  icellWeights[celli]/8.0
222  );
223  }
224 
225  if (returnReduce(newCellsToRefine.size(), sumOp<label>()) == 0)
226  {
227  break;
228  }
229 
230  // Mesh changing engine.
231  polyTopoChange meshMod(mesh_);
232 
233  // Play refinement commands into mesh changer.
234  meshCutter_.setRefinement(newCellsToRefine, meshMod);
235 
236  // Create mesh, return map from old to new mesh.
237  autoPtr<mapPolyMesh> map = meshMod.changeMesh
238  (
239  mesh_,
240  false, // inflate
241  true, // syncParallel
242  true, // orderCells (to reduce cell transfers)
243  false // orderPoints
244  );
245 
246  // Update fields
247  mesh_.updateMesh(map);
248 
249  // Update numbering of cells/vertices.
250  meshCutter_.updateMesh(map);
251 
252  {
253  // Map volumeStatus
254 
255  const labelList& cellMap = map().cellMap();
256 
257  List<volumeType> newVolumeStatus(cellMap.size());
258 
259  forAll(cellMap, newCelli)
260  {
261  label oldCelli = cellMap[newCelli];
262 
263  if (oldCelli == -1)
264  {
265  newVolumeStatus[newCelli] = volumeType::unknown;
266  }
267  else
268  {
269  newVolumeStatus[newCelli] = volumeStatus[oldCelli];
270  }
271  }
272 
273  volumeStatus.transfer(newVolumeStatus);
274  }
275 
276  Info<< " Background mesh refined from "
277  << returnReduce(map().nOldCells(), sumOp<label>())
278  << " to " << mesh_.globalData().nTotalCells()
279  << " cells." << endl;
280  }
281 
282  // Determine/update the status of each cell
283  forAll(volumeStatus, celli)
284  {
285  if (volumeStatus[celli] == volumeType::unknown)
286  {
287  treeBoundBox cellBb
288  (
289  mesh_.cells()[celli].points
290  (
291  mesh_.faces(),
292  mesh_.points()
293  )
294  );
295 
296  if (geometry.overlaps(cellBb))
297  {
298  volumeStatus[celli] = volumeType::mixed;
299  }
300  else if (geometry.inside(cellBb.midpoint()))
301  {
302  volumeStatus[celli] = volumeType::inside;
303  }
304  else
305  {
306  volumeStatus[celli] = volumeType::outside;
307  }
308  }
309  }
310 
311  // Hard code switch for this stage for testing
312  bool removeOutsideCells = false;
313 
314  if (removeOutsideCells)
315  {
316  DynamicList<label> cellsToRemove;
317 
318  forAll(volumeStatus, celli)
319  {
320  if (volumeStatus[celli] == volumeType::outside)
321  {
322  cellsToRemove.append(celli);
323  }
324  }
325 
326  removeCells cellRemover(mesh_);
327 
328  // Mesh changing engine.
329  polyTopoChange meshMod(mesh_);
330 
331  labelList exposedFaces = cellRemover.getExposedFaces
332  (
333  cellsToRemove
334  );
335 
336  // Play refinement commands into mesh changer.
337  cellRemover.setRefinement
338  (
339  cellsToRemove,
340  exposedFaces,
341  labelList(exposedFaces.size(), 0), // patchID dummy
342  meshMod
343  );
344 
345  // Create mesh, return map from old to new mesh.
346  autoPtr<mapPolyMesh> map = meshMod.changeMesh
347  (
348  mesh_,
349  false, // inflate
350  true, // syncParallel
351  true, // orderCells (to reduce cell transfers)
352  false // orderPoints
353  );
354 
355  // Update fields
356  mesh_.updateMesh(map);
357 
358  // Update numbering of cells/vertices.
359  meshCutter_.updateMesh(map);
360  cellRemover.updateMesh(map);
361 
362  {
363  // Map volumeStatus
364 
365  const labelList& cellMap = map().cellMap();
366 
367  List<volumeType> newVolumeStatus(cellMap.size());
368 
369  forAll(cellMap, newCelli)
370  {
371  label oldCelli = cellMap[newCelli];
372 
373  if (oldCelli == -1)
374  {
375  newVolumeStatus[newCelli] = volumeType::unknown;
376  }
377  else
378  {
379  newVolumeStatus[newCelli] = volumeStatus[oldCelli];
380  }
381  }
382 
383  volumeStatus.transfer(newVolumeStatus);
384  }
385 
386  Info<< "Removed "
387  << returnReduce(map().nOldCells(), sumOp<label>())
388  - mesh_.globalData().nTotalCells()
389  << " cells." << endl;
390  }
391 
392  if (debug)
393  {
394  // const_cast<Time&>(mesh_.time())++;
395  // Info<< "Time " << mesh_.time().timeName() << endl;
396  meshCutter_.write();
397  mesh_.write();
398  cellWeights.write();
399  }
400 
401  labelList newDecomp = decomposer.decompose
402  (
403  mesh_,
404  mesh_.cellCentres(),
405  icellWeights
406  );
407 
408  fvMeshDistribute distributor(mesh_, mergeDist_);
409 
410  autoPtr<mapDistributePolyMesh> mapDist = distributor.distribute
411  (
412  newDecomp
413  );
414 
415  meshCutter_.distribute(mapDist);
416 
417  mapDist().distributeCellData(volumeStatus);
418 
419  if (debug)
420  {
421  printMeshData(mesh_);
422 
423  // const_cast<Time&>(mesh_.time())++;
424  // Info<< "Time " << mesh_.time().timeName() << endl;
425  meshCutter_.write();
426  mesh_.write();
427  cellWeights.write();
428  }
429  }
430  }
431 
432  if (debug)
433  {
434  // const_cast<Time&>(mesh_.time())++;
435  // Info<< "Time " << mesh_.time().timeName() << endl;
436  cellWeights.write();
437  mesh_.write();
438  }
439 
440  buildPatchAndTree();
441 }
442 
443 
444 void Foam::backgroundMeshDecomposition::printMeshData
445 (
446  const polyMesh& mesh
447 ) const
448 {
449  // Collect all data on master
450 
451  globalIndex globalCells(mesh.nCells());
452  // labelListList patchNeiProcNo(Pstream::nProcs());
453  // labelListList patchSize(Pstream::nProcs());
454  // const labelList& pPatches = mesh.globalData().processorPatches();
455  // patchNeiProcNo[Pstream::myProcNo()].setSize(pPatches.size());
456  // patchSize[Pstream::myProcNo()].setSize(pPatches.size());
457  // forAll(pPatches, i)
458  // {
459  // const processorPolyPatch& ppp = refCast<const processorPolyPatch>
460  // (
461  // mesh.boundaryMesh()[pPatches[i]]
462  // );
463  // patchNeiProcNo[Pstream::myProcNo()][i] = ppp.neighbProcNo();
464  // patchSize[Pstream::myProcNo()][i] = ppp.size();
465  // }
466  // Pstream::gatherList(patchNeiProcNo);
467  // Pstream::gatherList(patchSize);
468 
469 
470  // // Print stats
471 
472  // globalIndex globalBoundaryFaces(mesh.nFaces()-mesh.nInternalFaces());
473 
474  for (label proci = 0; proci < Pstream::nProcs(); proci++)
475  {
476  Info<< "Processor " << proci << " "
477  << "Number of cells = " << globalCells.localSize(proci)
478  << endl;
479 
480  // label nProcFaces = 0;
481 
482  // const labelList& nei = patchNeiProcNo[proci];
483 
484  // forAll(patchNeiProcNo[proci], i)
485  // {
486  // Info<< " Number of faces shared with processor "
487  // << patchNeiProcNo[proci][i] << " = " << patchSize[proci][i]
488  // << endl;
489 
490  // nProcFaces += patchSize[proci][i];
491  // }
492 
493  // Info<< " Number of processor patches = " << nei.size() << nl
494  // << " Number of processor faces = " << nProcFaces << nl
495  // << " Number of boundary faces = "
496  // << globalBoundaryFaces.localSize(proci) << endl;
497  }
498 }
499 
500 
501 bool Foam::backgroundMeshDecomposition::refineCell
502 (
503  label celli,
504  volumeType volType,
505  scalar& weightEstimate
506 ) const
507 {
508  // Sample the box to find an estimate of the min size, and a volume
509  // estimate when overlapping == true.
510 
511 // const conformationSurfaces& geometry = geometryToConformTo_;
512 
513  treeBoundBox cellBb
514  (
515  mesh_.cells()[celli].points
516  (
517  mesh_.faces(),
518  mesh_.points()
519  )
520  );
521 
522  weightEstimate = 1.0;
523 
524  if (volType == volumeType::mixed)
525  {
526 // // Assess the cell size at the nearest point on the surface for the
527 // // MIXED cells, if the cell is large with respect to the cell size,
528 // // then refine it.
529 //
530 // pointField samplePoints
531 // (
532 // volRes_*volRes_*volRes_,
533 // Zero
534 // );
535 //
536 // // scalar sampleVol = cellBb.volume()/samplePoints.size();
537 //
538 // vector delta = cellBb.span()/volRes_;
539 //
540 // label pI = 0;
541 //
542 // for (label i = 0; i < volRes_; i++)
543 // {
544 // for (label j = 0; j < volRes_; j++)
545 // {
546 // for (label k = 0; k < volRes_; k++)
547 // {
548 // samplePoints[pI++] =
549 // cellBb.min()
550 // + vector
551 // (
552 // delta.x()*(i + 0.5),
553 // delta.y()*(j + 0.5),
554 // delta.z()*(k + 0.5)
555 // );
556 // }
557 // }
558 // }
559 //
560 // List<pointIndexHit> hitInfo;
561 // labelList hitSurfaces;
562 //
563 // geometry.findSurfaceNearest
564 // (
565 // samplePoints,
566 // scalarField(samplePoints.size(), sqr(great)),
567 // hitInfo,
568 // hitSurfaces
569 // );
570 //
571 // // weightEstimate = 0.0;
572 //
573 // scalar minCellSize = great;
574 //
575 // forAll(samplePoints, i)
576 // {
577 // scalar s = cellShapeControl_.cellSize
578 // (
579 // hitInfo[i].hitPoint()
580 // );
581 //
582 // // Info<< "cellBb.midpoint() " << cellBb.midpoint() << nl
583 // // << samplePoints[i] << nl
584 // // << hitInfo[i] << nl
585 // // << "cellBb.span() " << cellBb.span() << nl
586 // // << "cellBb.mag() " << cellBb.mag() << nl
587 // // << s << endl;
588 //
589 // if (s < minCellSize)
590 // {
591 // minCellSize = max(s, minCellSizeLimit_);
592 // }
593 //
594 // // Estimate the number of points in the cell by the surface size,
595 // // this is likely to be too small, so reduce.
596 // // weightEstimate += sampleVol/pow3(s);
597 // }
598 //
599 // if (sqr(spanScale_)*sqr(minCellSize) < magSqr(cellBb.span()))
600 // {
601 // return true;
602 // }
603  }
604  else if (volType == volumeType::inside)
605  {
606  // scalar s =
607  // foamyHexMesh_.cellShapeControl_.cellSize(cellBb.midpoint());
608 
609  // Estimate the number of points in the cell by the size at the cell
610  // midpoint
611  // weightEstimate = cellBb.volume()/pow3(s);
612 
613  return false;
614  }
615  // else
616  // {
617  // weightEstimate = 1.0;
618 
619  // return false;
620  // }
621 
622  return false;
623 }
624 
625 
626 Foam::labelList Foam::backgroundMeshDecomposition::selectRefinementCells
627 (
628  List<volumeType>& volumeStatus,
629  volScalarField& cellWeights
630 ) const
631 {
632  volScalarField::Internal& icellWeights = cellWeights;
633 
634  labelHashSet cellsToRefine;
635 
636  // Determine/update the status of each cell
637  forAll(volumeStatus, celli)
638  {
639  if (volumeStatus[celli] == volumeType::mixed)
640  {
641  if (meshCutter_.cellLevel()[celli] < minLevels_)
642  {
643  cellsToRefine.insert(celli);
644  }
645  }
646 
647  if (volumeStatus[celli] != volumeType::outside)
648  {
649  if
650  (
651  refineCell
652  (
653  celli,
654  volumeStatus[celli],
655  icellWeights[celli]
656  )
657  )
658  {
659  cellsToRefine.insert(celli);
660  }
661  }
662  }
663 
664  return cellsToRefine.toc();
665 }
666 
667 
668 void Foam::backgroundMeshDecomposition::buildPatchAndTree()
669 {
670  primitivePatch tmpBoundaryFaces
671  (
672  SubList<face>
673  (
674  mesh_.faces(),
675  mesh_.nFaces() - mesh_.nInternalFaces(),
676  mesh_.nInternalFaces()
677  ),
678  mesh_.points()
679  );
680 
681  boundaryFacesPtr_.reset
682  (
683  new bPatch
684  (
685  tmpBoundaryFaces.localFaces(),
686  tmpBoundaryFaces.localPoints()
687  )
688  );
689 
690  // Overall bb
691  treeBoundBox overallBb(boundaryFacesPtr_().localPoints());
692 
693  bFTreePtr_.reset
694  (
695  new indexedOctree<treeDataBPatch>
696  (
698  (
699  false,
700  boundaryFacesPtr_(),
702  ),
703  overallBb.extend(1e-4),
704  10, // maxLevel
705  10, // leafSize
706  3.0 // duplicity
707  )
708  );
709 
710  // Give the bounds of every processor to every other processor
711  allBackgroundMeshBounds_[Pstream::myProcNo()] = overallBb;
712 
713  Pstream::gatherList(allBackgroundMeshBounds_);
714  Pstream::scatterList(allBackgroundMeshBounds_);
715 
716  point bbMin(great, great, great);
717  point bbMax(-great, -great, -great);
718 
719  forAll(allBackgroundMeshBounds_, proci)
720  {
721  bbMin = min(bbMin, allBackgroundMeshBounds_[proci].min());
722  bbMax = max(bbMax, allBackgroundMeshBounds_[proci].max());
723  }
724 
725  globalBackgroundBounds_ = treeBoundBox(bbMin, bbMax);
726 
727  if (false)
728  {
729  OFstream fStr
730  (
731  mesh_.time().path()
732  /"backgroundMeshDecomposition_proc_"
734  + "_boundaryFaces.obj"
735  );
736 
737  const faceList& faces = boundaryFacesPtr_().localFaces();
738  const List<point>& points = boundaryFacesPtr_().localPoints();
739 
740  Map<label> foamToObj(points.size());
741 
742  label vertI = 0;
743 
744  forAll(faces, i)
745  {
746  const face& f = faces[i];
747 
748  forAll(f, fPI)
749  {
750  if (foamToObj.insert(f[fPI], vertI))
751  {
752  meshTools::writeOBJ(fStr, points[f[fPI]]);
753  vertI++;
754  }
755  }
756 
757  fStr<< 'f';
758 
759  forAll(f, fPI)
760  {
761  fStr<< ' ' << foamToObj[f[fPI]] + 1;
762  }
763 
764  fStr<< nl;
765  }
766  }
767 }
768 
769 
770 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
771 
773 (
774  const Time& runTime,
775  Random& rndGen,
776  const conformationSurfaces& geometryToConformTo,
777  const dictionary& coeffsDict
778 )
779 :
780  runTime_(runTime),
781  geometryToConformTo_(geometryToConformTo),
782  mesh_
783  (
784  IOobject
785  (
786  "backgroundMeshDecomposition",
787  runTime_.timeName(),
788  runTime_,
789  IOobject::MUST_READ,
790  IOobject::AUTO_WRITE,
791  false
792  )
793  ),
794  meshCutter_
795  (
796  mesh_,
797  labelList(mesh_.nCells(), 0),
798  labelList(mesh_.nPoints(), 0)
799  ),
800  boundaryFacesPtr_(),
801  bFTreePtr_(),
802  allBackgroundMeshBounds_(Pstream::nProcs()),
803  globalBackgroundBounds_(),
804  decomposeDict_
805  (
806  IOobject
807  (
808  "decomposeParDict",
809  runTime_.system(),
810  runTime_,
811  IOobject::MUST_READ_IF_MODIFIED,
812  IOobject::NO_WRITE
813  )
814  ),
815  decomposerPtr_(decompositionMethod::New(decomposeDict_)),
816  mergeDist_(1e-6*mesh_.bounds().mag()),
817  spanScale_(coeffsDict.lookup<scalar>("spanScale")),
818  minCellSizeLimit_
819  (
820  coeffsDict.lookupOrDefault<scalar>("minCellSizeLimit", 0.0)
821  ),
822  minLevels_(coeffsDict.lookup<label>("minLevels")),
823  volRes_(coeffsDict.lookup<label>("sampleResolution")),
824  maxCellWeightCoeff_(coeffsDict.lookup<scalar>("maxCellWeightCoeff"))
825 {
826  if (!Pstream::parRun())
827  {
829  << "This cannot be used when not running in parallel."
830  << exit(FatalError);
831  }
832 
833  if (!decomposerPtr_().parallelAware())
834  {
836  << "You have selected decomposition method "
837  << decomposerPtr_().typeName
838  << " which is not parallel aware." << endl
839  << exit(FatalError);
840  }
841 
842  Info<< nl << "Building initial background mesh decomposition" << endl;
843 
844  initialRefinement();
845 }
846 
847 
848 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
849 
851 {}
852 
853 
854 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
855 
858 (
859  volScalarField& cellWeights
860 )
861 {
862  if (debug)
863  {
864  // const_cast<Time&>(mesh_.time())++;
865  // Info<< "Time " << mesh_.time().timeName() << endl;
866  cellWeights.write();
867  mesh_.write();
868  }
869 
870  volScalarField::Internal& icellWeights = cellWeights;
871 
872  while (true)
873  {
874  // Refine large cells if necessary
875 
876  label nOccupiedCells = 0;
877 
878  forAll(icellWeights, cI)
879  {
880  if (icellWeights[cI] > 1 - small)
881  {
882  nOccupiedCells++;
883  }
884  }
885 
886  // Only look at occupied cells, as there is a possibility of runaway
887  // refinement if the number of cells grows too fast. Also, clip the
888  // minimum cellWeightLimit at maxCellWeightCoeff_
889 
890  scalar cellWeightLimit = max
891  (
892  maxCellWeightCoeff_
893  *sum(cellWeights).value()
894  /returnReduce(nOccupiedCells, sumOp<label>()),
895  maxCellWeightCoeff_
896  );
897 
898  if (debug)
899  {
900  Info<< " cellWeightLimit " << cellWeightLimit << endl;
901 
902  Pout<< " sum(cellWeights) " << sum(cellWeights.primitiveField())
903  << " max(cellWeights) " << max(cellWeights.primitiveField())
904  << endl;
905  }
906 
907  labelHashSet cellsToRefine;
908 
909  forAll(icellWeights, cWI)
910  {
911  if (icellWeights[cWI] > cellWeightLimit)
912  {
913  cellsToRefine.insert(cWI);
914  }
915  }
916 
917  if (returnReduce(cellsToRefine.size(), sumOp<label>()) == 0)
918  {
919  break;
920  }
921 
922  // Maintain 2:1 ratio
923  labelList newCellsToRefine
924  (
925  meshCutter_.consistentRefinement
926  (
927  cellsToRefine.toc(),
928  true // extend set
929  )
930  );
931 
932  if (debug && !cellsToRefine.empty())
933  {
934  Pout<< " cellWeights too large in " << cellsToRefine.size()
935  << " cells" << endl;
936  }
937 
938  forAll(newCellsToRefine, nCTRI)
939  {
940  label celli = newCellsToRefine[nCTRI];
941 
942  icellWeights[celli] /= 8.0;
943  }
944 
945  // Mesh changing engine.
946  polyTopoChange meshMod(mesh_);
947 
948  // Play refinement commands into mesh changer.
949  meshCutter_.setRefinement(newCellsToRefine, meshMod);
950 
951  // Create mesh, return map from old to new mesh.
952  autoPtr<mapPolyMesh> map = meshMod.changeMesh
953  (
954  mesh_,
955  false, // inflate
956  true, // syncParallel
957  true, // orderCells (to reduce cell motion)
958  false // orderPoints
959  );
960 
961  // Update fields
962  mesh_.updateMesh(map);
963 
964  // Update numbering of cells/vertices.
965  meshCutter_.updateMesh(map);
966 
967  Info<< " Background mesh refined from "
968  << returnReduce(map().nOldCells(), sumOp<label>())
969  << " to " << mesh_.globalData().nTotalCells()
970  << " cells." << endl;
971 
972  if (debug)
973  {
974  // const_cast<Time&>(mesh_.time())++;
975  // Info<< "Time " << mesh_.time().timeName() << endl;
976  cellWeights.write();
977  mesh_.write();
978  }
979  }
980 
981  if (debug)
982  {
983  printMeshData(mesh_);
984 
985  Pout<< " Pre distribute sum(cellWeights) "
986  << sum(icellWeights)
987  << " max(cellWeights) "
988  << max(icellWeights)
989  << endl;
990  }
991 
992  labelList newDecomp = decomposerPtr_().decompose
993  (
994  mesh_,
995  mesh_.cellCentres(),
996  icellWeights
997  );
998 
999  Info<< " Redistributing background mesh cells" << endl;
1000 
1001  fvMeshDistribute distributor(mesh_, mergeDist_);
1002 
1003  autoPtr<mapDistributePolyMesh> mapDist = distributor.distribute(newDecomp);
1004 
1005  meshCutter_.distribute(mapDist);
1006 
1007  if (debug)
1008  {
1009  printMeshData(mesh_);
1010 
1011  Pout<< " Post distribute sum(cellWeights) "
1012  << sum(icellWeights)
1013  << " max(cellWeights) "
1014  << max(icellWeights)
1015  << endl;
1016 
1017  // const_cast<Time&>(mesh_.time())++;
1018  // Info<< "Time " << mesh_.time().timeName() << endl;
1019  mesh_.write();
1020  cellWeights.write();
1021  }
1022 
1023  buildPatchAndTree();
1024 
1025  return mapDist;
1026 }
1027 
1028 
1030 (
1031  const point& pt
1032 ) const
1033 {
1034 // return bFTreePtr_().findAnyOverlap(pt, 0.0);
1035 
1036  return (bFTreePtr_().getVolumeType(pt) == volumeType::inside);
1037 }
1038 
1039 
1041 (
1042  const List<point>& pts
1043 ) const
1044 {
1045  boolList posProc(pts.size(), true);
1046 
1047  forAll(pts, pI)
1048  {
1049  posProc[pI] = positionOnThisProcessor(pts[pI]);
1050  }
1051 
1052  return posProc;
1053 }
1054 
1055 
1057 (
1058  const treeBoundBox& box
1059 ) const
1060 {
1061 // return !procBounds().contains(box);
1062  return !bFTreePtr_().findBox(box).empty();
1063 }
1064 
1065 
1067 (
1068  const point& centre,
1069  const scalar radiusSqr
1070 ) const
1071 {
1072  // return bFTreePtr_().findAnyOverlap(centre, radiusSqr);
1073 
1074  return bFTreePtr_().findNearest(centre, radiusSqr).hit();
1075 }
1076 
1077 
1079 (
1080  const point& start,
1081  const point& end
1082 ) const
1083 {
1084  return bFTreePtr_().findLine(start, end);
1085 }
1086 
1087 
1089 (
1090  const point& start,
1091  const point& end
1092 ) const
1093 {
1094  return bFTreePtr_().findLineAny(start, end);
1095 }
1096 
1097 
1099 (
1100  const List<point>& pts
1101 ) const
1102 {
1103  DynamicList<label> toCandidateProc;
1104  DynamicList<point> testPoints;
1105  labelList ptBlockStart(pts.size(), -1);
1106  labelList ptBlockSize(pts.size(), -1);
1107 
1108  label nTotalCandidates = 0;
1109 
1110  forAll(pts, pI)
1111  {
1112  const point& pt = pts[pI];
1113 
1114  label nCandidates = 0;
1115 
1116  forAll(allBackgroundMeshBounds_, proci)
1117  {
1118  // Candidate points may lie just outside a processor box, increase
1119  // test range by using overlaps rather than contains
1120  if (allBackgroundMeshBounds_[proci].overlaps(pt, sqr(small*100)))
1121  {
1122  toCandidateProc.append(proci);
1123  testPoints.append(pt);
1124 
1125  nCandidates++;
1126  }
1127  }
1128 
1129  ptBlockStart[pI] = nTotalCandidates;
1130  ptBlockSize[pI] = nCandidates;
1131 
1132  nTotalCandidates += nCandidates;
1133  }
1134 
1135  // Needed for reverseDistribute
1136  label preDistributionToCandidateProcSize = toCandidateProc.size();
1137 
1138  autoPtr<mapDistribute> map(buildMap(toCandidateProc));
1139 
1140  map().distribute(testPoints);
1141 
1142  List<scalar> distanceSqrToCandidate(testPoints.size(), sqr(great));
1143 
1144  // Test candidate points on candidate processors
1145  forAll(testPoints, tPI)
1146  {
1147  pointIndexHit info = bFTreePtr_().findNearest
1148  (
1149  testPoints[tPI],
1150  sqr(great)
1151  );
1152 
1153  if (info.hit())
1154  {
1155  distanceSqrToCandidate[tPI] = magSqr
1156  (
1157  testPoints[tPI] - info.hitPoint()
1158  );
1159  }
1160  }
1161 
1162  map().reverseDistribute
1163  (
1164  preDistributionToCandidateProcSize,
1165  distanceSqrToCandidate
1166  );
1167 
1168  labelList ptNearestProc(pts.size(), -1);
1169 
1170  forAll(pts, pI)
1171  {
1172  // Extract the sub list of results for this point
1173 
1174  SubList<scalar> ptNearestProcResults
1175  (
1176  distanceSqrToCandidate,
1177  ptBlockSize[pI],
1178  ptBlockStart[pI]
1179  );
1180 
1181  scalar nearestProcDistSqr = great;
1182 
1183  forAll(ptNearestProcResults, pPRI)
1184  {
1185  if (ptNearestProcResults[pPRI] < nearestProcDistSqr)
1186  {
1187  nearestProcDistSqr = ptNearestProcResults[pPRI];
1188 
1189  ptNearestProc[pI] = toCandidateProc[ptBlockStart[pI] + pPRI];
1190  }
1191  }
1192 
1193  if (debug)
1194  {
1195  Pout<< pts[pI] << " nearestProcDistSqr " << nearestProcDistSqr
1196  << " ptNearestProc[pI] " << ptNearestProc[pI] << endl;
1197  }
1198 
1199  if (ptNearestProc[pI] < 0)
1200  {
1202  << "The position " << pts[pI]
1203  << " did not find a nearest point on the background mesh."
1204  << exit(FatalError);
1205  }
1206  }
1207 
1208  return ptNearestProc;
1209 }
1210 
1211 
1212 
1215 (
1216  const List<point>& starts,
1217  const List<point>& ends,
1218  bool includeOwnProcessor
1219 ) const
1220 {
1221  DynamicList<label> toCandidateProc;
1222  DynamicList<point> testStarts;
1223  DynamicList<point> testEnds;
1224  labelList segmentBlockStart(starts.size(), -1);
1225  labelList segmentBlockSize(starts.size(), -1);
1226 
1227  label nTotalCandidates = 0;
1228 
1229  forAll(starts, sI)
1230  {
1231  const point& s = starts[sI];
1232  const point& e = ends[sI];
1233 
1234  // Dummy point for treeBoundBox::intersects
1235  point p(Zero);
1236 
1237  label nCandidates = 0;
1238 
1239  forAll(allBackgroundMeshBounds_, proci)
1240  {
1241  // It is assumed that the sphere in question overlaps the source
1242  // processor, so don't test it, unless includeOwnProcessor is true
1243  if
1244  (
1245  (includeOwnProcessor || proci != Pstream::myProcNo())
1246  && allBackgroundMeshBounds_[proci].intersects(s, e, p)
1247  )
1248  {
1249  toCandidateProc.append(proci);
1250  testStarts.append(s);
1251  testEnds.append(e);
1252 
1253  nCandidates++;
1254  }
1255  }
1256 
1257  segmentBlockStart[sI] = nTotalCandidates;
1258  segmentBlockSize[sI] = nCandidates;
1259 
1260  nTotalCandidates += nCandidates;
1261  }
1262 
1263  // Needed for reverseDistribute
1264  label preDistributionToCandidateProcSize = toCandidateProc.size();
1265 
1266  autoPtr<mapDistribute> map(buildMap(toCandidateProc));
1267 
1268  map().distribute(testStarts);
1269  map().distribute(testEnds);
1270 
1271  List<pointIndexHit> segmentIntersectsCandidate(testStarts.size());
1272 
1273  // Test candidate segments on candidate processors
1274  forAll(testStarts, sI)
1275  {
1276  const point& s = testStarts[sI];
1277  const point& e = testEnds[sI];
1278 
1279  // If the sphere finds some elements of the patch, then it overlaps
1280  segmentIntersectsCandidate[sI] = bFTreePtr_().findLine(s, e);
1281  }
1282 
1283  map().reverseDistribute
1284  (
1285  preDistributionToCandidateProcSize,
1286  segmentIntersectsCandidate
1287  );
1288 
1289  List<List<pointIndexHit>> segmentHitProcs(starts.size());
1290 
1291  // Working storage for assessing processors
1292  DynamicList<pointIndexHit> tmpProcHits;
1293 
1294  forAll(starts, sI)
1295  {
1296  tmpProcHits.clear();
1297 
1298  // Extract the sub list of results for this point
1299 
1300  SubList<pointIndexHit> segmentProcResults
1301  (
1302  segmentIntersectsCandidate,
1303  segmentBlockSize[sI],
1304  segmentBlockStart[sI]
1305  );
1306 
1307  forAll(segmentProcResults, sPRI)
1308  {
1309  if (segmentProcResults[sPRI].hit())
1310  {
1311  tmpProcHits.append(segmentProcResults[sPRI]);
1312 
1313  tmpProcHits.last().setIndex
1314  (
1315  toCandidateProc[segmentBlockStart[sI] + sPRI]
1316  );
1317  }
1318  }
1319 
1320  segmentHitProcs[sI] = tmpProcHits;
1321  }
1322 
1323  return segmentHitProcs;
1324 }
1325 
1326 
1328 (
1329  const point& centre,
1330  const scalar& radiusSqr
1331 ) const
1332 {
1333  forAll(allBackgroundMeshBounds_, proci)
1334  {
1335  if (bFTreePtr_().findNearest(centre, radiusSqr).hit())
1336  {
1337  return true;
1338  }
1339  }
1340 
1341  return false;
1342 }
1343 
1344 
1346 (
1347  const point& centre,
1348  const scalar radiusSqr
1349 ) const
1350 {
1351  DynamicList<label> toProc(Pstream::nProcs());
1352 
1353  forAll(allBackgroundMeshBounds_, proci)
1354  {
1355  // Test against the bounding box of the processor
1356  if
1357  (
1358  proci != Pstream::myProcNo()
1359  && allBackgroundMeshBounds_[proci].overlaps(centre, radiusSqr)
1360  )
1361  {
1362  // Expensive test
1363 // if (bFTreePtr_().findNearest(centre, radiusSqr).hit())
1364  {
1365  toProc.append(proci);
1366  }
1367  }
1368  }
1369 
1370  return Foam::move(toProc);
1371 }
1372 
1373 
1374 //Foam::labelListList Foam::backgroundMeshDecomposition::overlapsProcessors
1375 //(
1376 // const List<point>& centres,
1377 // const List<scalar>& radiusSqrs,
1378 // const Delaunay& T,
1379 // bool includeOwnProcessor
1380 //) const
1381 //{
1382 // DynamicList<label> toCandidateProc;
1383 // DynamicList<point> testCentres;
1384 // DynamicList<scalar> testRadiusSqrs;
1385 // labelList sphereBlockStart(centres.size(), -1);
1386 // labelList sphereBlockSize(centres.size(), -1);
1387 //
1388 // label nTotalCandidates = 0;
1389 //
1390 // forAll(centres, sI)
1391 // {
1392 // const point& c = centres[sI];
1393 // scalar rSqr = radiusSqrs[sI];
1394 //
1395 // label nCandidates = 0;
1396 //
1397 // forAll(allBackgroundMeshBounds_, proci)
1398 // {
1399 // // It is assumed that the sphere in question overlaps the source
1400 // // processor, so don't test it, unless includeOwnProcessor is true
1401 // if
1402 // (
1403 // (includeOwnProcessor || proci != Pstream::myProcNo())
1404 // && allBackgroundMeshBounds_[proci].overlaps(c, rSqr)
1405 // )
1406 // {
1407 // if (bFTreePtr_().findNearest(c, rSqr).hit())
1408 // {
1409 // toCandidateProc.append(proci);
1410 // testCentres.append(c);
1411 // testRadiusSqrs.append(rSqr);
1412 //
1413 // nCandidates++;
1414 // }
1415 // }
1416 // }
1417 //
1418 // sphereBlockStart[sI] = nTotalCandidates;
1419 // sphereBlockSize[sI] = nCandidates;
1420 //
1421 // nTotalCandidates += nCandidates;
1422 // }
1423 //
1424 // // Needed for reverseDistribute
1431 //
1432 // // TODO: This is faster, but results in more vertices being referred
1433 // boolList sphereOverlapsCandidate(testCentres.size(), true);
1470 //
1476 //
1477 // labelListList sphereProcs(centres.size());
1478 //
1479 // // Working storage for assessing processors
1480 // DynamicList<label> tmpProcs;
1481 //
1482 // forAll(centres, sI)
1483 // {
1484 // tmpProcs.clear();
1485 //
1486 // // Extract the sub list of results for this point
1487 //
1488 // SubList<bool> sphereProcResults
1489 // (
1490 // sphereOverlapsCandidate,
1491 // sphereBlockSize[sI],
1492 // sphereBlockStart[sI]
1493 // );
1494 //
1495 // forAll(sphereProcResults, sPRI)
1496 // {
1497 // if (sphereProcResults[sPRI])
1498 // {
1499 // tmpProcs.append(toCandidateProc[sphereBlockStart[sI] + sPRI]);
1500 // }
1501 // }
1502 //
1503 // sphereProcs[sI] = tmpProcs;
1504 // }
1505 //
1506 // return sphereProcs;
1507 //}
1508 
1509 
1510 //Foam::labelListList Foam::backgroundMeshDecomposition::overlapProcessors
1511 //(
1512 // const point& cc,
1513 // const scalar rSqr
1514 //) const
1515 //{
1516 // DynamicList<label> toCandidateProc;
1517 // label sphereBlockStart(-1);
1518 // label sphereBlockSize(-1);
1519 //
1520 // label nCandidates = 0;
1521 //
1522 // forAll(allBackgroundMeshBounds_, proci)
1523 // {
1524 // // It is assumed that the sphere in question overlaps the source
1525 // // processor, so don't test it, unless includeOwnProcessor is true
1526 // if
1527 // (
1528 // (includeOwnProcessor || proci != Pstream::myProcNo())
1529 // && allBackgroundMeshBounds_[proci].overlaps(cc, rSqr)
1530 // )
1531 // {
1532 // toCandidateProc.append(proci);
1533 //
1534 // nCandidates++;
1535 // }
1536 // }
1537 //
1538 // sphereBlockSize = nCandidates;
1539 // nTotalCandidates += nCandidates;
1540 //
1541 // // Needed for reverseDistribute
1542 // label preDistributionToCandidateProcSize = toCandidateProc.size();
1543 //
1544 // autoPtr<mapDistribute> map(buildMap(toCandidateProc));
1545 //
1546 // map().distribute(testCentres);
1547 // map().distribute(testRadiusSqrs);
1548 //
1549 // // TODO: This is faster, but results in more vertices being referred
1551 // boolList sphereOverlapsCandidate(testCentres.size(), false);
1552 //
1553 // // Test candidate spheres on candidate processors
1554 // forAll(testCentres, sI)
1555 // {
1556 // const point& c = testCentres[sI];
1557 // const scalar rSqr = testRadiusSqrs[sI];
1558 //
1559 // const bool flagOverlap = bFTreePtr_().findNearest(c, rSqr).hit();
1560 //
1561 // if (flagOverlap)
1562 // {
1563 // // if (vertexOctree.findAnyOverlap(c, rSqr))
1568 //
1583 // sphereOverlapsCandidate[sI] = true;
1585 // }
1586 // }
1587 //
1588 // map().reverseDistribute
1589 // (
1590 // preDistributionToCandidateProcSize,
1591 // sphereOverlapsCandidate
1592 // );
1593 //
1594 // labelListList sphereProcs(centres.size());
1595 //
1596 // // Working storage for assessing processors
1597 // DynamicList<label> tmpProcs;
1598 //
1599 // forAll(centres, sI)
1600 // {
1601 // tmpProcs.clear();
1602 //
1603 // // Extract the sub list of results for this point
1604 //
1605 // SubList<bool> sphereProcResults
1606 // (
1607 // sphereOverlapsCandidate,
1608 // sphereBlockSize[sI],
1609 // sphereBlockStart[sI]
1610 // );
1611 //
1612 // forAll(sphereProcResults, sPRI)
1613 // {
1614 // if (sphereProcResults[sPRI])
1615 // {
1616 // tmpProcs.append(toCandidateProc[sphereBlockStart[sI] + sPRI]);
1617 // }
1618 // }
1619 //
1620 // sphereProcs[sI] = tmpProcs;
1621 // }
1622 //
1623 // return sphereProcs;
1624 //}
1625 
1626 
1627 // ************************************************************************* //
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
bool positionOnThisProcessor(const point &pt) const
Is the given position inside the domain of this decomposition.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual Ostream & write(const char)=0
Write character.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
List< face > faceList
Definition: faceListFwd.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static scalar & perturbTol()
Get the perturbation tolerance.
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
labelList overlapProcessors(const point &centre, const scalar radiusSqr) const
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
static autoPtr< mapDistribute > buildMap(const List< label > &toProc)
Build a mapDistribute for the supplied destination processor data.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
treeDataPrimitivePatch< bPatch > treeDataBPatch
static void exchangeSizes(const Container &sendData, labelList &sizes, const label comm=UPstream::worldComm)
Helper: exchange sizes of sendData. sendData is the data per.
Definition: exchange.C:137
labelList processorNearestPosition(const List< point > &pts) const
What is the nearest processor to the given position?
stressControl lookup("compactNormalStress") >> compactNormalStress
bool overlapsOtherProcessors(const point &centre, const scalar &radiusSqr) const
label nPoints
pointIndexHit findLine(const point &start, const point &end) const
Find nearest intersection of line between start and end, (exposing.
List< List< pointIndexHit > > intersectsProcessors(const List< point > &starts, const List< point > &ends, bool includeOwnProcessor=false) const
Which processors are intersected by the line segment, returns all.
word timeName
Definition: getTimeIndex.H:3
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:97
dimensioned< scalar > magSqr(const dimensioned< Type > &)
static const char nl
Definition: Ostream.H:260
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
pointIndexHit findLineAny(const point &start, const point &end) const
Find any intersection of line between start and end, (exposing.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
Addressing for a faceList slice.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
vector point
Point is a vector.
Definition: point.H:41
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
PrimitivePatch< faceList, const pointField > bPatch
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
autoPtr< mapDistributePolyMesh > distribute(volScalarField &cellWeights)
Redistribute the background mesh based on a supplied weight field,.
volScalarField & p
backgroundMeshDecomposition(const Time &runTime, Random &rndGen, const conformationSurfaces &geometryToConformTo, const dictionary &coeffsDict)
Construct from components in foamyHexMesh operation.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
~backgroundMeshDecomposition()
Destructor.
Namespace for OpenFOAM.
int system(const std::string &command)
Execute the specified command.
Definition: POSIX.C:1234
bool overlapsThisProcessor(const treeBoundBox &box) const
Does the given box overlap the faces of the boundary of this.