snappyHexMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  snappyHexMesh
26 
27 Description
28  Automatic split hex mesher
29 
30  Refines, snaps to surface and adds surface layers.
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #include "argList.H"
35 #include "fvMesh.H"
36 #include "snappyRefineDriver.H"
37 #include "snappySnapDriver.H"
38 #include "snappyLayerDriver.H"
39 #include "searchableSurfaceList.H"
40 #include "refinementSurfaces.H"
41 #include "refinementFeatures.H"
42 #include "refinementRegions.H"
43 #include "decompositionMethod.H"
44 #include "fvMeshDistribute.H"
45 #include "wallPolyPatch.H"
46 #include "snapParameters.H"
47 #include "layerParameters.H"
48 #include "faceSet.H"
49 #include "meshCheck.H"
51 #include "MeshedSurface.H"
52 #include "IOmanip.H"
53 #include "fvMeshTools.H"
54 #include "systemDict.H"
55 
56 using namespace Foam;
57 
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 
60 // Convert size (as fraction of defaultCellSize) to refinement level
61 label sizeCoeffToRefinement
62 (
63  const scalar level0Coeff, // ratio of hex cell size v.s. defaultCellSize
64  const scalar sizeCoeff
65 )
66 {
67  return round(::log(level0Coeff/sizeCoeff)/::log(2));
68 }
69 
70 
71 autoPtr<refinementSurfaces> createRefinementSurfaces
72 (
73  const searchableSurfaceList& allGeometry,
74  const dictionary& surfacesDict,
75  const dictionary& shapeControlDict,
76  const label gapLevelIncrement,
77  const scalar level0Coeff
78 )
79 {
80  autoPtr<refinementSurfaces> surfacePtr;
81 
82  // Count number of surfaces.
83  label surfI = 0;
84  forAll(allGeometry.names(), geomI)
85  {
86  const word& geomName = allGeometry.names()[geomI];
87 
88  if (surfacesDict.found(geomName))
89  {
90  surfI++;
91  }
92  }
93 
94  labelList surfaces(surfI);
95  wordList names(surfI);
96  PtrList<surfaceZonesInfo> surfZones(surfI);
97 
98  labelList regionOffset(surfI);
99 
100  labelList globalMinLevel(surfI, 0);
101  labelList globalMaxLevel(surfI, 0);
102  labelList globalLevelIncr(surfI, 0);
103  PtrList<dictionary> globalPatchInfo(surfI);
104  List<Map<label>> regionMinLevel(surfI);
105  List<Map<label>> regionMaxLevel(surfI);
106  List<Map<label>> regionLevelIncr(surfI);
107  List<Map<scalar>> regionAngle(surfI);
108  List<Map<autoPtr<dictionary>>> regionPatchInfo(surfI);
109 
110  HashSet<word> unmatchedKeys(surfacesDict.toc());
111 
112  surfI = 0;
113  forAll(allGeometry.names(), geomI)
114  {
115  const word& geomName = allGeometry.names()[geomI];
116 
117  const entry* ePtr = surfacesDict.lookupEntryPtr(geomName, false, true);
118 
119  if (ePtr)
120  {
121  const dictionary& shapeDict = ePtr->dict();
122  unmatchedKeys.erase(ePtr->keyword());
123 
124  names[surfI] = geomName;
125  surfaces[surfI] = geomI;
126 
127  const searchableSurface& surface = allGeometry[geomI];
128 
129  // Find the index in shapeControlDict
130  // Invert surfaceCellSize to get the refinementLevel
131 
132  const word scsFuncName =
133  shapeDict.lookup("surfaceCellSizeFunction");
134  const dictionary& scsDict =
135  shapeDict.optionalSubDict(scsFuncName + "Coeffs");
136 
137  const scalar surfaceCellSize =
138  scsDict.lookup<scalar>("surfaceCellSizeCoeff");
139 
140  const label refLevel = sizeCoeffToRefinement
141  (
142  level0Coeff,
143  surfaceCellSize
144  );
145 
146  globalMinLevel[surfI] = refLevel;
147  globalMaxLevel[surfI] = refLevel;
148  globalLevelIncr[surfI] = gapLevelIncrement;
149 
150  // Surface zones
151  surfZones.set(surfI, new surfaceZonesInfo(surface, shapeDict));
152 
153 
154  // Global perpendicular angle
155  if (shapeDict.found("patchInfo"))
156  {
157  globalPatchInfo.set
158  (
159  surfI,
160  shapeDict.subDict("patchInfo").clone()
161  );
162  }
163 
164 
165  // Per region override of patchInfo
166 
167  if (shapeDict.found("regions"))
168  {
169  const dictionary& regionsDict = shapeDict.subDict("regions");
170  const wordList& regionNames =
171  allGeometry[surfaces[surfI]].regions();
172 
173  forAll(regionNames, regionI)
174  {
175  if (regionsDict.found(regionNames[regionI]))
176  {
177  // Get the dictionary for region
178  const dictionary& regionDict = regionsDict.subDict
179  (
180  regionNames[regionI]
181  );
182 
183  if (regionDict.found("patchInfo"))
184  {
185  regionPatchInfo[surfI].insert
186  (
187  regionI,
188  regionDict.subDict("patchInfo").clone()
189  );
190  }
191  }
192  }
193  }
194 
195  // Per region override of cellSize
196  if (shapeDict.found("regions"))
197  {
198  const dictionary& shapeControlRegionsDict =
199  shapeDict.subDict("regions");
200  const wordList& regionNames =
201  allGeometry[surfaces[surfI]].regions();
202 
203  forAll(regionNames, regionI)
204  {
205  if (shapeControlRegionsDict.found(regionNames[regionI]))
206  {
207  const dictionary& shapeControlRegionDict =
208  shapeControlRegionsDict.subDict
209  (
210  regionNames[regionI]
211  );
212 
213  const word scsFuncName =
214  shapeControlRegionDict.lookup
215  (
216  "surfaceCellSizeFunction"
217  );
218  const dictionary& scsDict =
219  shapeControlRegionDict.subDict
220  (
221  scsFuncName + "Coeffs"
222  );
223 
224  const scalar surfaceCellSize =
225  scsDict.lookup<scalar>("surfaceCellSizeCoeff");
226 
227  const label refLevel = sizeCoeffToRefinement
228  (
229  level0Coeff,
230  surfaceCellSize
231  );
232 
233  regionMinLevel[surfI].insert(regionI, refLevel);
234  regionMaxLevel[surfI].insert(regionI, refLevel);
235  regionLevelIncr[surfI].insert(regionI, 0);
236  }
237  }
238  }
239 
240  surfI++;
241  }
242  }
243 
244  // Calculate local to global region offset
245  label nRegions = 0;
246 
247  forAll(surfaces, surfI)
248  {
249  regionOffset[surfI] = nRegions;
250  nRegions += allGeometry[surfaces[surfI]].regions().size();
251  }
252 
253  // Rework surface specific information into information per global region
254  labelList minLevel(nRegions, 0);
255  labelList maxLevel(nRegions, 0);
256  labelList gapLevel(nRegions, -1);
257  PtrList<dictionary> patchInfo(nRegions);
258 
259  forAll(globalMinLevel, surfI)
260  {
261  label nRegions = allGeometry[surfaces[surfI]].regions().size();
262 
263  // Initialise to global (i.e. per surface)
264  for (label i = 0; i < nRegions; i++)
265  {
266  label globalRegionI = regionOffset[surfI] + i;
267  minLevel[globalRegionI] = globalMinLevel[surfI];
268  maxLevel[globalRegionI] = globalMaxLevel[surfI];
269  gapLevel[globalRegionI] =
270  maxLevel[globalRegionI]
271  + globalLevelIncr[surfI];
272 
273  if (globalPatchInfo.set(surfI))
274  {
275  patchInfo.set
276  (
277  globalRegionI,
278  globalPatchInfo[surfI].clone()
279  );
280  }
281  }
282 
283  // Overwrite with region specific information
284  forAllConstIter(Map<label>, regionMinLevel[surfI], iter)
285  {
286  label globalRegionI = regionOffset[surfI] + iter.key();
287 
288  minLevel[globalRegionI] = iter();
289  maxLevel[globalRegionI] = regionMaxLevel[surfI][iter.key()];
290  gapLevel[globalRegionI] =
291  maxLevel[globalRegionI]
292  + regionLevelIncr[surfI][iter.key()];
293  }
294 
295  const Map<autoPtr<dictionary>>& localInfo = regionPatchInfo[surfI];
296  forAllConstIter(Map<autoPtr<dictionary>>, localInfo, iter)
297  {
298  label globalRegionI = regionOffset[surfI] + iter.key();
299  patchInfo.set(globalRegionI, iter()().clone());
300  }
301  }
302 
303  surfacePtr.set
304  (
306  (
307  allGeometry,
308  surfaces,
309  names,
310  surfZones,
311  regionOffset,
312  minLevel,
313  maxLevel,
314  gapLevel,
315  scalarField(nRegions, -great), // perpendicularAngle,
316  patchInfo
317  )
318  );
319 
320 
321  const refinementSurfaces& rf = surfacePtr();
322 
323  // Determine maximum region name length
324  label maxLen = 0;
325  forAll(rf.surfaces(), surfI)
326  {
327  label geomI = rf.surfaces()[surfI];
328  const wordList& regionNames = allGeometry.regionNames()[geomI];
329  forAll(regionNames, regionI)
330  {
331  maxLen = Foam::max(maxLen, label(regionNames[regionI].size()));
332  }
333  }
334 
335 
336  Info<< setw(maxLen) << "Region"
337  << setw(10) << "Min Level"
338  << setw(10) << "Max Level"
339  << setw(10) << "Gap Level" << nl
340  << setw(maxLen) << "------"
341  << setw(10) << "---------"
342  << setw(10) << "---------"
343  << setw(10) << "---------" << endl;
344 
345  forAll(rf.surfaces(), surfI)
346  {
347  label geomI = rf.surfaces()[surfI];
348 
349  Info<< rf.names()[surfI] << ':' << nl;
350 
351  const wordList& regionNames = allGeometry.regionNames()[geomI];
352 
353  forAll(regionNames, regionI)
354  {
355  label globalI = rf.globalRegion(surfI, regionI);
356 
357  Info<< setw(maxLen) << regionNames[regionI]
358  << setw(10) << rf.minLevel()[globalI]
359  << setw(10) << rf.maxLevel()[globalI]
360  << setw(10) << rf.gapLevel()[globalI] << endl;
361  }
362  }
363 
364 
365  return surfacePtr;
366 }
367 
368 
369 void extractSurface
370 (
371  const polyMesh& mesh,
372  const Time& runTime,
374  const fileName& outFileName
375 )
376 {
378 
379  // Collect sizes. Hash on names to handle local-only patches (e.g.
380  // processor patches)
381  HashTable<label> patchSize(1000);
382  label nFaces = 0;
384  {
385  const polyPatch& pp = bMesh[iter.key()];
386  patchSize.insert(pp.name(), pp.size());
387  nFaces += pp.size();
388  }
390 
391 
392  // Allocate zone/patch for all patches
393  HashTable<label> compactZoneID(1000);
394  forAllConstIter(HashTable<label>, patchSize, iter)
395  {
396  label sz = compactZoneID.size();
397  compactZoneID.insert(iter.key(), sz);
398  }
399  Pstream::mapCombineScatter(compactZoneID);
400 
401 
402  // Rework HashTable into labelList just for speed of conversion
403  labelList patchToCompactZone(bMesh.size(), -1);
404  forAllConstIter(HashTable<label>, compactZoneID, iter)
405  {
406  label patchi = bMesh.findIndex(iter.key());
407  if (patchi != -1)
408  {
409  patchToCompactZone[patchi] = iter();
410  }
411  }
412 
413  // Collect faces on zones
414  DynamicList<label> faceLabels(nFaces);
415  DynamicList<label> compactZones(nFaces);
417  {
418  const polyPatch& pp = bMesh[iter.key()];
419  forAll(pp, i)
420  {
421  faceLabels.append(pp.start()+i);
422  compactZones.append(patchToCompactZone[pp.index()]);
423  }
424  }
425 
426  // Addressing engine for all faces
427  uindirectPrimitivePatch allBoundary
428  (
429  UIndirectList<face>(mesh.faces(), faceLabels),
430  mesh.points()
431  );
432 
433 
434  // Find correspondence to master points
435  labelList pointToGlobal;
436  labelList uniqueMeshPoints;
438  (
439  allBoundary.meshPoints(),
440  allBoundary.meshPointMap(),
441  pointToGlobal,
442  uniqueMeshPoints
443  );
444 
445  // Gather all unique points on master
446  List<pointField> gatheredPoints(Pstream::nProcs());
447  gatheredPoints[Pstream::myProcNo()] = pointField
448  (
449  mesh.points(),
450  uniqueMeshPoints
451  );
452  Pstream::gatherList(gatheredPoints);
453 
454  // Gather all faces
455  List<faceList> gatheredFaces(Pstream::nProcs());
456  gatheredFaces[Pstream::myProcNo()] = allBoundary.localFaces();
457  forAll(gatheredFaces[Pstream::myProcNo()], i)
458  {
459  inplaceRenumber(pointToGlobal, gatheredFaces[Pstream::myProcNo()][i]);
460  }
461  Pstream::gatherList(gatheredFaces);
462 
463  // Gather all ZoneIDs
464  List<labelList> gatheredZones(Pstream::nProcs());
465  gatheredZones[Pstream::myProcNo()] = move(compactZones);
466  Pstream::gatherList(gatheredZones);
467 
468  // On master combine all points, faces, zones
469  if (Pstream::master())
470  {
471  pointField allPoints = ListListOps::combine<pointField>
472  (
473  gatheredPoints,
475  );
476  gatheredPoints.clear();
477 
478  faceList allFaces = ListListOps::combine<faceList>
479  (
480  gatheredFaces,
482  );
483  gatheredFaces.clear();
484 
485  labelList allZones = ListListOps::combine<labelList>
486  (
487  gatheredZones,
489  );
490  gatheredZones.clear();
491 
492 
493  // Zones
494  surfZoneIdentifierList surfZones(compactZoneID.size());
495  forAllConstIter(HashTable<label>, compactZoneID, iter)
496  {
497  surfZones[iter()] = surfZoneIdentifier(iter.key(), iter());
498  Info<< "surfZone " << iter() << " : " << surfZones[iter()].name()
499  << endl;
500  }
501 
502  UnsortedMeshedSurface<face> unsortedFace
503  (
504  move(allPoints),
505  move(allFaces),
506  move(allZones),
507  move(surfZones)
508  );
509 
510 
511  MeshedSurface<face> sortedFace(unsortedFace);
512 
513  Info<< "Writing merged surface to "
514  << runTime.globalPath()/outFileName << endl;
515 
516  sortedFace.write(runTime.globalPath()/outFileName);
517  }
518 }
519 
520 
521 // Check writing tolerance before doing any serious work
522 scalar getMergeDistance(const polyMesh& mesh, const scalar mergeTol)
523 {
524  const boundBox& meshBb = mesh.bounds();
525  scalar mergeDist = mergeTol * meshBb.mag();
526 
527  Info<< nl
528  << "Overall mesh bounding box : " << meshBb << nl
529  << "Relative tolerance : " << mergeTol << nl
530  << "Absolute matching distance : " << mergeDist << nl
531  << endl;
532 
533  // check writing tolerance
535  {
536  const scalar writeTol = std::pow
537  (
538  scalar(10.0),
539  -scalar(IOstream::defaultPrecision())
540  );
541 
542  if (mergeTol < writeTol)
543  {
545  << "Your current settings specify ASCII writing with "
546  << IOstream::defaultPrecision() << " digits precision." << nl
547  << "Your merging tolerance (" << mergeTol
548  << ") is finer than this." << nl
549  << "Change to binary writeFormat, "
550  << "or increase the writePrecision" << endl
551  << "or adjust the merge tolerance (mergeTol)."
552  << exit(FatalError);
553  }
554  }
555 
556  return mergeDist;
557 }
558 
559 
560 void removeZeroSizedPatches(fvMesh& mesh)
561 {
562  // Remove non-constraint zero-sized patches
563 
564  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
565 
566  labelList oldToNew(pbm.size(), -1);
567  label newPatchi = 0;
568  forAll(pbm, patchi)
569  {
570  const polyPatch& pp = pbm[patchi];
571 
572  if
573  (
574  polyPatch::constraintType(pp.type())
575  || returnReduce(pp.size(), sumOp<label>())
576  )
577  {
578  oldToNew[patchi] = newPatchi++;
579  }
580  }
581 
582  const label nKeepPatches = newPatchi;
583 
584  // Shuffle unused ones to end
585  if (nKeepPatches != pbm.size())
586  {
587  Info<< endl
588  << "Removing zero-sized patches:" << endl << incrIndent;
589 
590  forAll(oldToNew, patchi)
591  {
592  if (oldToNew[patchi] == -1)
593  {
594  Info<< indent << pbm[patchi].name()
595  << " type " << pbm[patchi].type()
596  << " at position " << patchi << endl;
597  oldToNew[patchi] = newPatchi++;
598  }
599  }
600  Info<< decrIndent;
601 
602  fvMeshTools::reorderPatches(mesh, oldToNew, nKeepPatches, true);
603  Info<< endl;
604  }
605 }
606 
607 
608 // Write mesh and additional information
609 void writeMesh
610 (
611  const string& msg,
612  const meshRefinement& meshRefiner,
613  const meshRefinement::debugType debugLevel,
614  const meshRefinement::writeType writeLevel
615 )
616 {
617  const fvMesh& mesh = meshRefiner.mesh();
618 
619  meshRefiner.printMeshInfo(debugLevel, msg);
620  Info<< "Writing mesh to time " << meshRefiner.name() << endl;
621 
622  meshRefiner.write
623  (
624  debugLevel,
626  mesh.time().path()/meshRefiner.name()
627  );
628  Info<< "Wrote mesh in = "
629  << mesh.time().cpuTimeIncrement() << " s." << endl;
630 }
631 
632 
633 int main(int argc, char *argv[])
634 {
635  Foam::argList::removeOption("noFunctionObjects");
636  #include "addNoOverwriteOption.H"
638  (
639  "checkGeometry",
640  "check all surface geometry for quality"
641  );
643  (
644  "surfaceSimplify",
645  "boundBox",
646  "simplify the surface using snappyHexMesh starting from a boundBox"
647  );
649  (
650  "patches",
651  "(patch0 .. patchN)",
652  "only triangulate selected patches (wildcards supported)"
653  );
655  (
656  "outFile",
657  "file",
658  "name of the file to save the simplified surface to"
659  );
660  #include "addDictOption.H"
661  #include "addMeshOption.H"
662  #include "addRegionOption.H"
663 
664  #include "setRootCase.H"
667 
668  Info<< "Read mesh in = "
669  << runTime.cpuTimeIncrement() << " s" << endl;
670 
671  #include "setNoOverwrite.H"
672  const bool checkGeometry = args.optionFound("checkGeometry");
673  const bool surfaceSimplify = args.optionFound("surfaceSimplify");
674 
675  // Check that the read mesh is fully 3D
676  // as required for mesh relaxation after snapping
677  if (mesh.nSolutionD() != 3)
678  {
680  << "Mesh provided is not fully 3D"
681  " as required for mesh relaxation after snapping" << nl
682  << "Convert all empty patches to appropriate types for a 3D mesh,"
683  " current patch types are" << nl << mesh.boundaryMesh().types()
684  << exit(FatalError);
685  }
686 
687  // Check patches and faceZones are synchronised
690 
691 
692  // Read meshing dictionary
693  const IOdictionary meshDict(systemDict("snappyHexMeshDict", args, mesh));
694 
695 
696  // all surface geometry
697  const dictionary& geometryDict = meshDict.subDict("geometry");
698 
699  // refinement parameters
700  const dictionary& refineDict = meshDict.subDict("castellatedMeshControls");
701 
702  // mesh motion and mesh quality parameters
703  const dictionary& motionDict = meshDict.subDict("meshQualityControls");
704 
705  // snap-to-surface parameters
706  const dictionary& snapDict = meshDict.subDict("snapControls");
707 
708  // absolute merge distance
709  const scalar mergeDist = getMergeDistance
710  (
711  mesh,
712  meshDict.lookup<scalar>("mergeTolerance")
713  );
714 
715  const Switch keepPatches(meshDict.lookupOrDefault("keepPatches", false));
716 
717 
718  // Read decomposePar dictionary
719  dictionary decomposeDict;
720  {
721  if (Pstream::parRun())
722  {
723  decomposeDict = decompositionMethod::decomposeParDict(runTime);
724  }
725  else
726  {
727  decomposeDict.add("method", "none");
728  decomposeDict.add("numberOfSubdomains", 1);
729  }
730  }
731 
732 
733  // Debug
734  // ~~~~~
735 
736  // Set debug level
738  (
739  meshDict.lookupOrDefault<label>
740  (
741  "debug",
742  0
743  )
744  );
745  {
746  wordList flags;
747  if (meshDict.readIfPresent("debugFlags", flags))
748  {
749  debugLevel = meshRefinement::debugType
750  (
752  (
754  flags
755  )
756  );
757  }
758  }
759  if (debugLevel > 0)
760  {
761  meshRefinement::debug = debugLevel;
762  snappyRefineDriver::debug = debugLevel;
763  snappySnapDriver::debug = debugLevel;
764  snappyLayerDriver::debug = debugLevel;
765  }
766 
767  // Set file writing level
768  {
769  wordList flags;
770  if (meshDict.readIfPresent("writeFlags", flags))
771  {
773  (
775  (
777  (
779  flags
780  )
781  )
782  );
783  }
784  }
785 
786  // Set output level
787  {
788  wordList flags;
789  if (meshDict.readIfPresent("outputFlags", flags))
790  {
792  (
794  (
796  (
798  flags
799  )
800  )
801  );
802  }
803  }
804 
805 
806  // Read geometry
807  // ~~~~~~~~~~~~~
808 
809  searchableSurfaceList allGeometry
810  (
811  IOobject
812  (
813  "abc",
814  mesh.time().constant(),
816  mesh.time(),
819  ),
820  geometryDict,
821  meshDict.lookupOrDefault("singleRegionName", true)
822  );
823 
824 
825  // Read refinement surfaces
826  // ~~~~~~~~~~~~~~~~~~~~~~~~
827 
828  Info<< "Reading refinement surfaces..." << endl;
829 
830  refinementSurfaces surfaces
831  (
832  allGeometry,
833  refineDict.found("refinementSurfaces")
834  ? refineDict.subDict("refinementSurfaces")
836  refineDict.lookupOrDefault("gapLevelIncrement", 0)
837  );
838 
839  Info<< "Read refinement surfaces in = "
840  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
841 
842 
843  // Checking only?
844 
845  if (checkGeometry)
846  {
847  // Extract patchInfo
848  List<wordList> patchTypes(allGeometry.size());
849 
850  const PtrList<dictionary>& patchInfo = surfaces.patchInfo();
851  const labelList& surfaceGeometry = surfaces.surfaces();
852  forAll(surfaceGeometry, surfI)
853  {
854  label geomI = surfaceGeometry[surfI];
855  const wordList& regNames = allGeometry.regionNames()[geomI];
856 
857  patchTypes[geomI].setSize(regNames.size());
858  forAll(regNames, regionI)
859  {
860  label globalRegionI = surfaces.globalRegion(surfI, regionI);
861 
862  if (patchInfo.set(globalRegionI))
863  {
864  patchTypes[geomI][regionI] =
865  word(patchInfo[globalRegionI].lookup("type"));
866  }
867  else
868  {
869  patchTypes[geomI][regionI] = wallPolyPatch::typeName;
870  }
871  }
872  }
873 
874  // Write some stats
875  allGeometry.writeStats(patchTypes, Info);
876  // Check topology
877  allGeometry.checkTopology(true);
878  // Check geometry
879  allGeometry.checkGeometry
880  (
881  100.0, // max size ratio
882  1e-9, // intersection tolerance
883  0.01, // min triangle quality
884  true
885  );
886 
887  return 0;
888  }
889 
890 
891 
892  // Read refinement shells
893  // ~~~~~~~~~~~~~~~~~~~~~~
894 
895  Info<< "Reading refinement regions..." << endl;
896  refinementRegions shells
897  (
898  allGeometry,
899  refineDict.found("refinementRegions")
900  ? refineDict.subDict("refinementRegions")
902  );
903  Info<< "Read refinement regions in = "
904  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
905 
906 
907 
908  // Read feature meshes
909  // ~~~~~~~~~~~~~~~~~~~
910 
911  Info<< "Reading features..." << endl;
912  refinementFeatures features
913  (
914  mesh,
915  refineDict.found("features")
916  ? refineDict.lookup("features")
918  );
919  Info<< "Read features in = "
920  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
921 
922 
923 
924  // Refinement engine
925  // ~~~~~~~~~~~~~~~~~
926 
927  Info<< nl
928  << "Determining initial surface intersections" << nl
929  << "-----------------------------------------" << nl
930  << endl;
931 
932  // Main refinement engine
933  meshRefinement meshRefiner
934  (
935  mesh,
936  refineDict,
937  mergeDist, // tolerance used in sorting coordinates
938  overwrite, // overwrite mesh files?
939  surfaces, // for surface intersection refinement
940  features, // for feature edges/point based refinement
941  shells // for volume (inside/outside) refinement
942  );
943  Info<< "Calculated surface intersections in = "
944  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
945 
946  // Some stats
947  meshRefiner.printMeshInfo(debugLevel, "Initial mesh");
948 
949  meshRefiner.write
950  (
953  mesh.time().path()/meshRefiner.name()
954  );
955 
956 
957  // Add all the surface regions as patches
958  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
959 
960  //- Global surface region to patch (non faceZone surface) or patches
961  // (faceZone surfaces)
962  labelList globalToMasterPatch;
963  labelList globalToSlavePatch;
964  {
965  Info<< nl
966  << "Adding patches for surface regions" << nl
967  << "----------------------------------" << nl
968  << endl;
969 
970  // From global region number to mesh patch.
971  globalToMasterPatch.setSize(surfaces.nRegions(), -1);
972  globalToSlavePatch.setSize(surfaces.nRegions(), -1);
973 
974  Info<< setf(ios_base::left)
975  << setw(6) << "Patch"
976  << setw(20) << "Type"
977  << setw(30) << "Region" << nl
978  << setw(6) << "-----"
979  << setw(20) << "----"
980  << setw(30) << "------" << endl;
981 
982  const labelList& surfaceGeometry = surfaces.surfaces();
983  const PtrList<dictionary>& surfacePatchInfo = surfaces.patchInfo();
984 
985  forAll(surfaceGeometry, surfI)
986  {
987  label geomI = surfaceGeometry[surfI];
988 
989  const wordList& regNames = allGeometry.regionNames()[geomI];
990 
991  Info<< surfaces.names()[surfI] << ':' << nl << nl;
992 
993  if (surfaces.surfZones()[surfI].faceZoneName().empty())
994  {
995  // 'Normal' surface
996  forAll(regNames, i)
997  {
998  label globalRegionI = surfaces.globalRegion(surfI, i);
999 
1000  label patchi;
1001 
1002  if (surfacePatchInfo.set(globalRegionI))
1003  {
1004  patchi = meshRefiner.addMeshedPatch
1005  (
1006  regNames[i],
1007  surfacePatchInfo[globalRegionI]
1008  );
1009  }
1010  else
1011  {
1012  dictionary patchInfo;
1013  patchInfo.set("type", wallPolyPatch::typeName);
1014 
1015  patchi = meshRefiner.addMeshedPatch
1016  (
1017  regNames[i],
1018  patchInfo
1019  );
1020  }
1021 
1022  Info<< setf(ios_base::left)
1023  << setw(6) << patchi
1024  << setw(20) << mesh.boundaryMesh()[patchi].type()
1025  << setw(30) << regNames[i] << nl;
1026 
1027  globalToMasterPatch[globalRegionI] = patchi;
1028  globalToSlavePatch[globalRegionI] = patchi;
1029  }
1030  }
1031  else
1032  {
1033  // Zoned surface
1034  forAll(regNames, i)
1035  {
1036  label globalRegionI = surfaces.globalRegion(surfI, i);
1037 
1038  // Add master side patch
1039  {
1040  label patchi;
1041 
1042  if (surfacePatchInfo.set(globalRegionI))
1043  {
1044  patchi = meshRefiner.addMeshedPatch
1045  (
1046  regNames[i],
1047  surfacePatchInfo[globalRegionI]
1048  );
1049  }
1050  else
1051  {
1052  dictionary patchInfo;
1053  patchInfo.set("type", wallPolyPatch::typeName);
1054 
1055  patchi = meshRefiner.addMeshedPatch
1056  (
1057  regNames[i],
1058  patchInfo
1059  );
1060  }
1061 
1062  Info<< setf(ios_base::left)
1063  << setw(6) << patchi
1064  << setw(20) << mesh.boundaryMesh()[patchi].type()
1065  << setw(30) << regNames[i] << nl;
1066 
1067  globalToMasterPatch[globalRegionI] = patchi;
1068  }
1069  // Add slave side patch
1070  {
1071  const word slaveName = regNames[i] + "_slave";
1072  label patchi;
1073 
1074  if (surfacePatchInfo.set(globalRegionI))
1075  {
1076  patchi = meshRefiner.addMeshedPatch
1077  (
1078  slaveName,
1079  surfacePatchInfo[globalRegionI]
1080  );
1081  }
1082  else
1083  {
1084  dictionary patchInfo;
1085  patchInfo.set("type", wallPolyPatch::typeName);
1086 
1087  patchi = meshRefiner.addMeshedPatch
1088  (
1089  slaveName,
1090  patchInfo
1091  );
1092  }
1093 
1094  Info<< setf(ios_base::left)
1095  << setw(6) << patchi
1096  << setw(20) << mesh.boundaryMesh()[patchi].type()
1097  << setw(30) << slaveName << nl;
1098 
1099  globalToSlavePatch[globalRegionI] = patchi;
1100  }
1101  }
1102  }
1103 
1104  Info<< nl;
1105  }
1106 
1107  meshRefiner.addedMeshedPatches();
1108 
1109  Info<< "Added patches in = "
1110  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
1111  }
1112 
1113 
1114  // Parallel
1115  // ~~~~~~~~
1116 
1117  // Decomposition
1118  autoPtr<decompositionMethod> decomposerPtr
1119  (
1121  );
1122  decompositionMethod& decomposer = decomposerPtr();
1123 
1124  // Mesh distribution engine (uses tolerance to reconstruct meshes)
1125  fvMeshDistribute distributor(mesh);
1126 
1127 
1128  // Now do the real work -refinement -snapping -layers
1129  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1130 
1131  const Switch wantRefine(meshDict.lookup("castellatedMesh"));
1132  const Switch wantSnap(meshDict.lookup("snap"));
1133  const Switch wantLayers(meshDict.lookup("addLayers"));
1134 
1135  // Refinement parameters
1136  const refinementParameters refineParams(refineDict);
1137 
1138  // Snap parameters
1139  const snapParameters snapParams(snapDict);
1140 
1141  if (wantRefine)
1142  {
1143  cpuTime timer;
1144 
1145  snappyRefineDriver refineDriver
1146  (
1147  meshRefiner,
1148  decomposer,
1149  distributor,
1150  globalToMasterPatch,
1151  globalToSlavePatch
1152  );
1153 
1154  if (!overwrite && !debugLevel)
1155  {
1156  const_cast<Time&>(mesh.time())++;
1157  }
1158 
1159  refineDriver.doRefine
1160  (
1161  refineDict,
1162  refineParams,
1163  snapParams,
1164  refineParams.handleSnapProblems(),
1165  motionDict
1166  );
1167 
1168  if (!keepPatches && !wantSnap && !wantLayers)
1169  {
1170  removeZeroSizedPatches(mesh);
1171  }
1172 
1173  writeMesh
1174  (
1175  "Refined mesh",
1176  meshRefiner,
1177  debugLevel,
1179  );
1180 
1181  Info<< "Mesh refined in = "
1182  << timer.cpuTimeIncrement() << " s." << endl;
1183  }
1184 
1185  if (wantSnap)
1186  {
1187  cpuTime timer;
1188 
1189  snappySnapDriver snapDriver
1190  (
1191  meshRefiner,
1192  globalToMasterPatch,
1193  globalToSlavePatch
1194  );
1195 
1196  if (!overwrite && !debugLevel)
1197  {
1198  const_cast<Time&>(mesh.time())++;
1199  }
1200 
1201  // Use the resolveFeatureAngle from the refinement parameters
1202  scalar curvature = refineParams.curvature();
1203  scalar planarAngle = refineParams.planarAngle();
1204 
1205  snapDriver.doSnap
1206  (
1207  snapDict,
1208  motionDict,
1209  curvature,
1210  planarAngle,
1211  snapParams
1212  );
1213 
1214  if (!keepPatches && !wantLayers)
1215  {
1216  removeZeroSizedPatches(mesh);
1217  }
1218 
1219  writeMesh
1220  (
1221  "Snapped mesh",
1222  meshRefiner,
1223  debugLevel,
1225  );
1226 
1227  Info<< "Mesh snapped in = "
1228  << timer.cpuTimeIncrement() << " s." << endl;
1229  }
1230 
1231  if (wantLayers)
1232  {
1233  cpuTime timer;
1234 
1235  // Layer addition parameters dictionary
1236  const dictionary& layersDict = meshDict.subDict("addLayersControls");
1237 
1238  // Layer addition parameters
1239  const layerParameters layerParams(layersDict, mesh.boundaryMesh());
1240 
1241  snappyLayerDriver layerDriver
1242  (
1243  meshRefiner,
1244  globalToMasterPatch,
1245  globalToSlavePatch
1246  );
1247 
1248  // Use the maxLocalCells from the refinement parameters
1249  bool preBalance = returnReduce
1250  (
1251  (mesh.nCells() >= refineParams.maxLocalCells()),
1252  orOp<bool>()
1253  );
1254 
1255 
1256  if (!overwrite && !debugLevel)
1257  {
1258  const_cast<Time&>(mesh.time())++;
1259  }
1260 
1261  layerDriver.doLayers
1262  (
1263  layersDict,
1264  motionDict,
1265  layerParams,
1266  preBalance,
1267  decomposer,
1268  distributor
1269  );
1270 
1271  if (!keepPatches)
1272  {
1273  removeZeroSizedPatches(mesh);
1274  }
1275 
1276  writeMesh
1277  (
1278  "Layer mesh",
1279  meshRefiner,
1280  debugLevel,
1282  );
1283 
1284  Info<< "Layers added in = "
1285  << timer.cpuTimeIncrement() << " s." << endl;
1286  }
1287 
1288 
1289  {
1290  // Check final mesh
1291  Info<< "Checking final mesh ..." << endl;
1292  faceSet wrongFaces(mesh, "wrongFaces", mesh.nFaces()/100);
1293  meshCheck::checkMesh(false, mesh, motionDict, wrongFaces);
1294  const label nErrors = returnReduce
1295  (
1296  wrongFaces.size(),
1297  sumOp<label>()
1298  );
1299 
1300  if (nErrors > 0)
1301  {
1302  Info<< "Finished meshing with " << nErrors << " illegal faces"
1303  << " (concave, zero area or negative cell pyramid volume)"
1304  << endl;
1305  wrongFaces.write();
1306  }
1307  else
1308  {
1309  Info<< "Finished meshing without any errors" << endl;
1310  }
1311  }
1312 
1313 
1314  if (surfaceSimplify)
1315  {
1317 
1319 
1320  if (args.optionFound("patches"))
1321  {
1323  (
1324  wordReList(args.optionLookup("patches")())
1325  );
1326  }
1327  else
1328  {
1329  forAll(bMesh, patchi)
1330  {
1331  const polyPatch& patch = bMesh[patchi];
1332 
1333  if (!isA<processorPolyPatch>(patch))
1334  {
1336  }
1337  }
1338  }
1339 
1340  fileName outFileName
1341  (
1343  (
1344  "outFile",
1345  "constant/triSurface/simplifiedSurface.stl"
1346  )
1347  );
1348 
1349  extractSurface
1350  (
1351  mesh,
1352  runTime,
1354  outFileName
1355  );
1356 
1357  pointIOField cellCentres
1358  (
1359  IOobject
1360  (
1361  "internalCellCentres",
1362  runTime.name(),
1363  mesh,
1366  ),
1367  mesh.cellCentres()
1368  );
1369 
1370  cellCentres.write();
1371  }
1372 
1373 
1374  Info<< "Finished meshing in = "
1375  << runTime.elapsedCpuTime() << " s." << endl;
1376 
1377  Info<< "End\n" << endl;
1378 
1379  return 0;
1380 }
1381 
1382 
1383 // ************************************************************************* //
Istream and Ostream manipulators taking arguments.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
A HashTable with keys but without contents.
Definition: HashSet.H:62
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
bool erase(T *p)
Remove the specified element from the list and delete it.
Definition: ILList.C:105
A primitive field of type <Type> with automated input and output.
Definition: IOField.H:53
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:307
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:473
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
virtual Ostream & write(const char)=0
Write character.
A list of faces which address into the list of points.
static void mapCombineScatter(const List< commsStruct > &comms, Container &Values, const int tag, const label comm)
Scatter data. Reverse of combineGather.
static void mapCombineGather(const List< commsStruct > &comms, Container &Values, const CombineOp &cop, const int tag, const label comm)
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
static const word & constant()
Return constant name.
Definition: TimePaths.H:122
fileName globalPath() const
Return the global path.
Definition: TimePaths.H:132
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
fileName path() const
Explicitly inherit path from TimePaths to disambiguate from.
Definition: TimePaths.H:138
IOstream::streamFormat writeFormat() const
Default write format.
Definition: Time.H:281
A List with indirect addressing.
Definition: UIndirectList.H:60
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
A surface geometry mesh, in which the surface zone information is conveyed by the 'zoneId' associated...
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:128
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
const word & executable() const
Name of executable without the path.
Definition: argListI.H:36
IStringStream optionLookup(const word &opt) const
Return an IStringStream from the named option.
Definition: argListI.H:120
T optionLookupOrDefault(const word &opt, const T &deflt) const
Read a value from the named option if present.
Definition: argListI.H:294
static void removeOption(const word &opt)
Remove option from validOptions and from optionUsage.
Definition: argList.C:168
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
void set(T *)
Set pointer to that given.
Definition: autoPtrI.H:99
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
Starts timing CPU usage and return elapsed time from start.
Definition: cpuTime.H:55
double cpuTimeIncrement() const
Return CPU time (in seconds) since last call to cpuTimeIncrement()
Definition: cpuTime.C:60
double elapsedCpuTime() const
Return CPU time (in seconds) from the start.
Definition: cpuTime.C:53
Abstract base class for decomposition.
static IOdictionary decomposeParDict(const Time &time)
Read and return the decomposeParDict.
static autoPtr< decompositionMethod > NewDistributor(const dictionary &decompositionDict)
Return a reference to the selected decomposition method.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
autoPtr< dictionary > clone() const
Construct and return clone.
Definition: dictionary.C:421
const entry * lookupEntryPtr(const word &, bool recursive, bool patternMatch) const
Find and return an entry data stream pointer if present.
Definition: dictionary.C:578
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:927
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:849
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1020
T lookupOrDefault(const word &, const T &, const bool writeDefault=writeOptionalEntries > 0) const
Find and return a T, if not found return the given default.
wordList toc() const
Return the table of contents.
Definition: dictionary.C:982
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:539
static const dictionary null
Null dictionary.
Definition: dictionary.H:273
const word & name() const
Return const reference to name.
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:68
const keyType & keyword() const
Return keyword.
Definition: entry.H:136
virtual const dictionary & dict() const =0
Return dictionary if this entry is a dictionary.
A list of face labels.
Definition: faceSet.H:51
A class for handling file names.
Definition: fileName.H:82
Sends/receives parts of mesh+fvfields to neighbouring processors. Used in load balancing.
static void reorderPatches(fvMesh &, const labelList &oldToNew, const label nPatches, const bool validBoundary)
Reorder and remove trailing patches. If validBoundary call is.
Definition: fvMeshTools.C:104
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
autoPtr< globalIndex > mergePoints(labelList &pointToGlobal, labelList &uniquePoints) const
Helper for merging (collocated!) mesh point data.
Simple container to keep together layer specific information.
Helper class which maintains intersections of (changing) mesh with (static) surfaces.
static int readFlags(const Enum &namedEnum, const wordList &)
Helper: convert wordList into bit pattern using provided.
word name() const
Replacement for Time::name() : return oldInstance (if.
const fvMesh & mesh() const
Reference to mesh.
static const NamedEnum< IOwriteType, 5 > IOwriteTypeNames
void printMeshInfo(const bool, const string &) const
Print some mesh stats.
static const NamedEnum< IOdebugType, 5 > IOdebugTypeNames
static void checkCoupledFaceZones(const polyMesh &)
Helper function: check that face zones are synced.
label addMeshedPatch(const word &name, const dictionary &)
Add patch originating from meshing. Update meshedPatches_.
static const NamedEnum< IOoutputType, 1 > IOoutputTypeNames
bool write() const
Write mesh and all data.
static writeType writeLevel()
Get/set write level.
static outputType outputLevel()
Get/set output level.
void addedMeshedPatches()
Complete adding patches originating from meshing.
label index() const
Return the index of this patch in the boundaryMesh.
const word & name() const
Return name.
Foam::polyBoundaryMesh.
label findIndex(const word &patchName) const
Find patch index given a name.
wordList types() const
Return a list of patch types.
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the patch set corresponding to the given names.
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order. Return.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1344
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1521
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:1040
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:407
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: polyPatch.C:238
const vectorField & cellCentres() const
label nCells() const
label nFaces() const
Encapsulates queries for features.
Simple container to keep together refinement specific information.
Encapsulates queries for volume refinement ('refine all cells within shell').
Container for data on surfaces used for surface-driven refinement. Contains all the data about the le...
const wordList & names() const
Names of surfaces.
const labelList & minLevel() const
From global region number to refinement level.
const labelList & maxLevel() const
From global region number to refinement level.
const labelList & surfaces() const
label globalRegion(const label surfi, const label regioni) const
From surface and region on surface to global region.
const labelList & gapLevel() const
From global region number to small gap refinement level.
Container for searchableSurfaces.
label checkGeometry(const scalar maxRatio, const scalar tolerance, const scalar minQuality, const bool report) const
All geometric checks. Return number of failed checks.
const wordList & names() const
const List< wordList > & regionNames() const
label checkTopology(const bool report) const
All topological checks. Return number of failed checks.
void writeStats(const List< wordList > &, Ostream &) const
Write some stats.
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
static const word & geometryDir()
Return the geometry directory name.
Simple container to keep together snap specific information.
All to do with adding layers.
All to do with snapping to surface.
An identifier for a surface zone on a meshed surface.
Implements a timeout mechanism via sigalarm.
Definition: timer.H:82
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:329
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
int main(int argc, char *argv[])
Definition: financialFoam.C:44
const polyBoundaryMesh & bMesh
label patchi
Functions for checking mesh topology and geometry.
bool checkMesh(const bool report, const polyMesh &mesh, const dictionary &dict, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet &wrongFaces)
Check (subset of mesh including baffles) with mesh settings in dict.
Definition: checkMesh.C:33
label checkGeometry(const polyMesh &mesh, const bool allGeometry, const scalar nonOrthThreshold, const scalar skewThreshold, const autoPtr< surfaceWriter > &, const autoPtr< setWriter > &)
Check the geometry.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Smanip< ios_base::fmtflags > setf(const ios_base::fmtflags flags)
Definition: IOmanip.H:164
const doubleScalar e
Definition: doubleScalar.H:106
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:242
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:235
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
T clone(const T &t)
Definition: List.H:55
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
error FatalError
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion, const fileName &path=fileName::null)
Definition: systemDict.C:93
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:228
static const char nl
Definition: Ostream.H:267
wordList patchTypes(nPatches)
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(1e-3))
labelHashSet includePatches
const bool overwrite
Definition: setNoOverwrite.H:1
const Foam::wordList regionNames(args.optionFound("allRegions") ? runTime .controlDict().subDict("regionSolvers").toc() :wordList(1, args.optionFound("region") ? args.optionRead< word >("region") :polyMesh::defaultRegion))
Foam::argList args(argc, argv)