meshOctreeCreatorAdjustOctreeToSurface.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | cfMesh: A library for mesh generation
4  \\ / O peration |
5  \\ / A nd | Author: Franjo Juretic (franjo.juretic@c-fields.com)
6  \\/ M anipulation | Copyright (C) Creative Fields, Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of cfMesh.
10 
11  cfMesh is free software; you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by the
13  Free Software Foundation; either version 3 of the License, or (at your
14  option) any later version.
15 
16  cfMesh 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 cfMesh. If not, see <http://www.gnu.org/licenses/>.
23 
24 Description
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "meshOctreeCreator.H"
29 #include "triSurf.H"
30 #include "edgeMesh.H"
31 #include "boundBox.H"
32 #include "demandDrivenData.H"
33 #include "objectRefinementList.H"
34 #include "VRWGraph.H"
35 #include "meshOctreeModifier.H"
36 #include "helperFunctions.H"
37 #include "HashSet.H"
38 
39 #include "coordinateModifier.H"
42 
43 # ifdef USE_OMP
44 #include <omp.h>
45 # endif
46 
47 //#define OCTREETiming
48 //#define DEBUGSearch
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace Foam
53 {
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
58 {
59  meshOctreeModifier octreeModifier(octree_);
60  const LongList<meshOctreeCube*>& leaves = octreeModifier.leavesAccess();
61 
62  //- refine DATA boxes to the given level
63  Info << "Refining boundary boxes to the given size" << endl;
64 
65  bool changed;
66  do
67  {
68  changed = false;
69 
70  # ifdef OCTREETiming
71  const scalar startIter = omp_get_wtime();
72  # endif
73 
74  labelList refineCubes(leaves.size(), 0);
75  scalarList rThickness(leaves.size(), 0.0);
76  bool useNLayers(false);
77 
78  //- select boxes which need to be refined
79  # ifdef USE_OMP
80  # pragma omp parallel for schedule(dynamic, 100)
81  # endif
82  forAll(leaves, leafI)
83  {
84  const meshOctreeCube& oc = *leaves[leafI];
85  # ifdef DEBUGSearch
86  Info << "Checking leaf " << oc << endl;
87  Info << "Leaf has elements " << oc.hasContainedElements() << endl;
88  # endif
89 
90  if( oc.hasContainedElements() )
91  {
92  const label elRowI = oc.containedElements();
93  const VRWGraph& containedTriangles =
95 
96  bool refine(false);
97  direction minRequestedLevel(255);
98  scalar maxThicknessForLevel(0.0);
99  forAllRow(containedTriangles, elRowI, tI)
100  {
101  const label triI = containedTriangles(elRowI, tI);
102  DynList<std::pair<direction, scalar> >& refRequests =
103  surfRefLevel_[triI];
104 
105  forAllReverse(refRequests, i)
106  {
107  const std::pair<direction, scalar>& rp = refRequests[i];
108  if( rp.first <= oc.level() )
109  {
110  continue;
111  }
112 
113  if( rp.first < minRequestedLevel )
114  {
115  refine = true;
116  minRequestedLevel = rp.first;
117  maxThicknessForLevel = rp.second;
118  }
119  else if
120  (
121  (rp.first == minRequestedLevel) &&
122  (rp.second > maxThicknessForLevel)
123  )
124  {
125  maxThicknessForLevel = rp.second;
126  }
127  }
128  }
129 
130  if( refine )
131  {
132  refineCubes[leafI] = 1;
133 
134  if( maxThicknessForLevel > VSMALL )
135  {
136  rThickness[leafI] = maxThicknessForLevel;
137  useNLayers = true;
138  }
139 
140  changed = true;
141  }
142  }
143  }
144 
145  //- mark additional boxes for refinement to achieve
146  //- correct refinement distance
147  reduce(useNLayers, maxOp<label>());
148  reduce(changed, maxOp<bool>());
149  if( useNLayers && changed )
150  {
152  (
153  refineCubes,
154  rThickness
155  );
156  }
157  else if( changed )
158  {
159  //- refine boxes
160  octreeModifier.refineSelectedBoxes(refineCubes, hexRefinement_);
161  }
162 
163  # ifdef OCTREETiming
164  const scalar refTime = omp_get_wtime();
165  Info << "Time for refinement " << (refTime-startIter) << endl;
166  # endif
167 
168  if( Pstream::parRun() )
169  {
170  if( changed )
171  {
172  octreeModifier.distributeLeavesToProcessors();
173 
174  # ifdef OCTREETiming
175  const scalar distTime = omp_get_wtime();
176  Info << "Time for distributing data to processors "
177  << (distTime-refTime) << endl;
178  # endif
179 
181 
182  # ifdef OCTREETiming
183  Info << "Time for load distribution "
184  << (omp_get_wtime()-distTime) << endl;
185  # endif
186  }
187  }
188 
189  } while( changed );
190 
191  Info << "Finished refining boundary boxes" << endl;
192 }
193 
194 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195 
197 {
198  if( !meshDictPtr_ || !meshDictPtr_->found("objectRefinements") )
199  {
200  return;
201  }
202 
203  Info << "Refining boxes inside objects" << endl;
204  objectRefinementList refObjects;
205 
206  // Read objects
207  if( meshDictPtr_->isDict("objectRefinements") )
208  {
209  const dictionary& dict = meshDictPtr_->subDict("objectRefinements");
210  const wordList objectNames = dict.toc();
211 
212  refObjects.setSize(objectNames.size());
213 
214  forAll(refObjects, objectI)
215  {
216  const entry& objectEntry =
217  dict.lookupEntry(objectNames[objectI], false, false);
218 
219  refObjects.set
220  (
221  objectI,
223  (
224  objectEntry.keyword(),
225  objectEntry.dict()
226  )
227  );
228  }
229  }
230  else
231  {
232  Istream& is = meshDictPtr_->lookup("objectRefinements");
233 
234  PtrList<entry> objectEntries(is);
235  refObjects.setSize(objectEntries.size());
236 
237  forAll(refObjects, objectI)
238  {
239  refObjects.set
240  (
241  objectI,
243  (
244  objectEntries[objectI].keyword(),
245  objectEntries[objectI].dict()
246  )
247  );
248  }
249 
250  objectEntries.clear();
251  }
252 
253  coordinateModifier* cModPtr = NULL;
254 
255  if( meshDictPtr_->found("anisotropicSources") )
256  {
257  cModPtr =
259  (
260  meshDictPtr_->subDict("anisotropicSources")
261  );
262  }
263 
264  //- calculate refinement levels
265  const scalar s(readScalar(meshDictPtr_->lookup("maxCellSize")));
266 
267  List<direction> refLevels(refObjects.size(), globalRefLevel_);
268 
269  forAll(refObjects, oI)
270  {
271  //- calculate additional refinement levels from cell size
272  refObjects[oI].calculateAdditionalRefLevels(s);
273 
274  refLevels[oI] += refObjects[oI].additionalRefinementLevels();
275  }
276 
277  //- read refinement thickness
278  scalarList refThickness(refObjects.size(), 0.0);
279 
280  forAll(refThickness, oI)
281  refThickness[oI] = refObjects[oI].refinementThickness();
282 
283  forAll(refLevels, i)
284  Info << "Ref level for object " << refObjects[i].name()
285  << " is " << label(refLevels[i])
286  << " refinement thickness " << refThickness[i] << endl;
287 
288  if( octree_.neiProcs().size() )
289  forAll(refObjects, oI)
290  {
291  label l = refLevels[oI];
292  reduce(l, maxOp<label>());
293  refLevels[oI] = l;
294  }
295 
296  //- start refining boxes inside the objects
297  const boundBox& rootBox = octree_.rootBox();
298  meshOctreeModifier octreeModifier(octree_);
299  const LongList<meshOctreeCube*>& leaves = octreeModifier.leavesAccess();
300 
301  bool changed;
302  do
303  {
304  # ifdef OCTREETiming
305  const scalar startIter = omp_get_wtime();
306  # endif
307 
308  changed = false;
309 
310  labelList refineCubes(leaves.size(), 0);
311  scalarList rThickness(leaves.size(), 0.0);
312  bool useNLayers(false);
313 
314  //- select boxes which need to be refined
315  # ifdef USE_OMP
316  # pragma omp parallel for if( leaves.size() > 1000 ) \
317  schedule(dynamic, 20)
318  # endif
319  forAll(leaves, leafI)
320  {
321  const meshOctreeCube& oc = *leaves[leafI];
322 
324  continue;
325 
326  boundBox bb;
327  oc.cubeBox(rootBox, bb.min(), bb.max());
328 
329  if( cModPtr )
330  {
331  //- transform bounding boxes into non-modified coordinates
332  bb.min() = cModPtr->backwardModifiedPoint(bb.min());
333  bb.max() = cModPtr->backwardModifiedPoint(bb.max());
334  }
335 
336  bool refine(false);
337  forAll(refObjects, oI)
338  {
339  if( refObjects[oI].intersectsObject(bb) )
340  {
341  # ifdef DEBUGSearch
342  Info << "Marking leaf " << leafI
343  << " for refinement" << endl;
344  # endif
345 
346  if( oc.level() < refLevels[oI] )
347  refine = true;
348 
349  if( refThickness[oI] > VSMALL )
350  {
351  rThickness[leafI] =
352  Foam::max(rThickness[leafI], refThickness[oI]);
353 
354  useNLayers = true;
355  }
356  }
357  }
358 
359  if( refine )
360  {
361  refineCubes[leafI] = 1;
362  changed = true;
363  }
364  }
365 
366  //- mark additional boxes for refinement to achieve
367  //- correct refinement distance
368  reduce(useNLayers, maxOp<label>());
369  reduce(changed, maxOp<bool>());
370  if( useNLayers && changed )
371  {
373  (
374  refineCubes,
375  rThickness
376  );
377  }
378  else if( changed )
379  {
380  //- refine boxes
381  octreeModifier.refineSelectedBoxes(refineCubes, hexRefinement_);
382  }
383 
384  # ifdef OCTREETiming
385  const scalar refTime = omp_get_wtime();
386  Info << "Time for refinement " << (refTime-startIter) << endl;
387  # endif
388 
389  if( octree_.neiProcs().size() != 0 )
390  {
391  if( changed )
392  {
393  octreeModifier.distributeLeavesToProcessors();
394 
395  # ifdef OCTREETiming
396  const scalar distTime = omp_get_wtime();
397  Info << "Time for distributing data to processors "
398  << (distTime-refTime) << endl;
399  # endif
400 
401  loadDistribution(false);
402 
403  # ifdef OCTREETiming
404  Info << "Time for load distribution "
405  << (omp_get_wtime()-distTime) << endl;
406  # endif
407  }
408  }
409 
410  } while( changed );
411 
412  //- delete coordinate modifier if it exists
413  deleteDemandDrivenData(cModPtr);
414 
415  Info << "Finished refinement of boxes inside objects" << endl;
416 
417  //- set up inside-outside information
419 }
420 
422 {
423  if( !meshDictPtr_ || !meshDictPtr_->found("surfaceMeshRefinement") )
424  {
425  return;
426  }
427 
428  Info << "Refining boxes intersecting surface meshes" << endl;
429 
430  label nMarked;
431 
432  //- read surface meshes and calculate the refinement level for each
433  //- surface mesh
434  const dictionary& surfDict = meshDictPtr_->subDict("surfaceMeshRefinement");
435  const wordList surfaces = surfDict.toc();
436  PtrList<const triSurf> surfaceMeshesPtr(surfaces.size());
437  List<direction> refLevels(surfaces.size(), globalRefLevel_);
438  scalarList refThickness(surfaces.size(), 0.0);
439 
440  //- load surface meshes into memory
441  forAll(surfaces, surfI)
442  {
443  const dictionary& dict = surfDict.subDict(surfaces[surfI]);
444 
445  const fileName fName(dict.lookup("surfaceFile"));
446 
447  if( meshDictPtr_->found("anisotropicSources") )
448  {
449  triSurf origSurf(fName);
450  surfaceMeshGeometryModification surfMod(origSurf, *meshDictPtr_);
451 
452  surfaceMeshesPtr.set(surfI, surfMod.modifyGeometry());
453  }
454  else
455  {
456  surfaceMeshesPtr.set(surfI, new triSurf(fName));
457  }
458 
459  direction addLevel(0);
460  if( dict.found("cellSize") )
461  {
462  scalar s(readScalar(meshDictPtr_->lookup("maxCellSize")));
463 
464  const scalar cs = readScalar(dict.lookup("cellSize"));
465 
466  do
467  {
468  nMarked = 0;
469  if( cs <= s * (1.+SMALL) )
470  {
471  ++nMarked;
472  ++addLevel;
473  }
474 
475  s /= 2.0;
476 
477  } while( nMarked != 0 );
478  }
479  else if( dict.found("additionalRefinementLevels") )
480  {
481  addLevel =
482  readLabel(dict.lookup("additionalRefinementLevels"));
483  }
484 
485  if( dict.found("refinementThickness") )
486  {
487  refThickness[surfI] =
488  readScalar(dict.lookup("refinementThickness"));
489  }
490 
491  //- set the refinement level for the current surface
492  refLevels[surfI] += addLevel;
493  }
494 
495  if( octree_.neiProcs().size() )
496  forAll(refLevels, oI)
497  {
498  label l = refLevels[oI];
499  reduce(l, maxOp<label>());
500  refLevels[oI] = l;
501  }
502 
503  //- start refining boxes intersecting triangles in each refinement surface
504  const boundBox& rootBox = octree_.rootBox();
505  const vector tol = SMALL * rootBox.span();
506  meshOctreeModifier octreeModifier(octree_);
507  const LongList<meshOctreeCube*>& leaves = octreeModifier.leavesAccess();
508  DynList<label> leavesInBox;
509 
510  bool changed;
511  do
512  {
513  # ifdef OCTREETiming
514  const scalar startIter = omp_get_wtime();
515  # endif
516 
517  changed = false;
518 
519  labelList refineCubes(leaves.size(), 0);
520  labelList nLayers(leaves.size(), 0);
521  scalarField rThickness(leaves.size(), 0.0);
522  bool useNLayers(false);
523 
524  //- select boxes that need to be refined
525  forAll(surfaceMeshesPtr, surfI)
526  {
527  const triSurf& surf = surfaceMeshesPtr[surfI];
528  const pointField& points = surf.points();
529 
530  const direction surfLevel = refLevels[surfI];
531  const scalar sThickness = refThickness[surfI];
532 
533  # ifdef USE_OMP
534  # pragma omp parallel for \
535  reduction( + : nMarked) schedule(dynamic, 10) private(leavesInBox)
536  # endif
537  forAll(surf, triI)
538  {
539  //- find the bounding box of the current triangle
540  const labelledTri& tri = surf[triI];
541 
542  boundBox triBB(points[tri[0]], points[tri[0]]);
543  for(label pI=1;pI<3;++pI)
544  {
545  triBB.min() = Foam::min(triBB.min(), points[tri[pI]]);
546  triBB.max() = Foam::max(triBB.max(), points[tri[pI]]);
547  }
548 
549  triBB.min() -= tol;
550  triBB.max() += tol;
551 
552  //- find octree leaves inside the bounding box
553  leavesInBox.clear();
554  octree_.findLeavesContainedInBox(triBB, leavesInBox);
555 
556  //- check which of the leaves are intersected by the triangle
557  forAll(leavesInBox, i)
558  {
559  const label leafI = leavesInBox[i];
560 
561  const meshOctreeCube& oc = *leaves[leafI];
562 
563  if( oc.intersectsTriangleExact(surf, rootBox, triI) )
564  {
565  # ifdef DEBUGSearch
566  Info << "Marking leaf " << leafI
567  << " with coordinates " << oc
568  << " for refinement" << endl;
569  # endif
570 
571  if( oc.level() < surfLevel )
572  {
573  //- mark boxes for refinement
574  changed = true;
575  refineCubes[leafI] = 1;
576  }
577 
578  if( sThickness > VSMALL )
579  {
580  useNLayers = true;
581 
582  const scalar cs = oc.size(rootBox);
583  const label numLayers =
584  ceil(sThickness / cs);
585 
586  nLayers[leafI] =
587  Foam::max
588  (
589  nLayers[leafI],
590  Foam::max(numLayers, 1)
591  );
592 
593  rThickness[leafI] =
594  max(rThickness[leafI], sThickness);
595 
596  }
597  }
598  }
599  }
600  }
601 
602  //- mark additional boxes for refinement to achieve
603  //- correct refinement distance
604  reduce(useNLayers, maxOp<label>());
605  reduce(changed, maxOp<bool>());
606 
607  if( useNLayers && changed )
608  {
610  (
611  refineCubes,
612  rThickness
613  );
614  }
615  else if( changed )
616  {
617  //- refine boxes
618  octreeModifier.refineSelectedBoxes(refineCubes, hexRefinement_);
619  }
620 
621  # ifdef OCTREETiming
622  const scalar refTime = omp_get_wtime();
623  Info << "Time for refinement " << (refTime-startIter) << endl;
624  # endif
625 
626  if( octree_.neiProcs().size() != 0 )
627  {
628  reduce(changed, maxOp<bool>());
629  if( changed )
630  {
631  octreeModifier.distributeLeavesToProcessors();
632 
633  # ifdef OCTREETiming
634  const scalar distTime = omp_get_wtime();
635  Info << "Time for distributing data to processors "
636  << (distTime-refTime) << endl;
637  # endif
638 
639  loadDistribution(false);
640 
641  # ifdef OCTREETiming
642  Info << "Time for load distribution "
643  << (omp_get_wtime()-distTime) << endl;
644  # endif
645  }
646  }
647  } while( changed );
648 
649  Info << "Finished refinement of boxes intersecting surface meshes" << endl;
650 }
651 
653 {
654  if( !meshDictPtr_ || !meshDictPtr_->found("edgeMeshRefinement") )
655  {
656  return;
657  }
658 
659  Info << "Refining boxes intersecting edge meshes" << endl;
660 
661  //- read edge meshes and calculate the refinement level for each
662  //- edge mesh
663  const dictionary& edgeDict = meshDictPtr_->subDict("edgeMeshRefinement");
664  const wordList edgeMeshNames = edgeDict.toc();
665 
666  PtrList<const edgeMesh> edgeMeshesPtr(edgeMeshNames.size());
667  List<direction> refLevels(edgeMeshNames.size(), globalRefLevel_);
668  scalarList refThickness(edgeMeshNames.size(), 0.0);
669 
670  //- load edge meshes into memory
671  forAll(edgeMeshNames, emI)
672  {
673  const dictionary& dict = edgeDict.subDict(edgeMeshNames[emI]);
674 
675  const fileName fName(dict.lookup("edgeFile"));
676 
677  if( meshDictPtr_->found("anisotropicSources") )
678  {
679  edgeMesh origEdgeMesh(fName);
680  edgeMeshGeometryModification edgeMod(origEdgeMesh, *meshDictPtr_);
681 
682  edgeMeshesPtr.set(emI, edgeMod.modifyGeometry());
683  }
684  else
685  {
686  edgeMeshesPtr.set(emI, new edgeMesh(fName));
687  }
688 
689  label nMarked;
690  direction addLevel(0);
691  if( dict.found("cellSize") )
692  {
693  scalar s(readScalar(meshDictPtr_->lookup("maxCellSize")));
694 
695  const scalar cs = readScalar(dict.lookup("cellSize"));
696 
697  do
698  {
699  nMarked = 0;
700  if( cs <= s * (1.+SMALL) )
701  {
702  ++nMarked;
703  ++addLevel;
704  }
705 
706  s /= 2.0;
707 
708  } while( nMarked != 0 );
709  }
710  else if( dict.found("additionalRefinementLevels") )
711  {
712  addLevel =
713  readLabel(dict.lookup("additionalRefinementLevels"));
714  }
715 
716  if( dict.found("refinementThickness") )
717  {
718  refThickness[emI] =
719  readScalar(dict.lookup("refinementThickness"));
720  }
721 
722  //- set the refinement level for the current mesh
723  refLevels[emI] += addLevel;
724  }
725 
726  if( octree_.neiProcs().size() )
727  forAll(refLevels, oI)
728  {
729  label l = refLevels[oI];
730  reduce(l, maxOp<label>());
731  refLevels[oI] = l;
732  }
733 
734  //- start refining boxes intersecting triangles in each refinement surface
735  const boundBox& rootBox = octree_.rootBox();
736  const vector tol = SMALL * rootBox.span();
737  meshOctreeModifier octreeModifier(octree_);
738  const LongList<meshOctreeCube*>& leaves = octreeModifier.leavesAccess();
739  DynList<label> leavesInBox, intersectedLeaves;
740 
741  bool changed;
742  do
743  {
744  changed = false;
745 
746  # ifdef OCTREETiming
747  const scalar startIter = omp_get_wtime();
748  # endif
749 
750  labelList refineCubes(leaves.size(), 0);
751  scalarList rThickness(leaves.size(), 0.0);
752  bool useNLayers(false);
753 
754  //- select boxes which need to be refined
755  forAll(edgeMeshNames, emI)
756  {
757  const edgeMesh& edgeMesh = edgeMeshesPtr[emI];
758  const pointField& points = edgeMesh.points();
759  const edgeList& edges = edgeMesh.edges();
760 
761  # ifdef USE_OMP
762  # pragma omp parallel for schedule(dynamic, 10) \
763  private(leavesInBox,intersectedLeaves)
764  # endif
765  forAll(edges, edgeI)
766  {
767  //- find the bounding box of the current triangle
768  const edge& e = edges[edgeI];
769 
770  const point& sp = points[e.start()];
771  const point& ep = points[e.end()];
772 
773  boundBox edgeBB(min(sp, ep), max(sp, ep));
774 
775  edgeBB.min() -= tol;
776  edgeBB.max() += tol;
777 
778  //- find octree leaves inside the bounding box
779  leavesInBox.clear();
780  octree_.findLeavesContainedInBox(edgeBB, leavesInBox);
781 
782  //- check which of the leaves are intersected by the triangle
783  intersectedLeaves.clear();
784  forAll(leavesInBox, i)
785  {
786  const label leafI = leavesInBox[i];
787 
788  const meshOctreeCube& oc = *leaves[leafI];
789 
790  if( oc.intersectsLine(rootBox, sp, ep) )
791  {
792  intersectedLeaves.append(leafI);
793 
794  if( oc.level() < refLevels[emI] )
795  {
796  # ifdef DEBUGSearch
797  Info << "Marking leaf " << leafI
798  << " with coordinates " << oc
799  << " for refinement" << endl;
800  # endif
801 
802  refineCubes[leafI] = 1;
803  changed = true;
804  }
805  }
806  }
807 
808  if( refThickness[emI] > VSMALL )
809  {
810  useNLayers = true;
811 
812  forAll(intersectedLeaves, i)
813  {
814  const label leafI = intersectedLeaves[i];
815 
816  rThickness[leafI] =
817  Foam::max(rThickness[leafI], refThickness[emI]);
818  }
819  }
820  }
821  }
822 
823  //- mark additional boxes for refinement to achieve
824  //- correct refinement distance
825  reduce(useNLayers, maxOp<label>());
826  reduce(changed, maxOp<bool>());
827  if( useNLayers && changed )
828  {
830  (
831  refineCubes,
832  rThickness
833  );
834  }
835  else if( changed )
836  {
837  //- refine boxes
838  octreeModifier.refineSelectedBoxes(refineCubes, hexRefinement_);
839  }
840 
841  # ifdef OCTREETiming
842  const scalar refTime = omp_get_wtime();
843  Info << "Time for refinement " << (refTime-startIter) << endl;
844  # endif
845 
846  if( octree_.neiProcs().size() != 0 )
847  {
848  if( changed )
849  {
850  octreeModifier.distributeLeavesToProcessors();
851 
852  # ifdef OCTREETiming
853  const scalar distTime = omp_get_wtime();
854  Info << "Time for distributing data to processors "
855  << (distTime-refTime) << endl;
856  # endif
857 
858  loadDistribution(false);
859 
860  # ifdef OCTREETiming
861  Info << "Time for load distribution "
862  << (omp_get_wtime()-distTime) << endl;
863  # endif
864  }
865  }
866  } while( changed );
867 
868  Info << "Finished refinement of boxes intersecting edge meshes" << endl;
869 }
870 
872 {
873  # ifdef OCTREETiming
874  const scalar startTime = omp_get_wtime();
875  # endif
876 
879 
880  Info << "Refining boxes near DATA boxes" << endl;
881 
882  //- ensure one to one matching with unknown boxes
883  meshOctreeModifier octreeModifier(octree_);
884  const LongList<meshOctreeCube*>& leaves = octreeModifier.leavesAccess();
885 
886  # ifdef DEBUGSearch
887  forAll(leaves, leafI)
888  Info << "Leaf " << leafI << " is " << *leaves[leafI]
889  << " type " << label(leaves[leafI]->cubeType()) << endl;
890  # endif
891 
892  labelList refineBox(leaves.size(), 0);
893 
894  labelHashSet transferCoordinates;
895  LongList<meshOctreeCubeCoordinates> checkCoordinates;
896 
897  # ifdef USE_OMP
898  # pragma omp parallel for if( leaves.size() > 1000 ) \
899  schedule(dynamic, 20)
900  # endif
901  forAll(leaves, leafI)
902  {
903  if( leaves[leafI]->hasContainedElements() )
904  {
905  const meshOctreeCube& oc = *leaves[leafI];
906 
907  # ifdef DEBUGSearch
908  Info << "Refining boxes near box " << leafI
909  << " with coordinates " << oc << endl;
910  # endif
911 
912  for(label k=0;k<26;++k)
913  {
914  const label neiLabel =
916  (
917  oc.coordinates() + rp[k]
918  );
919 
920  if( neiLabel == meshOctreeCube::OTHERPROC )
921  {
922  # ifdef USE_OMP
923  # pragma omp critical
924  # endif
925  {
926  if( !transferCoordinates.found(leafI) )
927  {
928  transferCoordinates.insert(leafI);
929  checkCoordinates.append(oc.coordinates());
930  }
931  }
932 
933  continue;
934  }
935 
936  if( neiLabel == -1 )
937  continue;
938 
939  const meshOctreeCube* nei = leaves[neiLabel];
940  if( nei->level() == oc.level() )
941  continue;
943  continue;
944 
945  # ifdef DEBUGSearch
946  Info << "Adding neighbour " << *nei << " of type"
947  << label(nei->cubeType()) << endl;
948  # endif
949 
950  refineBox[nei->cubeLabel()] = 1;
951  }
952  }
953  }
954 
955  if( octree_.neiProcs().size() )
956  {
957  LongList<meshOctreeCubeCoordinates> receivedCoordinates;
959  (
960  checkCoordinates,
961  receivedCoordinates
962  );
963 
964  # ifdef USE_OMP
965  # pragma omp parallel for if( receivedCoordinates.size() > 1000 ) \
966  schedule(dynamic, 20)
967  # endif
968  forAll(receivedCoordinates, ccI)
969  {
970  const meshOctreeCubeCoordinates& cc = receivedCoordinates[ccI];
971  for(label k=0;k<26;++k)
972  {
973  const label neiLabel =
975 
976  if( neiLabel < 0 )
977  continue;
978 
979  const meshOctreeCube* nei = leaves[neiLabel];
980  if( nei->level() == cc.level() )
981  continue;
983  continue;
984 
985  # ifdef DEBUGSearch
986  Info << "Adding neighbour " << *nei << " of type"
987  << label(nei->cubeType()) << endl;
988  # endif
989 
990  refineBox[nei->cubeLabel()] = 1;
991  }
992  }
993  }
994 
995  for(label i=1;i<nLayers;i++)
996  {
997  if( Pstream::parRun() )
998  {
999  checkCoordinates.clear();
1000  transferCoordinates.clear();
1001  }
1002 
1003  # ifdef USE_OMP
1004  # pragma omp parallel for if( leaves.size() > 1000 ) \
1005  schedule(dynamic, 20)
1006  # endif
1007  forAll(leaves, leafI)
1008  {
1009  if( refineBox[leafI] == i )
1010  {
1011  const meshOctreeCube& oc = *leaves[leafI];
1012 
1013  # ifdef DEBUGSearch
1014  Info << "Refining boxes near box " << leafI
1015  << " with coordinates " << oc << endl;
1016  # endif
1017 
1018  for(label k=0;k<26;++k)
1019  {
1020  const label neiLabel =
1022  (
1023  oc.coordinates() + rp[k]
1024  );
1025 
1026  if( neiLabel == meshOctreeCube::OTHERPROC )
1027  {
1028  # ifdef USE_OMP
1029  # pragma omp critical
1030  # endif
1031  {
1032  if( !transferCoordinates.found(leafI) )
1033  {
1034  transferCoordinates.insert(leafI);
1035  checkCoordinates.append(oc.coordinates());
1036  }
1037  }
1038 
1039  continue;
1040  }
1041 
1042  if( neiLabel == -1 )
1043  continue;
1044 
1045  const meshOctreeCube* nei = leaves[neiLabel];
1046  if( nei->level() == oc.level() )
1047  continue;
1048  if( nei->cubeType() & meshOctreeCubeBasic::OUTSIDE )
1049  continue;
1050 
1051  # ifdef DEBUGSearch
1052  Info << "Layer " << label(i) << endl;
1053  Info << "Adding neighbour " << *nei << " of type"
1054  << label(nei->cubeType()) << endl;
1055  # endif
1056 
1057  refineBox[nei->cubeLabel()] = i + 1;
1058  }
1059  }
1060  }
1061 
1062  if( octree_.neiProcs().size() )
1063  {
1064  LongList<meshOctreeCubeCoordinates> receivedCoordinates;
1066  (
1067  checkCoordinates,
1068  receivedCoordinates
1069  );
1070 
1071  # ifdef USE_OMP
1072  # pragma omp parallel for if( receivedCoordinates.size() > 1000 ) \
1073  schedule(dynamic, 20)
1074  # endif
1075  forAll(receivedCoordinates, ccI)
1076  {
1077  const meshOctreeCubeCoordinates& cc = receivedCoordinates[ccI];
1078 
1079  for(label k=0;k<26;++k)
1080  {
1081  const label neiLabel =
1083 
1084  if( neiLabel < 0 )
1085  continue;
1086 
1087  const meshOctreeCube* nei = leaves[neiLabel];
1088  if( nei->level() == cc.level() )
1089  continue;
1090  if( nei->cubeType() & meshOctreeCubeBasic::OUTSIDE )
1091  continue;
1092 
1093  # ifdef DEBUGSearch
1094  Info << "Adding neighbour " << *nei << " of type"
1095  << label(nei->cubeType()) << endl;
1096  # endif
1097 
1098  refineBox[nei->cubeLabel()] = i + 1;
1099  }
1100  }
1101  }
1102  }
1103 
1104  //- refine cubes
1105  octreeModifier.refineSelectedBoxes(refineBox, hexRefinement_);
1106 
1107  # ifdef OCTREETiming
1108  const scalar refTime = omp_get_wtime();
1109  Info << "Time for refinement " << (refTime-startTime) << endl;
1110  # endif
1111 
1112  if( Pstream::parRun() )
1113  {
1114  octreeModifier.distributeLeavesToProcessors();
1115  loadDistribution();
1116  }
1117 
1118  # ifdef OCTREETiming
1119  const scalar distTime = omp_get_wtime();
1120  Info << "Time for load balancing " << (distTime-refTime) << endl;
1121  # endif
1122 
1123  //- set up inside-outside information
1125 
1126  # ifdef OCTREETiming
1127  Info << "Time for calculation of inside/outside information "
1128  << (omp_get_wtime()-distTime) << endl;
1129  # endif
1130 
1131  Info << "Finished refining boxes near DATA boxes" << endl;
1132 }
1133 
1136  const direction refLevel,
1137  const direction cubeType
1138 )
1139 {
1140  label nRefined;
1141  meshOctreeModifier octreeMod(octree_);
1142 
1143  do
1144  {
1145  # ifdef OCTREETiming
1146  const scalar startIter = omp_get_wtime();
1147  # endif
1148 
1149  nRefined = 0;
1150 
1151  const LongList<meshOctreeCube*>& leaves = octreeMod.leavesAccess();
1152 
1153  labelList refineCubes(leaves.size(), direction(0));
1154 
1155  # ifdef USE_OMP
1156  # pragma omp parallel for if( leaves.size() > 1000 ) \
1157  reduction(+ : nRefined) schedule(dynamic, 20)
1158  # endif
1159  forAll(leaves, leafI)
1160  {
1161  const meshOctreeCube& oc = *leaves[leafI];
1162 
1163  if( (oc.cubeType() & cubeType) && (oc.level() < refLevel) )
1164  {
1165  refineCubes[leafI] = 1;
1166  ++nRefined;
1167  }
1168  }
1169 
1170  //- refine selected cubes
1171  octreeMod.refineSelectedBoxes(refineCubes, hexRefinement_);
1172 
1173  # ifdef OCTREETiming
1174  const scalar refTime = omp_get_wtime();
1175  Info << "Time for refinement " << (refTime-startIter) << endl;
1176  # endif
1177 
1178  if( Pstream::parRun() )
1179  {
1180  reduce(nRefined, sumOp<label>());
1181  if( nRefined )
1182  {
1183  octreeMod.distributeLeavesToProcessors();
1184 
1185  # ifdef OCTREETiming
1186  const scalar distTime = omp_get_wtime();
1187  Info << "Time for distributing data to processors "
1188  << (distTime-refTime) << endl;
1189  # endif
1190 
1191  loadDistribution();
1192 
1193  # ifdef OCTREETiming
1194  Info << "Time for load distribution "
1195  << (omp_get_wtime()-distTime) << endl;
1196  # endif
1197  }
1198  }
1199 
1200  } while( nRefined != 0 );
1201 }
1202 
1203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1204 
1205 } // End namespace Foam
1206 
1207 // ************************************************************************* //
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:65
Foam::meshOctreeCreator::refineBoxesIntersectingSurfaces
void refineBoxesIntersectingSurfaces()
Definition: meshOctreeCreatorAdjustOctreeToSurface.C:421
Foam::maxOp
Definition: ops.H:172
Foam::LongList::append
void append(const T &e)
Append an element at the end of the list.
Definition: LongListI.H:265
triSurf.H
Foam::meshOctree::regularityPositions
const FixedList< meshOctreeCubeCoordinates, 26 > & regularityPositions() const
Definition: meshOctreeI.H:141
Foam::meshOctreeCubeBasic::OUTSIDE
@ OUTSIDE
Definition: meshOctreeCubeBasic.H:89
Foam::meshOctreeModifier::refineSelectedBoxesAndAdditionalLayers
void refineSelectedBoxesAndAdditionalLayers(labelList &refineBox, const scalarList &refThickness)
refine selected boxes and the boxes within the given range
Definition: meshOctreeModifierRefineSelectedBoxes.C:500
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
Foam::entry::keyword
const keyType & keyword() const
Return keyword.
Definition: entry.H:120
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::meshOctreeCubeCoordinates
Definition: meshOctreeCubeCoordinates.H:55
VRWGraph.H
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::triSurfPoints::points
const pointField & points() const
access to points
Definition: triSurfPointsI.H:44
objectRefinementList.H
Foam::meshOctreeCube::cubeLabel
label cubeLabel() const
position of the cube in the list of leaves
Definition: meshOctreeCubeI.H:76
Foam::meshOctreeCube::hasContainedElements
bool hasContainedElements() const
return true if the box contains some triangles
Definition: meshOctreeCubeI.H:81
Foam::LongList::clear
void clear()
Clear the list, i.e. set next free to zero.
Definition: LongListI.H:230
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::meshOctreeModifier::distributeLeavesToProcessors
void distributeLeavesToProcessors()
Definition: meshOctreeModifierDistributeLeavesToProcessors.C:39
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::meshOctree::neiProcs
const labelList & neiProcs() const
neighbour processors of the current one
Definition: meshOctreeI.H:152
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::meshOctreeCubeCoordinates::intersectsLine
bool intersectsLine(const boundBox &, const point &, const point &) const
check if the cube intersects a line
Definition: meshOctreeCubeCoordinatesIntersections.C:309
Foam::meshOctreeCubeBasic::cubeType
direction cubeType() const
return type
Definition: meshOctreeCubeBasicI.H:75
Foam::LongList::size
label size() const
Size of the active part of the list.
Definition: LongListI.H:203
Foam::HashSet< label, Hash< label > >
Foam::meshOctree::findLeavesContainedInBox
void findLeavesContainedInBox(const boundBox &, DynList< const meshOctreeCube *, 256 > &) const
find leaves contained in the box
Definition: meshOctreeInsideCalculations.C:131
Foam::meshOctreeCubeBasic::coordinates
const meshOctreeCubeCoordinates & coordinates() const
return coordinates in the octree
Definition: meshOctreeCubeBasicI.H:90
Foam::dictionary::isDict
bool isDict(const word &) const
Check if entry is a sub-dictionary.
Definition: dictionary.C:600
Foam::meshOctreeModifier
Definition: meshOctreeModifier.H:48
meshOctreeModifier.H
Foam::edgeMeshGeometryModification::modifyGeometry
const edgeMesh * modifyGeometry() const
modify coordinates
Definition: edgeMeshGeometryModification.C:77
Foam::meshOctreeModifier::refineSelectedBoxes
void refineSelectedBoxes(labelList &refineBox, const bool hexRefinement=false)
Definition: meshOctreeModifierRefineSelectedBoxes.C:428
Foam::boundBox::span
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:84
Foam::LongList
Definition: LongList.H:55
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:40
Foam::meshOctreeCreator::refineBoxes
void refineBoxes(const direction refLevel, const direction cubeType)
refine boxes of the given flag to the given size
Definition: meshOctreeCreatorAdjustOctreeToSurface.C:1135
Foam::meshOctree::findLeafLabelForPosition
label findLeafLabelForPosition(const meshOctreeCubeCoordinates &) const
return leaf cube for the given position
Definition: meshOctreeNeighbourSearches.C:516
Foam::meshOctreeCreator::createInsideOutsideInformation
void createInsideOutsideInformation()
Definition: meshOctreeCreatorFrontalMarking.C:40
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::PtrList::set
bool set(const label) const
Is element set.
Foam::meshOctreeCreator::octree_
meshOctree & octree_
Reference to meshOctree.
Definition: meshOctreeCreator.H:63
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::dictionary::found
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:304
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
edgeMesh.H
Foam::Info
messageStream Info
Foam::meshOctree::exchangeRequestsWithNeighbourProcessors
void exchangeRequestsWithNeighbourProcessors(const LongList< meshOctreeCubeCoordinates > &dataToSend, LongList< meshOctreeCubeCoordinates > &dataToReceive) const
exchange requests with other processors generating the octree
Definition: meshOctreeParallelCommunication.C:41
coordinateModifier.H
Foam::meshOctreeCreator::refineBoxesNearDataBoxes
void refineBoxesNearDataBoxes(const label nLayers=1)
refine boxes near DATA boxes to get a nice smooth surface
Definition: meshOctreeCreatorAdjustOctreeToSurface.C:871
Foam::meshOctreeCubeCoordinates::intersectsTriangleExact
bool intersectsTriangleExact(const triSurf &, const boundBox &, const label) const
Definition: meshOctreeCubeCoordinatesIntersections.C:187
Foam::edgeMesh::points
const pointField & points() const
Return points.
Definition: edgeMeshI.H:39
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
forAllRow
#define forAllRow(graph, rowI, index)
Definition: VRWGraph.H:277
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
Foam::meshOctreeCube::slotPtr
const meshOctreeSlot * slotPtr() const
return the pointer to the slot containing the cube
Definition: meshOctreeCubeI.H:58
Foam::edgeMesh::edges
const edgeList & edges() const
Return edges.
Definition: edgeMeshI.H:45
HashSet.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::surfaceMeshGeometryModification
Definition: surfaceMeshGeometryModification.H:55
forAllReverse
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:418
Foam::meshOctreeCreator::globalRefLevel_
direction globalRefLevel_
ref level to achieve max cell size
Definition: meshOctreeCreator.H:108
Foam::HashTable::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
Foam::PtrList::clear
void clear()
Clear the PtrList, i.e. set size to zero deleting all the.
Definition: PtrList.C:185
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::entry::dict
virtual const dictionary & dict() const =0
Return dictionary if this entry is a dictionary.
Foam::meshOctreeCreator::loadDistribution
void loadDistribution(const bool distributeUsed=false)
Definition: meshOctreeCreatorLoadDistribution.C:46
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::DynList
Definition: DynList.H:53
Foam::meshOctreeCubeBasic::OTHERPROC
@ OTHERPROC
Definition: meshOctreeCubeBasic.H:93
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::meshOctreeCubeCoordinates::size
scalar size(const boundBox &) const
return size
Definition: meshOctreeCubeCoordinatesI.H:213
Foam::meshOctreeCreator::surfRefLevel_
List< DynList< std::pair< direction, scalar > > > surfRefLevel_
this list contains ref level for each surface triangle
Definition: meshOctreeCreator.H:111
boundBox.H
Foam::meshOctreeSlot::containedTriangles_
VRWGraph containedTriangles_
surface triangles contained in an octree cube
Definition: meshOctreeSlot.H:62
Foam::meshOctreeCreator::meshDictPtr_
const IOdictionary * meshDictPtr_
Dictionary containing information necessary to perform refinement.
Definition: meshOctreeCreator.H:69
Foam::meshOctree::rootBox
const boundBox & rootBox() const
return rootBox
Definition: meshOctreeI.H:135
Foam::sumOp
Definition: ops.H:162
Foam::meshOctreeCubeCoordinates::level
direction level() const
return level
Definition: meshOctreeCubeCoordinatesI.H:74
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::HashTable::clear
void clear()
Clear all entries from table.
Definition: HashTable.C:473
helperFunctions.H
Foam::Vector< scalar >
Foam::surfaceMeshGeometryModification::modifyGeometry
const triSurf * modifyGeometry() const
modify coordinates
Definition: surfaceMeshGeometryModification.C:77
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
Foam::meshOctreeCreator::refineBoxesIntersectingEdgeMeshes
void refineBoxesIntersectingEdgeMeshes()
Definition: meshOctreeCreatorAdjustOctreeToSurface.C:652
Foam::objectRefinement::New
static autoPtr< objectRefinement > New(const word &name, const dictionary &dict)
Select constructed from dictionary.
Definition: newObjectRefinement.C:38
Foam::FixedList
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:53
Foam::meshOctreeCube
Definition: meshOctreeCube.H:56
points
const pointField & points
Definition: gmvOutputHeader.H:1
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
Foam::PtrList::setSize
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:142
Foam::readLabel
label readLabel(Istream &is)
Definition: label.H:64
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::labelledTri
Triangle with additional region number.
Definition: labelledTri.H:49
Foam::edgeMeshGeometryModification
Definition: edgeMeshGeometryModification.H:55
Foam::boundBox
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:55
edgeMeshGeometryModification.H
Foam::meshOctreeCreator::refineBoxesContainedInObjects
void refineBoxesContainedInObjects()
refine boxes contained inside the objects for refinement
Definition: meshOctreeCreatorAdjustOctreeToSurface.C:196
startTime
Foam::label startTime
Definition: checkTimeOptions.H:5
Foam::meshOctreeCube::containedElements
label containedElements() const
Definition: meshOctreeCubeI.H:89
meshOctreeCreator.H
Foam::dictionary::toc
wordList toc() const
Return the table of contents.
Definition: dictionary.C:697
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::meshOctreeCreator::refineBoundary
void refineBoundary()
Definition: meshOctreeCreatorAdjustOctreeToSurface.C:57
Foam::coordinateModifier
Definition: coordinateModifier.H:54
Foam::meshOctreeCubeCoordinates::cubeBox
void cubeBox(const boundBox &, point &, point &) const
return min and max points
Definition: meshOctreeCubeCoordinatesI.H:174
Foam::coordinateModifier::backwardModifiedPoint
point backwardModifiedPoint(const point &) const
calculate the displacement vector for the backward modification
Definition: coordinateModifier.C:205
Foam::dictionary::subDict
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:631
Foam::VRWGraph
Definition: VRWGraph.H:101
Foam::triSurf
Definition: triSurf.H:59
Foam::edgeMesh
Points connected by edges.
Definition: edgeMesh.H:69
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::DynList::clear
void clear()
Clear the list, i.e. set next free to zero.
Definition: DynListI.H:279
Foam::meshOctreeCreator::hexRefinement_
bool hexRefinement_
hex refinement flag
Definition: meshOctreeCreator.H:72
Foam::meshOctreeModifier::leavesAccess
LongList< meshOctreeCube * > & leavesAccess()
return leaves
Definition: meshOctreeModifierI.H:95
surfaceMeshGeometryModification.H
Foam::DynList::append
void append(const T &e)
Append an element at the end of the list.
Definition: DynListI.H:304