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