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-2024 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 "searchableSurfaces.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 searchableSurfaces& 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 {
377  const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
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;
437  autoPtr<globalIndex> globalNumbers = mesh.globalData().mergePoints
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
534  if (mesh.time().writeFormat() == IOstream::ASCII)
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 "addOverwriteOption.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 "addRegionOption.H"
662 
663  #include "setRootCase.H"
666 
667  Info<< "Read mesh in = "
668  << runTime.cpuTimeIncrement() << " s" << endl;
669 
670  const bool overwrite = args.optionFound("overwrite");
671  const bool checkGeometry = args.optionFound("checkGeometry");
672  const bool surfaceSimplify = args.optionFound("surfaceSimplify");
673 
674  // Check that the read mesh is fully 3D
675  // as required for mesh relaxation after snapping
676  if (mesh.nSolutionD() != 3)
677  {
679  << "Mesh provided is not fully 3D"
680  " as required for mesh relaxation after snapping" << nl
681  << "Convert all empty patches to appropriate types for a 3D mesh,"
682  " current patch types are" << nl << mesh.boundaryMesh().types()
683  << exit(FatalError);
684  }
685 
686  // Check patches and faceZones are synchronised
687  mesh.boundaryMesh().checkParallelSync(true);
689 
690 
691  // Read meshing dictionary
692  const IOdictionary meshDict(systemDict("snappyHexMeshDict", args, mesh));
693 
694 
695  // all surface geometry
696  const dictionary& geometryDict = meshDict.subDict("geometry");
697 
698  // refinement parameters
699  const dictionary& refineDict = meshDict.subDict("castellatedMeshControls");
700 
701  // mesh motion and mesh quality parameters
702  const dictionary& motionDict = meshDict.subDict("meshQualityControls");
703 
704  // snap-to-surface parameters
705  const dictionary& snapDict = meshDict.subDict("snapControls");
706 
707  // absolute merge distance
708  const scalar mergeDist = getMergeDistance
709  (
710  mesh,
711  meshDict.lookup<scalar>("mergeTolerance")
712  );
713 
714  const Switch keepPatches(meshDict.lookupOrDefault("keepPatches", false));
715 
716 
717  // Read decomposePar dictionary
718  dictionary decomposeDict;
719  {
720  if (Pstream::parRun())
721  {
722  decomposeDict = decompositionMethod::decomposeParDict(runTime);
723  }
724  else
725  {
726  decomposeDict.add("method", "none");
727  decomposeDict.add("numberOfSubdomains", 1);
728  }
729  }
730 
731 
732  // Debug
733  // ~~~~~
734 
735  // Set debug level
737  (
738  meshDict.lookupOrDefault<label>
739  (
740  "debug",
741  0
742  )
743  );
744  {
745  wordList flags;
746  if (meshDict.readIfPresent("debugFlags", flags))
747  {
748  debugLevel = meshRefinement::debugType
749  (
751  (
753  flags
754  )
755  );
756  }
757  }
758  if (debugLevel > 0)
759  {
760  meshRefinement::debug = debugLevel;
761  snappyRefineDriver::debug = debugLevel;
762  snappySnapDriver::debug = debugLevel;
763  snappyLayerDriver::debug = debugLevel;
764  }
765 
766  // Set file writing level
767  {
768  wordList flags;
769  if (meshDict.readIfPresent("writeFlags", flags))
770  {
772  (
774  (
776  (
778  flags
779  )
780  )
781  );
782  }
783  }
784 
785  // Set output level
786  {
787  wordList flags;
788  if (meshDict.readIfPresent("outputFlags", flags))
789  {
791  (
793  (
795  (
797  flags
798  )
799  )
800  );
801  }
802  }
803 
804 
805  // Read geometry
806  // ~~~~~~~~~~~~~
807 
808  searchableSurfaces allGeometry
809  (
810  IOobject
811  (
812  "abc",
813  mesh.time().constant(),
815  mesh.time(),
818  ),
819  geometryDict,
820  meshDict.lookupOrDefault("singleRegionName", true)
821  );
822 
823 
824  // Read refinement surfaces
825  // ~~~~~~~~~~~~~~~~~~~~~~~~
826 
827  Info<< "Reading refinement surfaces..." << endl;
828 
829  refinementSurfaces surfaces
830  (
831  allGeometry,
832  refineDict.found("refinementSurfaces")
833  ? refineDict.subDict("refinementSurfaces")
835  refineDict.lookupOrDefault("gapLevelIncrement", 0)
836  );
837 
838  Info<< "Read refinement surfaces in = "
839  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
840 
841 
842  // Checking only?
843 
844  if (checkGeometry)
845  {
846  // Extract patchInfo
847  List<wordList> patchTypes(allGeometry.size());
848 
849  const PtrList<dictionary>& patchInfo = surfaces.patchInfo();
850  const labelList& surfaceGeometry = surfaces.surfaces();
851  forAll(surfaceGeometry, surfI)
852  {
853  label geomI = surfaceGeometry[surfI];
854  const wordList& regNames = allGeometry.regionNames()[geomI];
855 
856  patchTypes[geomI].setSize(regNames.size());
857  forAll(regNames, regionI)
858  {
859  label globalRegionI = surfaces.globalRegion(surfI, regionI);
860 
861  if (patchInfo.set(globalRegionI))
862  {
863  patchTypes[geomI][regionI] =
864  word(patchInfo[globalRegionI].lookup("type"));
865  }
866  else
867  {
868  patchTypes[geomI][regionI] = wallPolyPatch::typeName;
869  }
870  }
871  }
872 
873  // Write some stats
874  allGeometry.writeStats(patchTypes, Info);
875  // Check topology
876  allGeometry.checkTopology(true);
877  // Check geometry
878  allGeometry.checkGeometry
879  (
880  100.0, // max size ratio
881  1e-9, // intersection tolerance
882  0.01, // min triangle quality
883  true
884  );
885 
886  return 0;
887  }
888 
889 
890 
891  // Read refinement shells
892  // ~~~~~~~~~~~~~~~~~~~~~~
893 
894  Info<< "Reading refinement regions..." << endl;
895  refinementRegions shells
896  (
897  allGeometry,
898  refineDict.found("refinementRegions")
899  ? refineDict.subDict("refinementRegions")
901  );
902  Info<< "Read refinement regions in = "
903  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
904 
905 
906 
907  // Read feature meshes
908  // ~~~~~~~~~~~~~~~~~~~
909 
910  Info<< "Reading features..." << endl;
911  refinementFeatures features
912  (
913  mesh,
914  refineDict.found("features")
915  ? refineDict.lookup("features")
917  );
918  Info<< "Read features in = "
919  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
920 
921 
922 
923  // Refinement engine
924  // ~~~~~~~~~~~~~~~~~
925 
926  Info<< nl
927  << "Determining initial surface intersections" << nl
928  << "-----------------------------------------" << nl
929  << endl;
930 
931  // Main refinement engine
932  meshRefinement meshRefiner
933  (
934  mesh,
935  refineDict,
936  mergeDist, // tolerance used in sorting coordinates
937  overwrite, // overwrite mesh files?
938  surfaces, // for surface intersection refinement
939  features, // for feature edges/point based refinement
940  shells // for volume (inside/outside) refinement
941  );
942  Info<< "Calculated surface intersections in = "
943  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
944 
945  // Some stats
946  meshRefiner.printMeshInfo(debugLevel, "Initial mesh");
947 
948  meshRefiner.write
949  (
952  mesh.time().path()/meshRefiner.name()
953  );
954 
955 
956  // Add all the surface regions as patches
957  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
958 
959  //- Global surface region to patch (non faceZone surface) or patches
960  // (faceZone surfaces)
961  labelList globalToMasterPatch;
962  labelList globalToSlavePatch;
963  {
964  Info<< nl
965  << "Adding patches for surface regions" << nl
966  << "----------------------------------" << nl
967  << endl;
968 
969  // From global region number to mesh patch.
970  globalToMasterPatch.setSize(surfaces.nRegions(), -1);
971  globalToSlavePatch.setSize(surfaces.nRegions(), -1);
972 
973  Info<< setf(ios_base::left)
974  << setw(6) << "Patch"
975  << setw(20) << "Type"
976  << setw(30) << "Region" << nl
977  << setw(6) << "-----"
978  << setw(20) << "----"
979  << setw(30) << "------" << endl;
980 
981  const labelList& surfaceGeometry = surfaces.surfaces();
982  const PtrList<dictionary>& surfacePatchInfo = surfaces.patchInfo();
983 
984  forAll(surfaceGeometry, surfI)
985  {
986  label geomI = surfaceGeometry[surfI];
987 
988  const wordList& regNames = allGeometry.regionNames()[geomI];
989 
990  Info<< surfaces.names()[surfI] << ':' << nl << nl;
991 
992  if (surfaces.surfZones()[surfI].faceZoneName().empty())
993  {
994  // 'Normal' surface
995  forAll(regNames, i)
996  {
997  label globalRegionI = surfaces.globalRegion(surfI, i);
998 
999  label patchi;
1000 
1001  if (surfacePatchInfo.set(globalRegionI))
1002  {
1003  patchi = meshRefiner.addMeshedPatch
1004  (
1005  regNames[i],
1006  surfacePatchInfo[globalRegionI]
1007  );
1008  }
1009  else
1010  {
1011  dictionary patchInfo;
1012  patchInfo.set("type", wallPolyPatch::typeName);
1013 
1014  patchi = meshRefiner.addMeshedPatch
1015  (
1016  regNames[i],
1017  patchInfo
1018  );
1019  }
1020 
1021  Info<< setf(ios_base::left)
1022  << setw(6) << patchi
1023  << setw(20) << mesh.boundaryMesh()[patchi].type()
1024  << setw(30) << regNames[i] << nl;
1025 
1026  globalToMasterPatch[globalRegionI] = patchi;
1027  globalToSlavePatch[globalRegionI] = patchi;
1028  }
1029  }
1030  else
1031  {
1032  // Zoned surface
1033  forAll(regNames, i)
1034  {
1035  label globalRegionI = surfaces.globalRegion(surfI, i);
1036 
1037  // Add master side patch
1038  {
1039  label patchi;
1040 
1041  if (surfacePatchInfo.set(globalRegionI))
1042  {
1043  patchi = meshRefiner.addMeshedPatch
1044  (
1045  regNames[i],
1046  surfacePatchInfo[globalRegionI]
1047  );
1048  }
1049  else
1050  {
1051  dictionary patchInfo;
1052  patchInfo.set("type", wallPolyPatch::typeName);
1053 
1054  patchi = meshRefiner.addMeshedPatch
1055  (
1056  regNames[i],
1057  patchInfo
1058  );
1059  }
1060 
1061  Info<< setf(ios_base::left)
1062  << setw(6) << patchi
1063  << setw(20) << mesh.boundaryMesh()[patchi].type()
1064  << setw(30) << regNames[i] << nl;
1065 
1066  globalToMasterPatch[globalRegionI] = patchi;
1067  }
1068  // Add slave side patch
1069  {
1070  const word slaveName = regNames[i] + "_slave";
1071  label patchi;
1072 
1073  if (surfacePatchInfo.set(globalRegionI))
1074  {
1075  patchi = meshRefiner.addMeshedPatch
1076  (
1077  slaveName,
1078  surfacePatchInfo[globalRegionI]
1079  );
1080  }
1081  else
1082  {
1083  dictionary patchInfo;
1084  patchInfo.set("type", wallPolyPatch::typeName);
1085 
1086  patchi = meshRefiner.addMeshedPatch
1087  (
1088  slaveName,
1089  patchInfo
1090  );
1091  }
1092 
1093  Info<< setf(ios_base::left)
1094  << setw(6) << patchi
1095  << setw(20) << mesh.boundaryMesh()[patchi].type()
1096  << setw(30) << slaveName << nl;
1097 
1098  globalToSlavePatch[globalRegionI] = patchi;
1099  }
1100  }
1101  }
1102 
1103  Info<< nl;
1104  }
1105  Info<< "Added patches in = "
1106  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
1107  }
1108 
1109 
1110  // Parallel
1111  // ~~~~~~~~
1112 
1113  // Decomposition
1114  autoPtr<decompositionMethod> decomposerPtr
1115  (
1117  );
1118  decompositionMethod& decomposer = decomposerPtr();
1119 
1120  // Mesh distribution engine (uses tolerance to reconstruct meshes)
1121  fvMeshDistribute distributor(mesh);
1122 
1123 
1124  // Now do the real work -refinement -snapping -layers
1125  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1126 
1127  const Switch wantRefine(meshDict.lookup("castellatedMesh"));
1128  const Switch wantSnap(meshDict.lookup("snap"));
1129  const Switch wantLayers(meshDict.lookup("addLayers"));
1130 
1131  // Refinement parameters
1132  const refinementParameters refineParams(refineDict);
1133 
1134  // Snap parameters
1135  const snapParameters snapParams(snapDict);
1136 
1137  if (wantRefine)
1138  {
1139  cpuTime timer;
1140 
1141  snappyRefineDriver refineDriver
1142  (
1143  meshRefiner,
1144  decomposer,
1145  distributor,
1146  globalToMasterPatch,
1147  globalToSlavePatch
1148  );
1149 
1150  if (!overwrite && !debugLevel)
1151  {
1152  const_cast<Time&>(mesh.time())++;
1153  }
1154 
1155  refineDriver.doRefine
1156  (
1157  refineDict,
1158  refineParams,
1159  snapParams,
1160  refineParams.handleSnapProblems(),
1161  motionDict
1162  );
1163 
1164  if (!keepPatches && !wantSnap && !wantLayers)
1165  {
1166  removeZeroSizedPatches(mesh);
1167  }
1168 
1169  writeMesh
1170  (
1171  "Refined mesh",
1172  meshRefiner,
1173  debugLevel,
1175  );
1176 
1177  Info<< "Mesh refined in = "
1178  << timer.cpuTimeIncrement() << " s." << endl;
1179  }
1180 
1181  if (wantSnap)
1182  {
1183  cpuTime timer;
1184 
1185  snappySnapDriver snapDriver
1186  (
1187  meshRefiner,
1188  globalToMasterPatch,
1189  globalToSlavePatch
1190  );
1191 
1192  if (!overwrite && !debugLevel)
1193  {
1194  const_cast<Time&>(mesh.time())++;
1195  }
1196 
1197  // Use the resolveFeatureAngle from the refinement parameters
1198  scalar curvature = refineParams.curvature();
1199  scalar planarAngle = refineParams.planarAngle();
1200 
1201  snapDriver.doSnap
1202  (
1203  snapDict,
1204  motionDict,
1205  curvature,
1206  planarAngle,
1207  snapParams
1208  );
1209 
1210  if (!keepPatches && !wantLayers)
1211  {
1212  removeZeroSizedPatches(mesh);
1213  }
1214 
1215  writeMesh
1216  (
1217  "Snapped mesh",
1218  meshRefiner,
1219  debugLevel,
1221  );
1222 
1223  Info<< "Mesh snapped in = "
1224  << timer.cpuTimeIncrement() << " s." << endl;
1225  }
1226 
1227  if (wantLayers)
1228  {
1229  cpuTime timer;
1230 
1231  // Layer addition parameters dictionary
1232  const dictionary& layersDict = meshDict.subDict("addLayersControls");
1233 
1234  // Layer addition parameters
1235  const layerParameters layerParams(layersDict, mesh.boundaryMesh());
1236 
1237  snappyLayerDriver layerDriver
1238  (
1239  meshRefiner,
1240  globalToMasterPatch,
1241  globalToSlavePatch
1242  );
1243 
1244  // Use the maxLocalCells from the refinement parameters
1245  bool preBalance = returnReduce
1246  (
1247  (mesh.nCells() >= refineParams.maxLocalCells()),
1248  orOp<bool>()
1249  );
1250 
1251 
1252  if (!overwrite && !debugLevel)
1253  {
1254  const_cast<Time&>(mesh.time())++;
1255  }
1256 
1257  layerDriver.doLayers
1258  (
1259  layersDict,
1260  motionDict,
1261  layerParams,
1262  preBalance,
1263  decomposer,
1264  distributor
1265  );
1266 
1267  if (!keepPatches)
1268  {
1269  removeZeroSizedPatches(mesh);
1270  }
1271 
1272  writeMesh
1273  (
1274  "Layer mesh",
1275  meshRefiner,
1276  debugLevel,
1278  );
1279 
1280  Info<< "Layers added in = "
1281  << timer.cpuTimeIncrement() << " s." << endl;
1282  }
1283 
1284 
1285  {
1286  // Check final mesh
1287  Info<< "Checking final mesh ..." << endl;
1288  faceSet wrongFaces(mesh, "wrongFaces", mesh.nFaces()/100);
1289  meshCheck::checkMesh(false, mesh, motionDict, wrongFaces);
1290  const label nErrors = returnReduce
1291  (
1292  wrongFaces.size(),
1293  sumOp<label>()
1294  );
1295 
1296  if (nErrors > 0)
1297  {
1298  Info<< "Finished meshing with " << nErrors << " illegal faces"
1299  << " (concave, zero area or negative cell pyramid volume)"
1300  << endl;
1301  wrongFaces.write();
1302  }
1303  else
1304  {
1305  Info<< "Finished meshing without any errors" << endl;
1306  }
1307  }
1308 
1309 
1310  if (surfaceSimplify)
1311  {
1312  const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
1313 
1315 
1316  if (args.optionFound("patches"))
1317  {
1319  (
1320  wordReList(args.optionLookup("patches")())
1321  );
1322  }
1323  else
1324  {
1325  forAll(bMesh, patchi)
1326  {
1327  const polyPatch& patch = bMesh[patchi];
1328 
1329  if (!isA<processorPolyPatch>(patch))
1330  {
1332  }
1333  }
1334  }
1335 
1336  fileName outFileName
1337  (
1339  (
1340  "outFile",
1341  "constant/triSurface/simplifiedSurface.stl"
1342  )
1343  );
1344 
1345  extractSurface
1346  (
1347  mesh,
1348  runTime,
1350  outFileName
1351  );
1352 
1353  pointIOField cellCentres
1354  (
1355  IOobject
1356  (
1357  "internalCellCentres",
1358  runTime.name(),
1359  mesh,
1362  ),
1363  mesh.cellCentres()
1364  );
1365 
1366  cellCentres.write();
1367  }
1368 
1369 
1370  Info<< "Finished meshing in = "
1371  << runTime.elapsedCpuTime() << " s." << endl;
1372 
1373  Info<< "End\n" << endl;
1374 
1375  return 0;
1376 }
1377 
1378 
1379 // ************************************************************************* //
Istream and Ostream manipulators taking arguments.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
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:310
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:458
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:287
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:243
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 keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
autoPtr< dictionary > clone() const
Construct and return clone.
Definition: dictionary.C:398
const entry * lookupEntryPtr(const word &, bool recursive, bool patternMatch) const
Find and return an entry data stream pointer if present.
Definition: dictionary.C:548
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:1152
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:921
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:843
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1014
wordList toc() const
Return the table of contents.
Definition: dictionary.C:976
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
static const dictionary null
Null dictionary.
Definition: dictionary.H:258
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:123
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 parallel.
Definition: fvMeshTools.C:141
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
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.
const Time & time() const
Return time.
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 set of patch indices 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:1326
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1515
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:1017
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1313
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:410
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
label nCells() const
const vectorField & cellCentres() 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.
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
static const word & geometryDir()
Return the geometry directory name.
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.
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
#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:241
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
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:92
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
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:234
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:227
static const char nl
Definition: Ostream.H:266
wordList patchTypes(nPatches)
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(1e-3))
labelHashSet includePatches
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)