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