snappyHexMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  snappyHexMesh
29 
30 Group
31  grpMeshGenerationUtilities
32 
33 Description
34  Automatic split hex mesher. Refines and snaps to surface.
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #include "argList.H"
39 #include "Time.H"
40 #include "fvMesh.H"
41 #include "snappyRefineDriver.H"
42 #include "snappySnapDriver.H"
43 #include "snappyLayerDriver.H"
44 #include "searchableSurfaces.H"
45 #include "refinementSurfaces.H"
46 #include "refinementFeatures.H"
47 #include "shellSurfaces.H"
48 #include "decompositionMethod.H"
49 #include "fvMeshDistribute.H"
50 #include "wallPolyPatch.H"
51 #include "refinementParameters.H"
52 #include "snapParameters.H"
53 #include "layerParameters.H"
54 #include "vtkSetWriter.H"
55 #include "faceSet.H"
56 #include "motionSmoother.H"
57 #include "polyTopoChange.H"
59 #include "surfZoneIdentifierList.H"
60 #include "UnsortedMeshedSurface.H"
61 #include "MeshedSurface.H"
62 #include "globalIndex.H"
63 #include "IOmanip.H"
64 #include "decompositionModel.H"
65 #include "fvMeshTools.H"
66 #include "profiling.H"
67 #include "processorMeshes.H"
68 #include "snappyVoxelMeshDriver.H"
69 
70 using namespace Foam;
71 
72 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73 
74 // Convert size (as fraction of defaultCellSize) to refinement level
75 label sizeCoeffToRefinement
76 (
77  const scalar level0Coeff, // ratio of hex cell size v.s. defaultCellSize
78  const scalar sizeCoeff
79 )
80 {
81  return round(::log(level0Coeff/sizeCoeff)/::log(2));
82 }
83 
84 
85 autoPtr<refinementSurfaces> createRefinementSurfaces
86 (
87  const searchableSurfaces& allGeometry,
88  const dictionary& surfacesDict,
89  const dictionary& shapeControlDict,
90  const label gapLevelIncrement,
91  const scalar level0Coeff
92 )
93 {
94  autoPtr<refinementSurfaces> surfacePtr;
95 
96  // Count number of surfaces.
97  label surfi = 0;
98  forAll(allGeometry.names(), geomi)
99  {
100  const word& geomName = allGeometry.names()[geomi];
101 
102  if (surfacesDict.found(geomName))
103  {
104  surfi++;
105  }
106  }
107 
108  labelList surfaces(surfi);
109  wordList names(surfi);
110  PtrList<surfaceZonesInfo> surfZones(surfi);
111 
112  labelList regionOffset(surfi);
113 
114  labelList globalMinLevel(surfi, Zero);
115  labelList globalMaxLevel(surfi, Zero);
116  labelList globalLevelIncr(surfi, Zero);
117  PtrList<dictionary> globalPatchInfo(surfi);
118  List<Map<label>> regionMinLevel(surfi);
119  List<Map<label>> regionMaxLevel(surfi);
120  List<Map<label>> regionLevelIncr(surfi);
121  List<Map<scalar>> regionAngle(surfi);
122  List<Map<autoPtr<dictionary>>> regionPatchInfo(surfi);
123 
124  wordHashSet unmatchedKeys(surfacesDict.toc());
125 
126  surfi = 0;
127  forAll(allGeometry.names(), geomi)
128  {
129  const word& geomName = allGeometry.names()[geomi];
130 
131  const entry* ePtr = surfacesDict.findEntry(geomName, keyType::REGEX);
132 
133  if (ePtr)
134  {
135  const dictionary& shapeDict = ePtr->dict();
136  unmatchedKeys.erase(ePtr->keyword());
137 
138  names[surfi] = geomName;
139  surfaces[surfi] = geomi;
140 
141  const searchableSurface& surface = allGeometry[geomi];
142 
143  // Find the index in shapeControlDict
144  // Invert surfaceCellSize to get the refinementLevel
145 
146  const word scsFuncName =
147  shapeDict.get<word>("surfaceCellSizeFunction");
148 
149  const dictionary& scsDict =
150  shapeDict.optionalSubDict(scsFuncName + "Coeffs");
151 
152  const scalar surfaceCellSize =
153  scsDict.get<scalar>("surfaceCellSizeCoeff");
154 
155  const label refLevel = sizeCoeffToRefinement
156  (
157  level0Coeff,
158  surfaceCellSize
159  );
160 
161  globalMinLevel[surfi] = refLevel;
162  globalMaxLevel[surfi] = refLevel;
163  globalLevelIncr[surfi] = gapLevelIncrement;
164 
165  // Surface zones
166  surfZones.set
167  (
168  surfi,
169  new surfaceZonesInfo
170  (
171  surface,
172  shapeDict,
173  allGeometry.regionNames()[surfaces[surfi]]
174  )
175  );
176 
177 
178  // Global perpendicular angle
179  if (shapeDict.found("patchInfo"))
180  {
181  globalPatchInfo.set
182  (
183  surfi,
184  shapeDict.subDict("patchInfo").clone()
185  );
186  }
187 
188 
189  // Per region override of patchInfo
190 
191  if (shapeDict.found("regions"))
192  {
193  const dictionary& regionsDict = shapeDict.subDict("regions");
194  const wordList& regionNames =
195  allGeometry[surfaces[surfi]].regions();
196 
197  forAll(regionNames, regioni)
198  {
199  if (regionsDict.found(regionNames[regioni]))
200  {
201  // Get the dictionary for region
202  const dictionary& regionDict = regionsDict.subDict
203  (
204  regionNames[regioni]
205  );
206 
207  if (regionDict.found("patchInfo"))
208  {
209  regionPatchInfo[surfi].insert
210  (
211  regioni,
212  regionDict.subDict("patchInfo").clone()
213  );
214  }
215  }
216  }
217  }
218 
219  // Per region override of cellSize
220  if (shapeDict.found("regions"))
221  {
222  const dictionary& shapeControlRegionsDict =
223  shapeDict.subDict("regions");
224  const wordList& regionNames =
225  allGeometry[surfaces[surfi]].regions();
226 
227  forAll(regionNames, regioni)
228  {
229  if (shapeControlRegionsDict.found(regionNames[regioni]))
230  {
231  const dictionary& shapeControlRegionDict =
232  shapeControlRegionsDict.subDict
233  (
234  regionNames[regioni]
235  );
236 
237  const word scsFuncName =
238  shapeControlRegionDict.get<word>
239  (
240  "surfaceCellSizeFunction"
241  );
242  const dictionary& scsDict =
243  shapeControlRegionDict.subDict
244  (
245  scsFuncName + "Coeffs"
246  );
247 
248  const scalar surfaceCellSize =
249  scsDict.get<scalar>("surfaceCellSizeCoeff");
250 
251  const label refLevel = sizeCoeffToRefinement
252  (
253  level0Coeff,
254  surfaceCellSize
255  );
256 
257  regionMinLevel[surfi].insert(regioni, refLevel);
258  regionMaxLevel[surfi].insert(regioni, refLevel);
259  regionLevelIncr[surfi].insert(regioni, 0);
260  }
261  }
262  }
263 
264  surfi++;
265  }
266  }
267 
268  // Calculate local to global region offset
269  label nRegions = 0;
270 
271  forAll(surfaces, surfi)
272  {
273  regionOffset[surfi] = nRegions;
274  nRegions += allGeometry[surfaces[surfi]].regions().size();
275  }
276 
277  // Rework surface specific information into information per global region
278  labelList minLevel(nRegions, Zero);
279  labelList maxLevel(nRegions, Zero);
280  labelList gapLevel(nRegions, -1);
281  PtrList<dictionary> patchInfo(nRegions);
282 
283  forAll(globalMinLevel, surfi)
284  {
285  label nRegions = allGeometry[surfaces[surfi]].regions().size();
286 
287  // Initialise to global (i.e. per surface)
288  for (label i = 0; i < nRegions; i++)
289  {
290  label globalRegioni = regionOffset[surfi] + i;
291  minLevel[globalRegioni] = globalMinLevel[surfi];
292  maxLevel[globalRegioni] = globalMaxLevel[surfi];
293  gapLevel[globalRegioni] =
294  maxLevel[globalRegioni]
295  + globalLevelIncr[surfi];
296 
297  if (globalPatchInfo.set(surfi))
298  {
299  patchInfo.set
300  (
301  globalRegioni,
302  globalPatchInfo[surfi].clone()
303  );
304  }
305  }
306 
307  // Overwrite with region specific information
308  forAllConstIters(regionMinLevel[surfi], iter)
309  {
310  label globalRegioni = regionOffset[surfi] + iter.key();
311 
312  minLevel[globalRegioni] = iter();
313  maxLevel[globalRegioni] = regionMaxLevel[surfi][iter.key()];
314  gapLevel[globalRegioni] =
315  maxLevel[globalRegioni]
316  + regionLevelIncr[surfi][iter.key()];
317  }
318 
319  const Map<autoPtr<dictionary>>& localInfo = regionPatchInfo[surfi];
320  forAllConstIters(localInfo, iter)
321  {
322  label globalRegioni = regionOffset[surfi] + iter.key();
323  patchInfo.set(globalRegioni, iter()().clone());
324  }
325  }
326 
327  surfacePtr.reset
328  (
330  (
331  allGeometry,
332  surfaces,
333  names,
334  surfZones,
335  regionOffset,
336  minLevel,
337  maxLevel,
338  gapLevel,
339  scalarField(nRegions, -GREAT), //perpendicularAngle,
340  patchInfo,
341  false //dryRun
342  )
343  );
344 
345 
346  const refinementSurfaces& rf = surfacePtr();
347 
348  // Determine maximum region name length
349  label maxLen = 0;
350  forAll(rf.surfaces(), surfi)
351  {
352  label geomi = rf.surfaces()[surfi];
353  const wordList& regionNames = allGeometry.regionNames()[geomi];
354  forAll(regionNames, regioni)
355  {
356  maxLen = Foam::max(maxLen, label(regionNames[regioni].size()));
357  }
358  }
359 
360 
361  Info<< setw(maxLen) << "Region"
362  << setw(10) << "Min Level"
363  << setw(10) << "Max Level"
364  << setw(10) << "Gap Level" << nl
365  << setw(maxLen) << "------"
366  << setw(10) << "---------"
367  << setw(10) << "---------"
368  << setw(10) << "---------" << endl;
369 
370  forAll(rf.surfaces(), surfi)
371  {
372  label geomi = rf.surfaces()[surfi];
373 
374  Info<< rf.names()[surfi] << ':' << nl;
375 
376  const wordList& regionNames = allGeometry.regionNames()[geomi];
377 
378  forAll(regionNames, regioni)
379  {
380  label globali = rf.globalRegion(surfi, regioni);
381 
382  Info<< setw(maxLen) << regionNames[regioni]
383  << setw(10) << rf.minLevel()[globali]
384  << setw(10) << rf.maxLevel()[globali]
385  << setw(10) << rf.gapLevel()[globali] << endl;
386  }
387  }
388 
389 
390  return surfacePtr;
391 }
392 
393 
394 void extractSurface
395 (
396  const polyMesh& mesh,
397  const Time& runTime,
398  const labelHashSet& includePatches,
399  const fileName& outFileName
400 )
401 {
403 
404  // Collect sizes. Hash on names to handle local-only patches (e.g.
405  // processor patches)
406  HashTable<label> patchSize(1024);
407  label nFaces = 0;
408  for (const label patchi : includePatches)
409  {
410  const polyPatch& pp = bMesh[patchi];
411  patchSize.insert(pp.name(), pp.size());
412  nFaces += pp.size();
413  }
415 
416 
417  // Allocate zone/patch for all patches
418  HashTable<label> compactZoneID(1024);
419  forAllConstIters(patchSize, iter)
420  {
421  label sz = compactZoneID.size();
422  compactZoneID.insert(iter.key(), sz);
423  }
424  Pstream::mapCombineScatter(compactZoneID);
425 
426 
427  // Rework HashTable into labelList just for speed of conversion
428  labelList patchToCompactZone(bMesh.size(), -1);
429  forAllConstIters(compactZoneID, iter)
430  {
431  label patchi = bMesh.findPatchID(iter.key());
432  if (patchi != -1)
433  {
434  patchToCompactZone[patchi] = iter();
435  }
436  }
437 
438  // Collect faces on zones
439  DynamicList<label> faceLabels(nFaces);
440  DynamicList<label> compactZones(nFaces);
441  for (const label patchi : includePatches)
442  {
443  const polyPatch& pp = bMesh[patchi];
444  forAll(pp, i)
445  {
446  faceLabels.append(pp.start()+i);
447  compactZones.append(patchToCompactZone[pp.index()]);
448  }
449  }
450 
451  // Addressing engine for all faces
452  uindirectPrimitivePatch allBoundary
453  (
454  UIndirectList<face>(mesh.faces(), faceLabels),
455  mesh.points()
456  );
457 
458 
459  // Find correspondence to master points
460  labelList pointToGlobal;
461  labelList uniqueMeshPoints;
463  (
464  allBoundary.meshPoints(),
465  allBoundary.meshPointMap(),
466  pointToGlobal,
467  uniqueMeshPoints
468  );
469 
470  // Gather all unique points on master
471  List<pointField> gatheredPoints(Pstream::nProcs());
472  gatheredPoints[Pstream::myProcNo()] = pointField
473  (
474  mesh.points(),
475  uniqueMeshPoints
476  );
477  Pstream::gatherList(gatheredPoints);
478 
479  // Gather all faces
480  List<faceList> gatheredFaces(Pstream::nProcs());
481  gatheredFaces[Pstream::myProcNo()] = allBoundary.localFaces();
482  forAll(gatheredFaces[Pstream::myProcNo()], i)
483  {
484  inplaceRenumber(pointToGlobal, gatheredFaces[Pstream::myProcNo()][i]);
485  }
486  Pstream::gatherList(gatheredFaces);
487 
488  // Gather all ZoneIDs
489  List<labelList> gatheredZones(Pstream::nProcs());
490  gatheredZones[Pstream::myProcNo()].transfer(compactZones);
491  Pstream::gatherList(gatheredZones);
492 
493  // On master combine all points, faces, zones
494  if (Pstream::master())
495  {
496  pointField allPoints = ListListOps::combine<pointField>
497  (
498  gatheredPoints,
500  );
501  gatheredPoints.clear();
502 
503  faceList allFaces = ListListOps::combine<faceList>
504  (
505  gatheredFaces,
507  );
508  gatheredFaces.clear();
509 
510  labelList allZones = ListListOps::combine<labelList>
511  (
512  gatheredZones,
514  );
515  gatheredZones.clear();
516 
517 
518  // Zones
519  surfZoneIdentifierList surfZones(compactZoneID.size());
520  forAllConstIters(compactZoneID, iter)
521  {
522  surfZones[iter()] = surfZoneIdentifier(iter.key(), iter());
523  Info<< "surfZone " << iter() << " : " << surfZones[iter()].name()
524  << endl;
525  }
526 
527  UnsortedMeshedSurface<face> unsortedFace
528  (
529  std::move(allPoints),
530  std::move(allFaces),
531  std::move(allZones),
532  surfZones
533  );
534 
535 
536  MeshedSurface<face> sortedFace(unsortedFace);
537 
538  fileName globalCasePath
539  (
541  ? runTime.globalPath()/outFileName
542  : runTime.path()/outFileName
543  );
544  globalCasePath.clean(); // Remove unneeded ".."
545 
546  Info<< "Writing merged surface to " << globalCasePath << endl;
547 
548  sortedFace.write(globalCasePath);
549  }
550 }
551 
552 
553 label checkAlignment(const polyMesh& mesh, const scalar tol, Ostream& os)
554 {
555  // Check all edges aligned with one of the coordinate axes
556  const faceList& faces = mesh.faces();
557  const pointField& points = mesh.points();
558 
559  label nUnaligned = 0;
560 
561  forAll(faces, facei)
562  {
563  const face& f = faces[facei];
564  forAll(f, fp)
565  {
566  label fp1 = f.fcIndex(fp);
567  const linePointRef e(edge(f[fp], f[fp1]).line(points));
568  const vector v(e.vec());
569  const scalar magV(mag(v));
570  if (magV > ROOTVSMALL)
571  {
572  for
573  (
574  direction dir = 0;
575  dir < pTraits<vector>::nComponents;
576  ++dir
577  )
578  {
579  const scalar s(mag(v[dir]));
580  if (s > magV*tol && s < magV*(1-tol))
581  {
582  ++nUnaligned;
583  break;
584  }
585  }
586  }
587  }
588  }
589 
590  reduce(nUnaligned, sumOp<label>());
591 
592  if (nUnaligned)
593  {
594  os << "Initial mesh has " << nUnaligned
595  << " edges unaligned with any of the coordinate axes" << nl << endl;
596  }
597  return nUnaligned;
598 }
599 
600 
601 // Check writing tolerance before doing any serious work
602 scalar getMergeDistance
603 (
604  const polyMesh& mesh,
605  const scalar mergeTol,
606  const bool dryRun
607 )
608 {
609  const boundBox& meshBb = mesh.bounds();
610  scalar mergeDist = mergeTol * meshBb.mag();
611 
612  Info<< nl
613  << "Overall mesh bounding box : " << meshBb << nl
614  << "Relative tolerance : " << mergeTol << nl
615  << "Absolute matching distance : " << mergeDist << nl
616  << endl;
617 
618  // check writing tolerance
619  if (mesh.time().writeFormat() == IOstream::ASCII && !dryRun)
620  {
621  const scalar writeTol = std::pow
622  (
623  scalar(10),
624  -scalar(IOstream::defaultPrecision())
625  );
626 
627  if (mergeTol < writeTol)
628  {
630  << "Your current settings specify ASCII writing with "
631  << IOstream::defaultPrecision() << " digits precision." << nl
632  << "Your merging tolerance (" << mergeTol
633  << ") is finer than this." << nl
634  << "Change to binary writeFormat, "
635  << "or increase the writePrecision" << endl
636  << "or adjust the merge tolerance (mergeTol)."
637  << exit(FatalError);
638  }
639  }
640 
641  return mergeDist;
642 }
643 
644 
645 void removeZeroSizedPatches(fvMesh& mesh)
646 {
647  // Remove any zero-sized ones. Assumes
648  // - processor patches are already only there if needed
649  // - all other patches are available on all processors
650  // - but coupled ones might still be needed, even if zero-size
651  // (e.g. processorCyclic)
652  // See also logic in createPatch.
653  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
654 
655  labelList oldToNew(pbm.size(), -1);
656  label newPatchi = 0;
657  forAll(pbm, patchi)
658  {
659  const polyPatch& pp = pbm[patchi];
660 
661  if (!isA<processorPolyPatch>(pp))
662  {
663  if
664  (
665  isA<coupledPolyPatch>(pp)
666  || returnReduce(pp.size(), sumOp<label>())
667  )
668  {
669  // Coupled (and unknown size) or uncoupled and used
670  oldToNew[patchi] = newPatchi++;
671  }
672  }
673  }
674 
675  forAll(pbm, patchi)
676  {
677  const polyPatch& pp = pbm[patchi];
678 
679  if (isA<processorPolyPatch>(pp))
680  {
681  oldToNew[patchi] = newPatchi++;
682  }
683  }
684 
685 
686  const label nKeepPatches = newPatchi;
687 
688  // Shuffle unused ones to end
689  if (nKeepPatches != pbm.size())
690  {
691  Info<< endl
692  << "Removing zero-sized patches:" << endl << incrIndent;
693 
694  forAll(oldToNew, patchi)
695  {
696  if (oldToNew[patchi] == -1)
697  {
698  Info<< indent << pbm[patchi].name()
699  << " type " << pbm[patchi].type()
700  << " at position " << patchi << endl;
701  oldToNew[patchi] = newPatchi++;
702  }
703  }
704  Info<< decrIndent;
705 
706  fvMeshTools::reorderPatches(mesh, oldToNew, nKeepPatches, true);
707  Info<< endl;
708  }
709 }
710 
711 
712 // Write mesh and additional information
713 void writeMesh
714 (
715  const string& msg,
716  const meshRefinement& meshRefiner,
717  const meshRefinement::debugType debugLevel,
718  const meshRefinement::writeType writeLevel
719 )
720 {
721  const fvMesh& mesh = meshRefiner.mesh();
722 
723  meshRefiner.printMeshInfo(debugLevel, msg);
724  Info<< "Writing mesh to time " << meshRefiner.timeName() << endl;
725 
727  if (!debugLevel && !(writeLevel&meshRefinement::WRITELAYERSETS))
728  {
730  }
732 
733  meshRefiner.write
734  (
735  debugLevel,
737  mesh.time().path()/meshRefiner.timeName()
738  );
739  Info<< "Wrote mesh in = "
740  << mesh.time().cpuTimeIncrement() << " s." << endl;
741 }
742 
743 
744 int main(int argc, char *argv[])
745 {
747  (
748  "Automatic split hex mesher. Refines and snaps to surface"
749  );
750 
751  #include "addRegionOption.H"
752  #include "addOverwriteOption.H"
753  #include "addProfilingOption.H"
755  (
756  "checkGeometry",
757  "Check all surface geometry for quality"
758  );
760  (
761  "Check case set-up only using a single time step"
762  );
764  (
765  "surfaceSimplify",
766  "boundBox",
767  "Simplify the surface using snappyHexMesh starting from a boundBox"
768  );
770  (
771  "patches",
772  "(patch0 .. patchN)",
773  "Only triangulate selected patches (wildcards supported)"
774  );
776  (
777  "outFile",
778  "file",
779  "Name of the file to save the simplified surface to"
780  );
781  argList::addOption("dict", "file", "Alternative snappyHexMeshDict");
782 
783  argList::noFunctionObjects(); // Never use function objects
784 
785  #include "setRootCase.H"
786  #include "createTime.H"
787 
788  const bool overwrite = args.found("overwrite");
789  const bool checkGeometry = args.found("checkGeometry");
790  const bool surfaceSimplify = args.found("surfaceSimplify");
791  const bool dryRun = args.dryRun();
792 
793  if (dryRun)
794  {
795  Info<< "Operating in dry-run mode to detect set-up errors"
796  << nl << endl;
797  }
798 
799 
800  #include "createNamedMesh.H"
801  Info<< "Read mesh in = "
802  << runTime.cpuTimeIncrement() << " s" << endl;
803 
804  // Check patches and faceZones are synchronised
807 
808  if (dryRun)
809  {
810  // Check if mesh aligned with cartesian axes
811  checkAlignment(mesh, 1e-6, Pout); //FatalIOError);
812  }
813 
814 
815 
816  // Read meshing dictionary
817  const word dictName("snappyHexMeshDict");
818  #include "setSystemMeshDictionaryIO.H"
820 
821 
822  // all surface geometry
823  const dictionary& geometryDict =
824  meshRefinement::subDict(meshDict, "geometry", dryRun);
825 
826  // refinement parameters
827  const dictionary& refineDict =
828  meshRefinement::subDict(meshDict, "castellatedMeshControls", dryRun);
829 
830  // mesh motion and mesh quality parameters
831  const dictionary& motionDict =
832  meshRefinement::subDict(meshDict, "meshQualityControls", dryRun);
833 
834  // snap-to-surface parameters
835  const dictionary& snapDict =
836  meshRefinement::subDict(meshDict, "snapControls", dryRun);
837 
838  // layer addition parameters
839  const dictionary& layerDict =
840  meshRefinement::subDict(meshDict, "addLayersControls", dryRun);
841 
842  // absolute merge distance
843  const scalar mergeDist = getMergeDistance
844  (
845  mesh,
846  meshRefinement::get<scalar>
847  (
848  meshDict,
849  "mergeTolerance",
850  dryRun
851  ),
852  dryRun
853  );
854 
855  const bool keepPatches(meshDict.getOrDefault("keepPatches", false));
856 
857  // format to be used for writing lines
858  const word setFormat
859  (
861  (
862  "setFormat",
864  )
865  );
866  const autoPtr<writer<scalar>> setFormatter
867  (
869  );
870 
871  const scalar maxSizeRatio
872  (
873  meshDict.getOrDefault<scalar>("maxSizeRatio", 100)
874  );
875 
876 
877  // Read decomposePar dictionary
878  dictionary decomposeDict;
879  if (Pstream::parRun())
880  {
881  // Ensure demand-driven decompositionMethod finds alternative
882  // decomposeParDict location properly.
883 
884  IOdictionary* dictPtr = new IOdictionary
885  (
887  (
888  IOobject
889  (
891  runTime.system(),
892  runTime,
895  ),
896  args.getOrDefault<fileName>("decomposeParDict", "")
897  )
898  );
899 
900  // Store it on the object registry, but to be found it must also
901  // have the expected "decomposeParDict" name.
902 
904  runTime.store(dictPtr);
905 
906  decomposeDict = *dictPtr;
907  }
908  else
909  {
910  decomposeDict.add("method", "none");
911  decomposeDict.add("numberOfSubdomains", 1);
912  }
913 
914 
915  // Debug
916  // ~~~~~
917 
918  // Set debug level
920  (
921  meshDict.getOrDefault<label>
922  (
923  "debug",
924  0
925  )
926  );
927  {
928  wordList flags;
929  if (meshDict.readIfPresent("debugFlags", flags))
930  {
931  debugLevel = meshRefinement::debugType
932  (
934  (
936  flags
937  )
938  );
939  }
940  }
941  if (debugLevel > 0)
942  {
943  meshRefinement::debug = debugLevel;
944  snappyRefineDriver::debug = debugLevel;
945  snappySnapDriver::debug = debugLevel;
946  snappyLayerDriver::debug = debugLevel;
947  }
948 
949  // Set file writing level
950  {
951  wordList flags;
952  if (meshDict.readIfPresent("writeFlags", flags))
953  {
955  (
957  (
959  (
961  flags
962  )
963  )
964  );
965  }
966  }
967 
969  //{
970  // wordList flags;
971  // if (meshDict.readIfPresent("outputFlags", flags))
972  // {
973  // meshRefinement::outputLevel
974  // (
975  // meshRefinement::outputType
976  // (
977  // meshRefinement::readFlags
978  // (
979  // meshRefinement::outputTypeNames,
980  // flags
981  // )
982  // )
983  // );
984  // }
985  //}
986 
987  // for the impatient who want to see some output files:
989 
990  // Read geometry
991  // ~~~~~~~~~~~~~
992 
993  searchableSurfaces allGeometry
994  (
995  IOobject
996  (
997  "abc", // dummy name
998  mesh.time().constant(), // instance
999  //mesh.time().findInstance("triSurface", word::null),// instance
1000  "triSurface", // local
1001  mesh.time(), // registry
1004  ),
1005  geometryDict,
1006  meshDict.getOrDefault("singleRegionName", true)
1007  );
1008 
1009 
1010  // Read refinement surfaces
1011  // ~~~~~~~~~~~~~~~~~~~~~~~~
1012 
1013  autoPtr<refinementSurfaces> surfacesPtr;
1014 
1015  Info<< "Reading refinement surfaces." << endl;
1016 
1017  if (surfaceSimplify)
1018  {
1019  addProfiling(surfaceSimplify, "snappyHexMesh::surfaceSimplify");
1020  IOdictionary foamyHexMeshDict
1021  (
1022  IOobject
1023  (
1024  "foamyHexMeshDict",
1025  runTime.system(),
1026  runTime,
1029  )
1030  );
1031 
1032  const dictionary& conformationDict =
1033  foamyHexMeshDict.subDict("surfaceConformation").subDict
1034  (
1035  "geometryToConformTo"
1036  );
1037 
1038  const dictionary& motionDict =
1039  foamyHexMeshDict.subDict("motionControl");
1040 
1041  const dictionary& shapeControlDict =
1042  motionDict.subDict("shapeControlFunctions");
1043 
1044  // Calculate current ratio of hex cells v.s. wanted cell size
1045  const scalar defaultCellSize =
1046  motionDict.get<scalar>("defaultCellSize");
1047 
1048  const scalar initialCellSize = ::pow(mesh.V()[0], 1.0/3.0);
1049 
1050  //Info<< "Wanted cell size = " << defaultCellSize << endl;
1051  //Info<< "Current cell size = " << initialCellSize << endl;
1052  //Info<< "Fraction = " << initialCellSize/defaultCellSize
1053  // << endl;
1054 
1055  surfacesPtr =
1056  createRefinementSurfaces
1057  (
1058  allGeometry,
1059  conformationDict,
1060  shapeControlDict,
1061  refineDict.getOrDefault("gapLevelIncrement", 0),
1062  initialCellSize/defaultCellSize
1063  );
1064 
1066  }
1067  else
1068  {
1069  surfacesPtr.reset
1070  (
1071  new refinementSurfaces
1072  (
1073  allGeometry,
1075  (
1076  refineDict,
1077  "refinementSurfaces",
1078  dryRun
1079  ),
1080  refineDict.getOrDefault("gapLevelIncrement", 0),
1081  dryRun
1082  )
1083  );
1084 
1085  Info<< "Read refinement surfaces in = "
1086  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
1087  }
1088 
1089  refinementSurfaces& surfaces = surfacesPtr();
1090 
1091 
1092  // Checking only?
1093 
1094  if (checkGeometry)
1095  {
1096  // Check geometry amongst itself (e.g. intersection, size differences)
1097 
1098  // Extract patchInfo
1099  List<wordList> patchTypes(allGeometry.size());
1100 
1101  const PtrList<dictionary>& patchInfo = surfaces.patchInfo();
1102  const labelList& surfaceGeometry = surfaces.surfaces();
1103  forAll(surfaceGeometry, surfi)
1104  {
1105  label geomi = surfaceGeometry[surfi];
1106  const wordList& regNames = allGeometry.regionNames()[geomi];
1107 
1108  patchTypes[geomi].setSize(regNames.size());
1109  forAll(regNames, regioni)
1110  {
1111  label globalRegioni = surfaces.globalRegion(surfi, regioni);
1112 
1113  if (patchInfo.set(globalRegioni))
1114  {
1115  patchTypes[geomi][regioni] =
1116  meshRefinement::get<word>
1117  (
1118  patchInfo[globalRegioni],
1119  "type",
1120  dryRun,
1122  word::null
1123  );
1124  }
1125  else
1126  {
1127  patchTypes[geomi][regioni] = wallPolyPatch::typeName;
1128  }
1129  }
1130  }
1131 
1132  // Write some stats
1133  allGeometry.writeStats(patchTypes, Info);
1134  // Check topology
1135  allGeometry.checkTopology(true);
1136  // Check geometry
1137  allGeometry.checkGeometry
1138  (
1139  maxSizeRatio, // max size ratio
1140  1e-9, // intersection tolerance
1141  setFormatter,
1142  0.01, // min triangle quality
1143  true
1144  );
1145 
1146  if (!dryRun)
1147  {
1148  return 0;
1149  }
1150  }
1151 
1152 
1153  if (dryRun)
1154  {
1155  // Check geometry to mesh bounding box
1156  Info<< "Checking for geometry size relative to mesh." << endl;
1157  const boundBox& meshBb = mesh.bounds();
1158  forAll(allGeometry, geomi)
1159  {
1160  const searchableSurface& s = allGeometry[geomi];
1161  const boundBox& bb = s.bounds();
1162 
1163  scalar ratio = bb.mag() / meshBb.mag();
1164  if (ratio > maxSizeRatio || ratio < 1.0/maxSizeRatio)
1165  {
1166  Warning
1167  << " " << allGeometry.names()[geomi]
1168  << " bounds differ from mesh"
1169  << " by more than a factor " << maxSizeRatio << ":" << nl
1170  << " bounding box : " << bb << nl
1171  << " mesh bounding box : " << meshBb
1172  << endl;
1173  }
1174  if (!meshBb.contains(bb))
1175  {
1176  Warning
1177  << " " << allGeometry.names()[geomi]
1178  << " bounds not fully contained in mesh" << nl
1179  << " bounding box : " << bb << nl
1180  << " mesh bounding box : " << meshBb
1181  << endl;
1182  }
1183  }
1184  Info<< endl;
1185  }
1186 
1187 
1188 
1189 
1190  // Read refinement shells
1191  // ~~~~~~~~~~~~~~~~~~~~~~
1192 
1193  Info<< "Reading refinement shells." << endl;
1194  shellSurfaces shells
1195  (
1196  allGeometry,
1197  meshRefinement::subDict(refineDict, "refinementRegions", dryRun),
1198  dryRun
1199  );
1200  Info<< "Read refinement shells in = "
1201  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
1202 
1203 
1204  Info<< "Setting refinement level of surface to be consistent"
1205  << " with shells." << endl;
1206  surfaces.setMinLevelFields(shells);
1207  Info<< "Checked shell refinement in = "
1208  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
1209 
1210 
1211  // Optionally read limit shells
1212  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1213 
1214  const dictionary limitDict(refineDict.subOrEmptyDict("limitRegions"));
1215 
1216  if (!limitDict.empty())
1217  {
1218  Info<< "Reading limit shells." << endl;
1219  }
1220 
1221  shellSurfaces limitShells(allGeometry, limitDict, dryRun);
1222 
1223  if (!limitDict.empty())
1224  {
1225  Info<< "Read limit shells in = "
1226  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
1227  }
1228 
1229  if (dryRun)
1230  {
1231  // Check for use of all geometry
1232  const wordList& allGeomNames = allGeometry.names();
1233 
1234  labelHashSet unusedGeometries(identity(allGeomNames.size()));
1235  unusedGeometries.erase(surfaces.surfaces());
1236  unusedGeometries.erase(shells.shells());
1237  unusedGeometries.erase(limitShells.shells());
1238 
1239  if (unusedGeometries.size())
1240  {
1241  IOWarningInFunction(geometryDict)
1242  << "The following geometry entries are not used:" << nl;
1243  for (const label geomi : unusedGeometries)
1244  {
1245  Info<< " " << allGeomNames[geomi] << nl;
1246  }
1247  Info<< endl;
1248  }
1249  }
1250 
1251 
1252 
1253 
1254  // Read feature meshes
1255  // ~~~~~~~~~~~~~~~~~~~
1256 
1257  Info<< "Reading features." << endl;
1258  refinementFeatures features
1259  (
1260  mesh,
1262  (
1263  meshRefinement::lookup(refineDict, "features", dryRun)
1264  ),
1265  dryRun
1266  );
1267  Info<< "Read features in = "
1268  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
1269 
1270 
1271  if (dryRun)
1272  {
1273  // Check geometry to mesh bounding box
1274  Info<< "Checking for line geometry size relative to surface geometry."
1275  << endl;
1276 
1277  OStringStream os;
1278  bool hasErrors = features.checkSizes
1279  (
1280  maxSizeRatio, //const scalar maxRatio,
1281  mesh.bounds(),
1282  true, //const bool report,
1283  os //FatalIOError
1284  );
1285  if (hasErrors)
1286  {
1287  Warning<< os.str() << endl;
1288  }
1289  }
1290 
1291 
1292  // Refinement engine
1293  // ~~~~~~~~~~~~~~~~~
1294 
1295  Info<< nl
1296  << "Determining initial surface intersections" << nl
1297  << "-----------------------------------------" << nl
1298  << endl;
1299 
1300  // Main refinement engine
1301  meshRefinement meshRefiner
1302  (
1303  mesh,
1304  mergeDist, // tolerance used in sorting coordinates
1305  overwrite, // overwrite mesh files?
1306  surfaces, // for surface intersection refinement
1307  features, // for feature edges/point based refinement
1308  shells, // for volume (inside/outside) refinement
1309  limitShells, // limit of volume refinement
1310  labelList(), // initial faces to test
1311  dryRun
1312  );
1313 
1314  if (!dryRun)
1315  {
1316  meshRefiner.updateIntersections(identity(mesh.nFaces()));
1317  Info<< "Calculated surface intersections in = "
1318  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
1319  }
1320 
1321  // Some stats
1322  meshRefiner.printMeshInfo(debugLevel, "Initial mesh");
1323 
1324  meshRefiner.write
1325  (
1328  mesh.time().path()/meshRefiner.timeName()
1329  );
1330 
1331 
1332  // Refinement parameters
1333  const refinementParameters refineParams(refineDict, dryRun);
1334 
1335  // Snap parameters
1336  const snapParameters snapParams(snapDict, dryRun);
1337 
1338 
1339 
1340  // Add all the cellZones and faceZones
1341  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1342 
1343  // 1. cellZones relating to surface (faceZones added later)
1344 
1345  const labelList namedSurfaces
1346  (
1348  );
1349 
1351  (
1352  surfaces.surfZones(),
1353  namedSurfaces,
1354  mesh
1355  );
1356 
1357 
1358  // 2. cellZones relating to locations
1359 
1360  refineParams.addCellZonesToMesh(mesh);
1361 
1362 
1363 
1364  // Add all the surface regions as patches
1365  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1366 
1367  //- Global surface region to patch (non faceZone surface) or patches
1368  // (faceZone surfaces)
1369  labelList globalToMasterPatch;
1370  labelList globalToSlavePatch;
1371 
1372 
1373  {
1374  Info<< nl
1375  << "Adding patches for surface regions" << nl
1376  << "----------------------------------" << nl
1377  << endl;
1378 
1379  // From global region number to mesh patch.
1380  globalToMasterPatch.setSize(surfaces.nRegions(), -1);
1381  globalToSlavePatch.setSize(surfaces.nRegions(), -1);
1382 
1383  if (!dryRun)
1384  {
1385  Info<< setf(ios_base::left)
1386  << setw(6) << "Patch"
1387  << setw(20) << "Type"
1388  << setw(30) << "Region" << nl
1389  << setw(6) << "-----"
1390  << setw(20) << "----"
1391  << setw(30) << "------" << endl;
1392  }
1393 
1394  const labelList& surfaceGeometry = surfaces.surfaces();
1395  const PtrList<dictionary>& surfacePatchInfo = surfaces.patchInfo();
1396  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
1397 
1398  forAll(surfaceGeometry, surfi)
1399  {
1400  label geomi = surfaceGeometry[surfi];
1401 
1402  const wordList& regNames = allGeometry.regionNames()[geomi];
1403 
1404  if (!dryRun)
1405  {
1406  Info<< surfaces.names()[surfi] << ':' << nl << nl;
1407  }
1408 
1409  const wordList& fzNames =
1410  surfaces.surfZones()[surfi].faceZoneNames();
1411 
1412  if (fzNames.empty())
1413  {
1414  // 'Normal' surface
1415  forAll(regNames, i)
1416  {
1417  label globalRegioni = surfaces.globalRegion(surfi, i);
1418 
1419  label patchi;
1420 
1421  if (surfacePatchInfo.set(globalRegioni))
1422  {
1423  patchi = meshRefiner.addMeshedPatch
1424  (
1425  regNames[i],
1426  surfacePatchInfo[globalRegioni]
1427  );
1428  }
1429  else
1430  {
1431  dictionary patchInfo;
1432  patchInfo.set("type", wallPolyPatch::typeName);
1433 
1434  patchi = meshRefiner.addMeshedPatch
1435  (
1436  regNames[i],
1437  patchInfo
1438  );
1439  }
1440 
1441  if (!dryRun)
1442  {
1443  Info<< setf(ios_base::left)
1444  << setw(6) << patchi
1445  << setw(20) << pbm[patchi].type()
1446  << setw(30) << regNames[i] << nl;
1447  }
1448 
1449  globalToMasterPatch[globalRegioni] = patchi;
1450  globalToSlavePatch[globalRegioni] = patchi;
1451  }
1452  }
1453  else
1454  {
1455  // Zoned surface
1456  forAll(regNames, i)
1457  {
1458  label globalRegioni = surfaces.globalRegion(surfi, i);
1459 
1460  // Add master side patch
1461  {
1462  label patchi;
1463 
1464  if (surfacePatchInfo.set(globalRegioni))
1465  {
1466  patchi = meshRefiner.addMeshedPatch
1467  (
1468  regNames[i],
1469  surfacePatchInfo[globalRegioni]
1470  );
1471  }
1472  else
1473  {
1474  dictionary patchInfo;
1475  patchInfo.set("type", wallPolyPatch::typeName);
1476 
1477  patchi = meshRefiner.addMeshedPatch
1478  (
1479  regNames[i],
1480  patchInfo
1481  );
1482  }
1483 
1484  if (!dryRun)
1485  {
1486  Info<< setf(ios_base::left)
1487  << setw(6) << patchi
1488  << setw(20) << pbm[patchi].type()
1489  << setw(30) << regNames[i] << nl;
1490  }
1491 
1492  globalToMasterPatch[globalRegioni] = patchi;
1493  }
1494  // Add slave side patch
1495  {
1496  const word slaveName = regNames[i] + "_slave";
1497  label patchi;
1498 
1499  if (surfacePatchInfo.set(globalRegioni))
1500  {
1501  patchi = meshRefiner.addMeshedPatch
1502  (
1503  slaveName,
1504  surfacePatchInfo[globalRegioni]
1505  );
1506  }
1507  else
1508  {
1509  dictionary patchInfo;
1510  patchInfo.set("type", wallPolyPatch::typeName);
1511 
1512  patchi = meshRefiner.addMeshedPatch
1513  (
1514  slaveName,
1515  patchInfo
1516  );
1517  }
1518 
1519  if (!dryRun)
1520  {
1521  Info<< setf(ios_base::left)
1522  << setw(6) << patchi
1523  << setw(20) << pbm[patchi].type()
1524  << setw(30) << slaveName << nl;
1525  }
1526 
1527  globalToSlavePatch[globalRegioni] = patchi;
1528  }
1529  }
1530 
1531  // For now: have single faceZone per surface. Use first
1532  // region in surface for patch for zoning
1533  if (regNames.size())
1534  {
1535  forAll(fzNames, fzi)
1536  {
1537  const word& fzName = fzNames[fzi];
1538  label globalRegioni = surfaces.globalRegion(surfi, fzi);
1539 
1540  meshRefiner.addFaceZone
1541  (
1542  fzName,
1543  pbm[globalToMasterPatch[globalRegioni]].name(),
1544  pbm[globalToSlavePatch[globalRegioni]].name(),
1545  surfaces.surfZones()[surfi].faceType()
1546  );
1547  }
1548  }
1549  }
1550 
1551  if (!dryRun)
1552  {
1553  Info<< nl;
1554  }
1555  }
1556  Info<< "Added patches in = "
1557  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
1558  }
1559 
1560 
1561 
1562  // Add all information for all the remaining faceZones
1563  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1564 
1565  HashTable<Pair<word>> faceZoneToPatches;
1566  forAll(mesh.faceZones(), zonei)
1567  {
1568  const word& fzName = mesh.faceZones()[zonei].name();
1569 
1570  label mpI, spI;
1572  bool hasInfo = meshRefiner.getFaceZoneInfo(fzName, mpI, spI, fzType);
1573 
1574  if (!hasInfo)
1575  {
1576  // faceZone does not originate from a surface but presumably
1577  // from a cellZone pair instead
1578  string::size_type i = fzName.find("_to_");
1579  if (i != string::npos)
1580  {
1581  word cz0 = fzName.substr(0, i);
1582  word cz1 = fzName.substr(i+4, fzName.size()-i+4);
1583  word slaveName(cz1 + "_to_" + cz0);
1584  faceZoneToPatches.insert(fzName, Pair<word>(fzName, slaveName));
1585  }
1586  else
1587  {
1588  // Add as fzName + fzName_slave
1589  const word slaveName = fzName + "_slave";
1590  faceZoneToPatches.insert(fzName, Pair<word>(fzName, slaveName));
1591  }
1592  }
1593  }
1594 
1595  if (faceZoneToPatches.size())
1596  {
1598  (
1599  meshRefiner,
1600  refineParams,
1601  faceZoneToPatches
1602  );
1603  }
1604 
1605 
1606 
1607  // Re-do intersections on meshed boundaries since they use an extrapolated
1608  // other side
1609  {
1610  const labelList adaptPatchIDs(meshRefiner.meshedPatches());
1611 
1612  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
1613 
1614  label nFaces = 0;
1615  forAll(adaptPatchIDs, i)
1616  {
1617  nFaces += pbm[adaptPatchIDs[i]].size();
1618  }
1619 
1620  labelList faceLabels(nFaces);
1621  nFaces = 0;
1622  forAll(adaptPatchIDs, i)
1623  {
1624  const polyPatch& pp = pbm[adaptPatchIDs[i]];
1625  forAll(pp, i)
1626  {
1627  faceLabels[nFaces++] = pp.start()+i;
1628  }
1629  }
1630  meshRefiner.updateIntersections(faceLabels);
1631  }
1632 
1633 
1634 
1635  // Parallel
1636  // ~~~~~~~~
1637 
1638  // Construct decomposition engine. Note: cannot use decompositionModel
1639  // MeshObject since we're clearing out the mesh inside the mesh generation.
1640  autoPtr<decompositionMethod> decomposerPtr
1641  (
1643  (
1644  decomposeDict
1645  )
1646  );
1647  decompositionMethod& decomposer = *decomposerPtr;
1648 
1649  if (Pstream::parRun() && !decomposer.parallelAware())
1650  {
1652  << "You have selected decomposition method "
1653  << decomposer.typeName
1654  << " which is not parallel aware." << endl
1655  << "Please select one that is (hierarchical, ptscotch)"
1656  << exit(FatalError);
1657  }
1658 
1659  // Mesh distribution engine (uses tolerance to reconstruct meshes)
1660  fvMeshDistribute distributor(mesh);
1661 
1662 
1663 
1664 
1665 
1666  // Now do the real work -refinement -snapping -layers
1667  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1668 
1669  const bool wantRefine
1670  (
1671  meshRefinement::get<bool>(meshDict, "castellatedMesh", dryRun)
1672  );
1673  const bool wantSnap
1674  (
1675  meshRefinement::get<bool>(meshDict, "snap", dryRun)
1676  );
1677  const bool wantLayers
1678  (
1679  meshRefinement::get<bool>(meshDict, "addLayers", dryRun)
1680  );
1681 
1682  if (dryRun)
1683  {
1684  string errorMsg(FatalError.message());
1685  string IOerrorMsg(FatalIOError.message());
1686 
1687  if (errorMsg.size() || IOerrorMsg.size())
1688  {
1689  //errorMsg = "[dryRun] " + errorMsg;
1690  //errorMsg.replaceAll("\n", "\n[dryRun] ");
1691  //IOerrorMsg = "[dryRun] " + IOerrorMsg;
1692  //IOerrorMsg.replaceAll("\n", "\n[dryRun] ");
1693 
1694  Warning
1695  << nl
1696  << "Missing/incorrect required dictionary entries:" << nl
1697  << nl
1698  << IOerrorMsg.c_str() << nl
1699  << errorMsg.c_str() << nl << nl
1700  << "Exiting dry-run" << nl << endl;
1701 
1702  FatalError.clear();
1703  FatalIOError.clear();
1704 
1705  return 0;
1706  }
1707  }
1708 
1709 
1710  // How to treat co-planar faces
1711  meshRefinement::FaceMergeType mergeType =
1712  meshRefinement::FaceMergeType::GEOMETRIC;
1713  {
1714  const bool mergePatchFaces
1715  (
1716  meshDict.getOrDefault("mergePatchFaces", true)
1717  );
1718 
1719  if (!mergePatchFaces)
1720  {
1721  Info<< "Not merging patch-faces of cell to preserve"
1722  << " (split)hex cell shape."
1723  << nl << endl;
1724  mergeType = meshRefinement::FaceMergeType::NONE;
1725  }
1726  else
1727  {
1728  const bool mergeAcrossPatches
1729  (
1730  meshDict.getOrDefault("mergeAcrossPatches", false)
1731  );
1732 
1733  if (mergeAcrossPatches)
1734  {
1735  Info<< "Merging co-planar patch-faces of cells"
1736  << ", regardless of patch assignment"
1737  << nl << endl;
1738  mergeType = meshRefinement::FaceMergeType::IGNOREPATCH;
1739  }
1740  }
1741  }
1742 
1743 
1744 
1745  if (wantRefine)
1746  {
1747  cpuTime timer;
1748 
1749  snappyRefineDriver refineDriver
1750  (
1751  meshRefiner,
1752  decomposer,
1753  distributor,
1754  globalToMasterPatch,
1755  globalToSlavePatch,
1756  setFormatter(),
1757  dryRun
1758  );
1759 
1760 
1761  if (!overwrite && !debugLevel)
1762  {
1763  const_cast<Time&>(mesh.time())++;
1764  }
1765 
1766 
1767  refineDriver.doRefine
1768  (
1769  refineDict,
1770  refineParams,
1771  snapParams,
1772  refineParams.handleSnapProblems(),
1773  mergeType,
1774  motionDict
1775  );
1776 
1777  // Remove zero sized patches originating from faceZones
1778  if (!keepPatches && !wantSnap && !wantLayers)
1779  {
1781  }
1782 
1783  if (!dryRun)
1784  {
1785  writeMesh
1786  (
1787  "Refined mesh",
1788  meshRefiner,
1789  debugLevel,
1791  );
1792  }
1793 
1794  Info<< "Mesh refined in = "
1795  << timer.cpuTimeIncrement() << " s." << endl;
1796 
1798  }
1799 
1800  if (wantSnap)
1801  {
1802  cpuTime timer;
1803 
1804  snappySnapDriver snapDriver
1805  (
1806  meshRefiner,
1807  globalToMasterPatch,
1808  globalToSlavePatch,
1809  dryRun
1810  );
1811 
1812  if (!overwrite && !debugLevel)
1813  {
1814  const_cast<Time&>(mesh.time())++;
1815  }
1816 
1817  // Use the resolveFeatureAngle from the refinement parameters
1818  scalar curvature = refineParams.curvature();
1819  scalar planarAngle = refineParams.planarAngle();
1820 
1821  snapDriver.doSnap
1822  (
1823  snapDict,
1824  motionDict,
1825  mergeType,
1826  curvature,
1827  planarAngle,
1828  snapParams
1829  );
1830 
1831  // Remove zero sized patches originating from faceZones
1832  if (!keepPatches && !wantLayers)
1833  {
1835  }
1836 
1837  if (!dryRun)
1838  {
1839  writeMesh
1840  (
1841  "Snapped mesh",
1842  meshRefiner,
1843  debugLevel,
1845  );
1846  }
1847 
1848  Info<< "Mesh snapped in = "
1849  << timer.cpuTimeIncrement() << " s." << endl;
1850 
1852  }
1853 
1854  if (wantLayers)
1855  {
1856  cpuTime timer;
1857 
1858  // Layer addition parameters
1859  const layerParameters layerParams
1860  (
1861  layerDict,
1862  mesh.boundaryMesh(),
1863  dryRun
1864  );
1865 
1866  snappyLayerDriver layerDriver
1867  (
1868  meshRefiner,
1869  globalToMasterPatch,
1870  globalToSlavePatch,
1871  dryRun
1872  );
1873 
1874  // Use the maxLocalCells from the refinement parameters
1875  bool preBalance = returnReduce
1876  (
1877  (mesh.nCells() >= refineParams.maxLocalCells()),
1878  orOp<bool>()
1879  );
1880 
1881 
1882  if (!overwrite && !debugLevel)
1883  {
1884  const_cast<Time&>(mesh.time())++;
1885  }
1886 
1887  layerDriver.doLayers
1888  (
1889  layerDict,
1890  motionDict,
1891  layerParams,
1892  mergeType,
1893  preBalance,
1894  decomposer,
1895  distributor
1896  );
1897 
1898  // Remove zero sized patches originating from faceZones
1899  if (!keepPatches)
1900  {
1902  }
1903 
1904  if (!dryRun)
1905  {
1906  writeMesh
1907  (
1908  "Layer mesh",
1909  meshRefiner,
1910  debugLevel,
1912  );
1913  }
1914 
1915  Info<< "Layers added in = "
1916  << timer.cpuTimeIncrement() << " s." << endl;
1917 
1919  }
1920 
1921 
1922  {
1923  addProfiling(checkMesh, "snappyHexMesh::checkMesh");
1924 
1925  // Check final mesh
1926  Info<< "Checking final mesh ..." << endl;
1927  faceSet wrongFaces(mesh, "wrongFaces", mesh.nFaces()/100);
1928  motionSmoother::checkMesh(false, mesh, motionDict, wrongFaces, dryRun);
1929  const label nErrors = returnReduce
1930  (
1931  wrongFaces.size(),
1932  sumOp<label>()
1933  );
1934 
1935  if (nErrors > 0)
1936  {
1937  Info<< "Finished meshing with " << nErrors << " illegal faces"
1938  << " (concave, zero area or negative cell pyramid volume)"
1939  << endl;
1940  wrongFaces.write();
1941  }
1942  else
1943  {
1944  Info<< "Finished meshing without any errors" << endl;
1945  }
1946 
1948  }
1949 
1950 
1951  if (surfaceSimplify)
1952  {
1953  addProfiling(surfaceSimplify, "snappyHexMesh::surfaceSimplify");
1954 
1956 
1957  labelHashSet includePatches(bMesh.size());
1958 
1959  if (args.found("patches"))
1960  {
1961  includePatches = bMesh.patchSet
1962  (
1963  args.getList<wordRe>("patches")
1964  );
1965  }
1966  else
1967  {
1968  forAll(bMesh, patchi)
1969  {
1970  const polyPatch& patch = bMesh[patchi];
1971 
1972  if (!isA<processorPolyPatch>(patch))
1973  {
1974  includePatches.insert(patchi);
1975  }
1976  }
1977  }
1978 
1979  fileName outFileName
1980  (
1982  (
1983  "outFile",
1984  "constant/triSurface/simplifiedSurface.stl"
1985  )
1986  );
1987 
1988  extractSurface
1989  (
1990  mesh,
1991  runTime,
1992  includePatches,
1993  outFileName
1994  );
1995 
1996  pointIOField cellCentres
1997  (
1998  IOobject
1999  (
2000  "internalCellCentres",
2001  runTime.timeName(),
2002  mesh,
2005  ),
2006  mesh.cellCentres()
2007  );
2008 
2009  cellCentres.write();
2010  }
2011 
2013 
2014  Info<< "Finished meshing in = "
2015  << runTime.elapsedCpuTime() << " s." << endl;
2016 
2017 
2018  if (dryRun)
2019  {
2020  string errorMsg(FatalError.message());
2021  string IOerrorMsg(FatalIOError.message());
2022 
2023  if (errorMsg.size() || IOerrorMsg.size())
2024  {
2025  //errorMsg = "[dryRun] " + errorMsg;
2026  //errorMsg.replaceAll("\n", "\n[dryRun] ");
2027  //IOerrorMsg = "[dryRun] " + IOerrorMsg;
2028  //IOerrorMsg.replaceAll("\n", "\n[dryRun] ");
2029 
2030  Perr<< nl
2031  << "Missing/incorrect required dictionary entries:" << nl
2032  << nl
2033  << IOerrorMsg.c_str() << nl
2034  << errorMsg.c_str() << nl << nl
2035  << "Exiting dry-run" << nl << endl;
2036 
2037  FatalError.clear();
2038  FatalIOError.clear();
2039 
2040  return 0;
2041  }
2042  }
2043 
2044 
2045  Info<< "End\n" << endl;
2046 
2047  return 0;
2048 }
2049 
2050 
2051 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Foam::Pstream::mapCombineGather
static void mapCombineGather(const List< commsStruct > &comms, Container &Values, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:545
shellSurfaces.H
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:63
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:191
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:63
Foam::setf
Smanip< ios_base::fmtflags > setf(const ios_base::fmtflags flags)
Definition: IOmanip.H:169
profiling.H
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
Foam::HashTable::size
label size() const noexcept
Definition: HashTableI.H:45
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::polyMesh::points
virtual const pointField & points() const
Definition: polyMesh.C:1062
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:46
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:165
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:190
Foam::boundBox::mag
scalar mag() const
Definition: boundBoxI.H:126
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Foam::refinementSurfaces::maxLevel
const labelList & maxLevel() const
Definition: refinementSurfaces.H:208
Foam::meshRefinement::WRITEMESH
@ WRITEMESH
Definition: meshRefinement.H:110
Foam::surfaceZonesInfo::faceZoneType
faceZoneType
Definition: surfaceZonesInfo.H:82
Foam::accessOp
Definition: UList.H:683
snappySnapDriver.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::fileName
A class for handling file names.
Definition: fileName.H:71
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:57
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:88
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputSpray.H:25
uindirectPrimitivePatch.H
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:59
Foam::argList::getOrDefault
T getOrDefault(const word &optName, const T &deflt) const
Definition: argListI.H:300
Foam::meshRefinement::lookup
static ITstream & lookup(const dictionary &dict, const word &keyword, const bool noExit, enum keyType::option matchOpt=keyType::REGEX)
Definition: meshRefinement.C:3477
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
Foam::DynamicList< label >
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryI.H:80
Foam::refinementSurfaces::surfaces
const labelList & surfaces() const
Definition: refinementSurfaces.H:179
setFormat
const word setFormat(propsDict.getOrDefault< word >("setFormat", "vtk"))
Foam::Time::writeFormat
IOstream::streamFormat writeFormat() const
Definition: Time.H:383
addOverwriteOption.H
globalIndex.H
Foam::TimePaths::processorCase
bool processorCase() const noexcept
Definition: TimePathsI.H:29
Foam::Warning
messageStream Warning
Foam::meshRefinement::writeLevel
static writeType writeLevel()
Definition: meshRefinement.C:3421
Foam::fvMeshTools::removeEmptyPatches
static labelList removeEmptyPatches(fvMesh &, const bool validBoundary)
Definition: fvMeshTools.C:354
dictName
const word dictName("faMeshDefinition")
Foam::meshRefinement::printMeshInfo
void printMeshInfo(const bool, const string &) const
Definition: meshRefinement.C:3170
Foam::meshRefinement::addFaceZone
label addFaceZone(const word &fzName, const word &masterPatch, const word &slavePatch, const surfaceZonesInfo::faceZoneType &fzType)
Definition: meshRefinement.C:2385
Foam::searchableSurfaces::checkTopology
label checkTopology(const bool report) const
Definition: searchableSurfaces.C:813
polyTopoChange.H
motionSmoother.H
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:773
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:59
wallPolyPatch.H
Foam::meshRefinement::subDict
static const dictionary & subDict(const dictionary &dict, const word &keyword, const bool noExit, enum keyType::option matchOpt=keyType::REGEX)
Definition: meshRefinement.C:3446
Foam::argList::addNote
static void addNote(const string &note)
Definition: argList.C:405
Foam::HashTable::insert
bool insert(const Key &key, const T &obj)
Definition: HashTableI.H:173
Foam::Map
A HashTable to objects of type <T> with a label key.
Definition: lumpedPointController.H:65
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Definition: UPstream.H:453
Foam::faceSet
A list of face labels.
Definition: faceSet.H:47
Foam::surfaceZonesInfo::getNamedSurfaces
static labelList getNamedSurfaces(const PtrList< surfaceZonesInfo > &surfList)
Definition: surfaceZonesInfo.C:256
Foam::vtkSetWriter
Definition: vtkSetWriter.H:50
Foam::entry::keyword
const keyType & keyword() const noexcept
Definition: entry.H:191
Foam::FatalIOError
IOerror FatalIOError
Foam::cpuTimePosix::elapsedCpuTime
double elapsedCpuTime() const
Definition: cpuTimePosix.C:73
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Definition: polyMesh.H:440
Foam::motionSmootherAlgo::checkMesh
static bool checkMesh(const bool report, const polyMesh &mesh, const dictionary &dict, labelHashSet &wrongFaces, const bool dryRun=false)
Definition: motionSmootherAlgoCheck.C:455
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::searchableSurfaces::writeStats
void writeStats(const List< wordList > &, Ostream &) const
Definition: searchableSurfaces.C:863
Foam::PtrList::set
const T * set(const label i) const
Definition: PtrList.H:134
Foam::snappySnapDriver
All to do with snapping to surface.
Definition: snappySnapDriver.H:54
decompositionMethod.H
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:100
Foam::dictionary::set
entry * set(entry *entryPtr)
Definition: dictionary.C:773
meshBb
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(rndGen, 1e-3))
refinementFeatures.H
Foam::cpuTimePosix
Starts timing CPU usage and return elapsed time from start.
Definition: cpuTimePosix.H:52
Foam::Pout
prefixOSstream Pout
Foam::HashSet
A HashTable with keys but without contents that is similar to std::unordered_set.
Definition: HashSet.H:73
Foam::meshRefinement::mesh
const fvMesh & mesh() const
Definition: meshRefinement.H:940
processorMeshes.H
Foam::incrIndent
Ostream & incrIndent(Ostream &os)
Definition: Ostream.H:352
Foam::cpuTimePosix::cpuTimeIncrement
double cpuTimeIncrement() const
Definition: cpuTimePosix.C:80
Foam::wordRe
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition: wordRe.H:78
Foam::refinementSurfaces::gapLevel
const labelList & gapLevel() const
Definition: refinementSurfaces.H:214
setSystemMeshDictionaryIO.H
Foam::inplaceRenumber
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Definition: ListOpsTemplates.C:54
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:73
regionNames
wordList regionNames
Definition: getAllRegionOptions.H:31
Foam::regIOobject::store
bool store()
Definition: regIOobjectI.H:30
Foam::sumOp
Definition: ops.H:207
Foam::surfZoneIdentifier
Identifies a surface patch/zone by name and index, with optional geometric type.
Definition: surfZoneIdentifier.H:55
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
patchTypes
wordList patchTypes(nPatches)
Foam::meshRefinement::FaceMergeType
FaceMergeType
Definition: meshRefinement.H:129
Foam::globalMeshData::mergePoints
autoPtr< globalIndex > mergePoints(labelList &pointToGlobal, labelList &uniquePoints) const
Definition: globalMeshData.C:2367
Foam::decompositionModel::canonicalName
static const word canonicalName
Definition: decompositionModel.H:75
UnsortedMeshedSurface.H
searchableSurfaces.H
Foam::snapParameters
Simple container to keep together snap specific information.
Definition: snapParameters.H:48
Foam::decompositionMethod::parallelAware
virtual bool parallelAware() const =0
Foam::refinementSurfaces::surfZones
const PtrList< surfaceZonesInfo > & surfZones() const
Definition: refinementSurfaces.H:190
Foam::polyBoundaryMesh::checkParallelSync
bool checkParallelSync(const bool report=false) const
Definition: polyBoundaryMesh.C:965
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:45
Foam::argList::noFunctionObjects
static void noFunctionObjects(bool addWithOption=false)
Definition: argList.C:466
Foam::meshRefinement::OBJINTERSECTIONS
@ OBJINTERSECTIONS
Definition: meshRefinement.H:91
Foam::primitiveMesh::nCells
label nCells() const noexcept
Definition: primitiveMeshI.H:89
Foam::shellSurfaces
Encapsulates queries for volume refinement ('refine all cells within shell').
Definition: shellSurfaces.H:53
snappyVoxelMeshDriver.H
Foam::meshRefinement::timeName
word timeName() const
Definition: meshRefinement.C:3223
Foam::refinementSurfaces::names
const wordList & names() const
Definition: refinementSurfaces.H:185
Foam::argList::dryRun
int dryRun() const noexcept
Definition: argListI.H:109
Foam::Field
Generic templated field type.
Definition: Field.H:59
Foam::meshRefinement::checkCoupledFaceZones
static void checkCoupledFaceZones(const polyMesh &)
Definition: meshRefinement.C:1988
Foam::meshRefinement::readFlags
static int readFlags(const EnumContainer &namedEnum, const wordList &words)
Definition: meshRefinementTemplates.C:246
Foam::UnsortedMeshedSurface
A surface geometry mesh, in which the surface zone information is conveyed by the 'zoneId' associated...
Definition: MeshedSurface.H:79
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:64
refinementSurfaces.H
addProfilingOption.H
Foam::List::setSize
void setSize(const label n)
Definition: List.H:218
argList.H
Foam::searchableSurface
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
Definition: searchableSurface.H:65
Foam::meshRefinement::meshedPatches
labelList meshedPatches() const
Definition: meshRefinement.C:2358
faceSet.H
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:453
Foam::refinementSurfaces::minLevel
const labelList & minLevel() const
Definition: refinementSurfaces.H:202
addRegionOption.H
IOmanip.H
Istream and Ostream manipulators taking arguments.
Foam::refinementParameters
Simple container to keep together refinement specific information.
Definition: refinementParameters.H:54
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const noexcept
Definition: polyMesh.H:482
Foam::IOobject::selectIO
static IOobject selectIO(const IOobject &io, const fileName &altFile, const word &ioName="")
Definition: IOobject.C:231
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:55
size_type
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:69
Foam::fvMeshTools::reorderPatches
static void reorderPatches(fvMesh &, const labelList &oldToNew, const label nPatches, const bool validBoundary)
Definition: fvMeshTools.C:322
Foam::surfaceZonesInfo
Definition: surfaceZonesInfo.H:56
Foam::checkGeometry
label checkGeometry(const polyMesh &mesh, const bool allGeometry, autoPtr< surfaceWriter > &surfWriter, const autoPtr< writer< scalar >> &setWriter)
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:68
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Definition: hashSets.C:40
Foam::layerParameters
Simple container to keep together layer specific information.
Definition: layerParameters.H:54
Foam::meshRefinement::writeTypeNames
static const Enum< writeType > writeTypeNames
Definition: meshRefinement.H:117
addProfiling
#define addProfiling(name, descr)
Definition: profilingTrigger.H:112
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:119
createNamedMesh.H
Required Variables.
Foam::refinementSurfaces::globalRegion
label globalRegion(const label surfI, const label regionI) const
Definition: refinementSurfaces.H:266
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:255
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::decompositionMethod
Abstract base class for domain decomposition.
Definition: decompositionMethod.H:47
Foam::writer
Base class for graphics format writing. Entry points are.
Definition: writer.H:77
dictIO
IOobject dictIO
Definition: setConstantMeshDictionaryIO.H:1
snappyRefineDriver.H
Foam::Ostream::write
virtual bool write(const token &tok)=0
Foam::error::message
string message() const
Definition: error.C:312
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:81
fvMesh.H
Foam
Definition: atmBoundaryLayer.C:26
Foam::decrIndent
Ostream & decrIndent(Ostream &os)
Definition: Ostream.H:361
Foam::patchIdentifier::index
label index() const noexcept
Definition: patchIdentifier.H:143
Foam::refinementHistory::removeFiles
static void removeFiles(const polyMesh &)
Definition: refinementHistory.C:1720
Foam::setw
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Foam::entry::dict
virtual const dictionary & dict() const =0
Foam::indent
Ostream & indent(Ostream &os)
Definition: Ostream.H:343
Foam::polyPatch::start
label start() const
Definition: polyPatch.H:357
Foam::searchableSurfaces::regionNames
const List< wordList > & regionNames() const
Definition: searchableSurfaces.H:164
Foam::argList::addDryRunOption
static void addDryRunOption(const string &usage, bool advanced=false)
Definition: argList.C:447
Foam::dictionary::subOrEmptyDict
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Definition: dictionary.C:533
Foam::refinementSurfaces::patchInfo
const PtrList< dictionary > & patchInfo() const
Definition: refinementSurfaces.H:257
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Foam::snappyLayerDriver
All to do with adding layers.
Definition: snappyLayerDriver.H:55
Foam::HashTable< label >
Foam::polyMesh::bounds
const boundBox & bounds() const
Definition: polyMesh.H:446
Foam::autoPtr::reset
void reset(autoPtr< T > &&other) noexcept
Foam::Perr
prefixOSstream Perr
Foam::IOobject::name
const word & name() const noexcept
Definition: IOobjectI.H:58
Foam::argList::addBoolOption
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Definition: argList.C:317
Foam::IOstream::defaultPrecision
static unsigned int defaultPrecision() noexcept
Definition: IOstream.H:338
Time.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:49
setRootCase.H
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:78
vtkSetWriter.H
Foam::searchableSurfaces::checkGeometry
label checkGeometry(const scalar maxRatio, const scalar tolerance, const autoPtr< writer< scalar >> &setWriter, const scalar minQuality, const bool report) const
Definition: searchableSurfaces.C:833
Foam::TimePaths::system
const word & system() const
Definition: TimePathsI.H:95
Foam::IOstreamOption::ASCII
@ ASCII
"ascii" (normal default)
Definition: IOstreamOption.H:68
Foam::polyMesh::faces
virtual const faceList & faces() const
Definition: polyMesh.C:1087
FatalErrorInFunction
#define FatalErrorInFunction
Definition: error.H:465
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Definition: gatherScatterList.C:46
fvMeshDistribute.H
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Definition: UPstream.H:459
Foam::nl
constexpr char nl
Definition: Ostream.H:424
Foam::Pair< word >
Foam::Time::path
fileName path() const
Definition: Time.H:354
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::Pstream::mapCombineScatter
static void mapCombineScatter(const List< commsStruct > &comms, Container &Values, const int tag, const label comm)
Definition: combineGatherScatter.C:660
fvMeshTools.H
Foam::ILList::erase
bool erase(T *item)
Definition: ILList.C:90
Foam::refinementFeatures
Encapsulates queries for features.
Definition: refinementFeatures.H:49
f
labelList f(nPoints)
Foam::foamVersion::patch
const std::string patch
layerParameters.H
refinementParameters.H
Foam::Vector< scalar >
Foam::regIOobject::rename
virtual void rename(const word &newName)
Definition: regIOobject.C:408
Foam::meshRefinement::writeType
writeType
Definition: meshRefinement.H:108
Foam::UPstream::parRun
static bool & parRun() noexcept
Definition: UPstream.H:429
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:58
Foam::OStringStream
Definition: StringStream.H:229
Foam::searchableSurfaces
Container for searchableSurfaces. The collection is specified as a dictionary. For example,...
Definition: searchableSurfaces.H:88
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::refinementSurfaces::nRegions
label nRegions() const
Definition: refinementSurfaces.H:286
Foam::timer
Implements a timeout mechanism via sigalarm.
Definition: timer.H:82
Foam::meshRefinement
Helper class which maintains intersections of (changing) mesh with (static) surfaces.
Definition: meshRefinement.H:81
Foam::processorMeshes::removeFiles
static void removeFiles(const polyMesh &mesh)
Definition: processorMeshes.C:267
Foam::dictionary::findEntry
entry * findEntry(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Definition: dictionaryI.H:90
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Definition: createFields.H:11
Foam::dictionary::clone
autoPtr< dictionary > clone() const
Definition: dictionary.C:165
Foam::meshRefinement::write
bool write() const
Definition: meshRefinement.C:3063
Foam::keyType::REGEX
@ REGEX
Regular expression.
Definition: keyType.H:80
Foam::identity
labelList identity(const label len, label start=0)
Definition: labelList.C:31
Foam::decompositionMethod::New
static autoPtr< decompositionMethod > New(const dictionary &decompDict, const word &regionName="")
Definition: decompositionMethod.C:337
Foam::direction
uint8_t direction
Definition: direction.H:46
Foam::word::null
static const word null
Definition: word.H:78
Foam::error::clear
void clear() const
Definition: error.C:318
createTime.H
Foam::meshRefinement::addMeshedPatch
label addMeshedPatch(const word &name, const dictionary &)
Definition: meshRefinement.C:2308
Foam::IOobject::MUST_READ_IF_MODIFIED
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:182
snapParameters.H
Foam::line
A line primitive.
Definition: line.H:49
Foam::argList::getList
List< T > getList(const label index) const
Foam::snappyRefineDriver::addFaceZones
static void addFaceZones(meshRefinement &meshRefiner, const refinementParameters &refineParams, const HashTable< Pair< word >> &faceZoneToPatches)
Definition: snappyRefineDriver.C:2986
Foam::dictionary::optionalSubDict
const dictionary & optionalSubDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:570
Foam::refinementSurfaces::setMinLevelFields
void setMinLevelFields(const shellSurfaces &shells)
Definition: refinementSurfaces.C:700
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:57
Foam::TimePaths::globalPath
fileName globalPath() const
Definition: TimePathsI.H:73
Foam::name
word name(const expressions::valueTypeCode typeCode)
Definition: exprTraits.C:52
Foam::fvMesh::time
const Time & time() const
Definition: fvMesh.H:276
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:56
Foam::dictionary::add
entry * add(entry *entryPtr, bool mergeEntry=false)
Definition: dictionary.C:633
Foam::primitiveMesh::nFaces
label nFaces() const noexcept
Definition: primitiveMeshI.H:83
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
decompositionModel.H
Foam::DelaunayMeshTools::allPoints
tmp< pointField > allPoints(const Triangulation &t)
Foam::meshRefinement::updateIntersections
void updateIntersections(const labelList &changedFaces)
Definition: meshRefinement.C:291
Foam::patchIdentifier::name
const word & name() const noexcept
Definition: patchIdentifier.H:131
Foam::plusEqOp
Definition: ops.H:66
Foam::dictionary::toc
wordList toc() const
Definition: dictionary.C:595
Foam::searchableSurfaces::names
const wordList & names() const
Definition: searchableSurfaces.H:152
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:52
IOWarningInFunction
#define IOWarningInFunction(ios)
Definition: messageStream.H:371
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Definition: polyMesh.C:1288
Foam::meshRefinement::debugTypeNames
static const Enum< debugType > debugTypeNames
Definition: meshRefinement.H:97
Foam::TimePaths::constant
const word & constant() const
Definition: TimePathsI.H:89
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:141
Foam::MeshedSurface< face >
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:184
Foam::argList::addOption
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Definition: argList.C:328
args
Foam::argList args(argc, argv)
Foam::meshRefinement::WRITELAYERSETS
@ WRITELAYERSETS
Definition: meshRefinement.H:113
Foam::orOp
Definition: ops.H:228
Foam::refinementSurfaces
Container for data on surfaces used for surface-driven refinement. Contains all the data about the le...
Definition: refinementSurfaces.H:59
Foam::PtrListOps::names
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
Foam::topoSet::removeFiles
static void removeFiles(const polyMesh &)
Definition: topoSet.C:628
Foam::snappyRefineDriver
Definition: snappyRefineDriver.H:61
meshDict
const IOdictionary & meshDict
Definition: findBlockMeshDict.H:92
Foam::fvMeshDistribute
Sends/receives parts of mesh+fvfields to neighbouring processors. Used in load balancing.
Definition: fvMeshDistribute.H:66
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Definition: UPstream.H:441
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:398
Foam::meshRefinement::getFaceZoneInfo
bool getFaceZoneInfo(const word &fzName, label &masterPatchID, label &slavePatchID, surfaceZonesInfo::faceZoneType &fzType) const
Definition: meshRefinement.C:2409
Foam::argList::found
bool found(const word &optName) const
Definition: argListI.H:171
MeshedSurface.H
snappyLayerDriver.H
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:75
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Definition: fvMeshGeometry.C:172
Foam::surfaceZonesInfo::addCellZonesToMesh
static labelList addCellZonesToMesh(const PtrList< surfaceZonesInfo > &surfList, const labelList &namedSurfaces, polyMesh &mesh)
Definition: surfaceZonesInfo.C:451
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:181
Foam::meshRefinement::debugType
debugType
Definition: meshRefinement.H:88
surfZoneIdentifierList.H
Foam::profiling::writeNow
static bool writeNow()
Definition: profiling.C:127