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