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-2023 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. Refines and snaps to surface.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "argList.H"
33 #include "Time.H"
34 #include "fvMesh.H"
35 #include "snappyRefineDriver.H"
36 #include "snappySnapDriver.H"
37 #include "snappyLayerDriver.H"
38 #include "searchableSurfaces.H"
39 #include "refinementSurfaces.H"
40 #include "refinementFeatures.H"
41 #include "refinementRegions.H"
42 #include "decompositionMethod.H"
43 #include "noDecomp.H"
44 #include "fvMeshDistribute.H"
45 #include "wallPolyPatch.H"
46 #include "refinementParameters.H"
47 #include "snapParameters.H"
48 #include "layerParameters.H"
49 #include "faceSet.H"
50 #include "motionSmoother.H"
51 #include "polyTopoChange.H"
52 #include "cellModeller.H"
54 #include "surfZoneIdentifierList.H"
55 #include "UnsortedMeshedSurface.H"
56 #include "MeshedSurface.H"
57 #include "globalIndex.H"
58 #include "IOmanip.H"
59 #include "fvMeshTools.H"
60 #include "systemDict.H"
61 
62 using namespace Foam;
63 
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66 // Convert size (as fraction of defaultCellSize) to refinement level
67 label sizeCoeffToRefinement
68 (
69  const scalar level0Coeff, // ratio of hex cell size v.s. defaultCellSize
70  const scalar sizeCoeff
71 )
72 {
73  return round(::log(level0Coeff/sizeCoeff)/::log(2));
74 }
75 
76 
77 autoPtr<refinementSurfaces> createRefinementSurfaces
78 (
79  const searchableSurfaces& allGeometry,
80  const dictionary& surfacesDict,
81  const dictionary& shapeControlDict,
82  const label gapLevelIncrement,
83  const scalar level0Coeff
84 )
85 {
86  autoPtr<refinementSurfaces> surfacePtr;
87 
88  // Count number of surfaces.
89  label surfI = 0;
90  forAll(allGeometry.names(), geomI)
91  {
92  const word& geomName = allGeometry.names()[geomI];
93 
94  if (surfacesDict.found(geomName))
95  {
96  surfI++;
97  }
98  }
99 
100  labelList surfaces(surfI);
101  wordList names(surfI);
102  PtrList<surfaceZonesInfo> surfZones(surfI);
103 
104  labelList regionOffset(surfI);
105 
106  labelList globalMinLevel(surfI, 0);
107  labelList globalMaxLevel(surfI, 0);
108  labelList globalLevelIncr(surfI, 0);
109  PtrList<dictionary> globalPatchInfo(surfI);
110  List<Map<label>> regionMinLevel(surfI);
111  List<Map<label>> regionMaxLevel(surfI);
112  List<Map<label>> regionLevelIncr(surfI);
113  List<Map<scalar>> regionAngle(surfI);
114  List<Map<autoPtr<dictionary>>> regionPatchInfo(surfI);
115 
116  HashSet<word> unmatchedKeys(surfacesDict.toc());
117 
118  surfI = 0;
119  forAll(allGeometry.names(), geomI)
120  {
121  const word& geomName = allGeometry.names()[geomI];
122 
123  const entry* ePtr = surfacesDict.lookupEntryPtr(geomName, false, true);
124 
125  if (ePtr)
126  {
127  const dictionary& shapeDict = ePtr->dict();
128  unmatchedKeys.erase(ePtr->keyword());
129 
130  names[surfI] = geomName;
131  surfaces[surfI] = geomI;
132 
133  const searchableSurface& surface = allGeometry[geomI];
134 
135  // Find the index in shapeControlDict
136  // Invert surfaceCellSize to get the refinementLevel
137 
138  const word scsFuncName =
139  shapeDict.lookup("surfaceCellSizeFunction");
140  const dictionary& scsDict =
141  shapeDict.optionalSubDict(scsFuncName + "Coeffs");
142 
143  const scalar surfaceCellSize =
144  scsDict.lookup<scalar>("surfaceCellSizeCoeff");
145 
146  const label refLevel = sizeCoeffToRefinement
147  (
148  level0Coeff,
149  surfaceCellSize
150  );
151 
152  globalMinLevel[surfI] = refLevel;
153  globalMaxLevel[surfI] = refLevel;
154  globalLevelIncr[surfI] = gapLevelIncrement;
155 
156  // Surface zones
157  surfZones.set(surfI, new surfaceZonesInfo(surface, shapeDict));
158 
159 
160  // Global perpendicular angle
161  if (shapeDict.found("patchInfo"))
162  {
163  globalPatchInfo.set
164  (
165  surfI,
166  shapeDict.subDict("patchInfo").clone()
167  );
168  }
169 
170 
171  // Per region override of patchInfo
172 
173  if (shapeDict.found("regions"))
174  {
175  const dictionary& regionsDict = shapeDict.subDict("regions");
176  const wordList& regionNames =
177  allGeometry[surfaces[surfI]].regions();
178 
179  forAll(regionNames, regionI)
180  {
181  if (regionsDict.found(regionNames[regionI]))
182  {
183  // Get the dictionary for region
184  const dictionary& regionDict = regionsDict.subDict
185  (
186  regionNames[regionI]
187  );
188 
189  if (regionDict.found("patchInfo"))
190  {
191  regionPatchInfo[surfI].insert
192  (
193  regionI,
194  regionDict.subDict("patchInfo").clone()
195  );
196  }
197  }
198  }
199  }
200 
201  // Per region override of cellSize
202  if (shapeDict.found("regions"))
203  {
204  const dictionary& shapeControlRegionsDict =
205  shapeDict.subDict("regions");
206  const wordList& regionNames =
207  allGeometry[surfaces[surfI]].regions();
208 
209  forAll(regionNames, regionI)
210  {
211  if (shapeControlRegionsDict.found(regionNames[regionI]))
212  {
213  const dictionary& shapeControlRegionDict =
214  shapeControlRegionsDict.subDict
215  (
216  regionNames[regionI]
217  );
218 
219  const word scsFuncName =
220  shapeControlRegionDict.lookup
221  (
222  "surfaceCellSizeFunction"
223  );
224  const dictionary& scsDict =
225  shapeControlRegionDict.subDict
226  (
227  scsFuncName + "Coeffs"
228  );
229 
230  const scalar surfaceCellSize =
231  scsDict.lookup<scalar>("surfaceCellSizeCoeff");
232 
233  const label refLevel = sizeCoeffToRefinement
234  (
235  level0Coeff,
236  surfaceCellSize
237  );
238 
239  regionMinLevel[surfI].insert(regionI, refLevel);
240  regionMaxLevel[surfI].insert(regionI, refLevel);
241  regionLevelIncr[surfI].insert(regionI, 0);
242  }
243  }
244  }
245 
246  surfI++;
247  }
248  }
249 
250  // Calculate local to global region offset
251  label nRegions = 0;
252 
253  forAll(surfaces, surfI)
254  {
255  regionOffset[surfI] = nRegions;
256  nRegions += allGeometry[surfaces[surfI]].regions().size();
257  }
258 
259  // Rework surface specific information into information per global region
260  labelList minLevel(nRegions, 0);
261  labelList maxLevel(nRegions, 0);
262  labelList gapLevel(nRegions, -1);
263  PtrList<dictionary> patchInfo(nRegions);
264 
265  forAll(globalMinLevel, surfI)
266  {
267  label nRegions = allGeometry[surfaces[surfI]].regions().size();
268 
269  // Initialise to global (i.e. per surface)
270  for (label i = 0; i < nRegions; i++)
271  {
272  label globalRegionI = regionOffset[surfI] + i;
273  minLevel[globalRegionI] = globalMinLevel[surfI];
274  maxLevel[globalRegionI] = globalMaxLevel[surfI];
275  gapLevel[globalRegionI] =
276  maxLevel[globalRegionI]
277  + globalLevelIncr[surfI];
278 
279  if (globalPatchInfo.set(surfI))
280  {
281  patchInfo.set
282  (
283  globalRegionI,
284  globalPatchInfo[surfI].clone()
285  );
286  }
287  }
288 
289  // Overwrite with region specific information
290  forAllConstIter(Map<label>, regionMinLevel[surfI], iter)
291  {
292  label globalRegionI = regionOffset[surfI] + iter.key();
293 
294  minLevel[globalRegionI] = iter();
295  maxLevel[globalRegionI] = regionMaxLevel[surfI][iter.key()];
296  gapLevel[globalRegionI] =
297  maxLevel[globalRegionI]
298  + regionLevelIncr[surfI][iter.key()];
299  }
300 
301  const Map<autoPtr<dictionary>>& localInfo = regionPatchInfo[surfI];
302  forAllConstIter(Map<autoPtr<dictionary>>, localInfo, iter)
303  {
304  label globalRegionI = regionOffset[surfI] + iter.key();
305  patchInfo.set(globalRegionI, iter()().clone());
306  }
307  }
308 
309  surfacePtr.set
310  (
312  (
313  allGeometry,
314  surfaces,
315  names,
316  surfZones,
317  regionOffset,
318  minLevel,
319  maxLevel,
320  gapLevel,
321  scalarField(nRegions, -great), // perpendicularAngle,
322  patchInfo
323  )
324  );
325 
326 
327  const refinementSurfaces& rf = surfacePtr();
328 
329  // Determine maximum region name length
330  label maxLen = 0;
331  forAll(rf.surfaces(), surfI)
332  {
333  label geomI = rf.surfaces()[surfI];
334  const wordList& regionNames = allGeometry.regionNames()[geomI];
335  forAll(regionNames, regionI)
336  {
337  maxLen = Foam::max(maxLen, label(regionNames[regionI].size()));
338  }
339  }
340 
341 
342  Info<< setw(maxLen) << "Region"
343  << setw(10) << "Min Level"
344  << setw(10) << "Max Level"
345  << setw(10) << "Gap Level" << nl
346  << setw(maxLen) << "------"
347  << setw(10) << "---------"
348  << setw(10) << "---------"
349  << setw(10) << "---------" << endl;
350 
351  forAll(rf.surfaces(), surfI)
352  {
353  label geomI = rf.surfaces()[surfI];
354 
355  Info<< rf.names()[surfI] << ':' << nl;
356 
357  const wordList& regionNames = allGeometry.regionNames()[geomI];
358 
359  forAll(regionNames, regionI)
360  {
361  label globalI = rf.globalRegion(surfI, regionI);
362 
363  Info<< setw(maxLen) << regionNames[regionI]
364  << setw(10) << rf.minLevel()[globalI]
365  << setw(10) << rf.maxLevel()[globalI]
366  << setw(10) << rf.gapLevel()[globalI] << endl;
367  }
368  }
369 
370 
371  return surfacePtr;
372 }
373 
374 
375 void extractSurface
376 (
377  const polyMesh& mesh,
378  const Time& runTime,
380  const fileName& outFileName
381 )
382 {
383  const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
384 
385  // Collect sizes. Hash on names to handle local-only patches (e.g.
386  // processor patches)
387  HashTable<label> patchSize(1000);
388  label nFaces = 0;
390  {
391  const polyPatch& pp = bMesh[iter.key()];
392  patchSize.insert(pp.name(), pp.size());
393  nFaces += pp.size();
394  }
396 
397 
398  // Allocate zone/patch for all patches
399  HashTable<label> compactZoneID(1000);
400  forAllConstIter(HashTable<label>, patchSize, iter)
401  {
402  label sz = compactZoneID.size();
403  compactZoneID.insert(iter.key(), sz);
404  }
405  Pstream::mapCombineScatter(compactZoneID);
406 
407 
408  // Rework HashTable into labelList just for speed of conversion
409  labelList patchToCompactZone(bMesh.size(), -1);
410  forAllConstIter(HashTable<label>, compactZoneID, iter)
411  {
412  label patchi = bMesh.findPatchID(iter.key());
413  if (patchi != -1)
414  {
415  patchToCompactZone[patchi] = iter();
416  }
417  }
418 
419  // Collect faces on zones
420  DynamicList<label> faceLabels(nFaces);
421  DynamicList<label> compactZones(nFaces);
423  {
424  const polyPatch& pp = bMesh[iter.key()];
425  forAll(pp, i)
426  {
427  faceLabels.append(pp.start()+i);
428  compactZones.append(patchToCompactZone[pp.index()]);
429  }
430  }
431 
432  // Addressing engine for all faces
433  uindirectPrimitivePatch allBoundary
434  (
435  UIndirectList<face>(mesh.faces(), faceLabels),
436  mesh.points()
437  );
438 
439 
440  // Find correspondence to master points
441  labelList pointToGlobal;
442  labelList uniqueMeshPoints;
443  autoPtr<globalIndex> globalNumbers = mesh.globalData().mergePoints
444  (
445  allBoundary.meshPoints(),
446  allBoundary.meshPointMap(),
447  pointToGlobal,
448  uniqueMeshPoints
449  );
450 
451  // Gather all unique points on master
452  List<pointField> gatheredPoints(Pstream::nProcs());
453  gatheredPoints[Pstream::myProcNo()] = pointField
454  (
455  mesh.points(),
456  uniqueMeshPoints
457  );
458  Pstream::gatherList(gatheredPoints);
459 
460  // Gather all faces
461  List<faceList> gatheredFaces(Pstream::nProcs());
462  gatheredFaces[Pstream::myProcNo()] = allBoundary.localFaces();
463  forAll(gatheredFaces[Pstream::myProcNo()], i)
464  {
465  inplaceRenumber(pointToGlobal, gatheredFaces[Pstream::myProcNo()][i]);
466  }
467  Pstream::gatherList(gatheredFaces);
468 
469  // Gather all ZoneIDs
470  List<labelList> gatheredZones(Pstream::nProcs());
471  gatheredZones[Pstream::myProcNo()] = move(compactZones);
472  Pstream::gatherList(gatheredZones);
473 
474  // On master combine all points, faces, zones
475  if (Pstream::master())
476  {
477  pointField allPoints = ListListOps::combine<pointField>
478  (
479  gatheredPoints,
481  );
482  gatheredPoints.clear();
483 
484  faceList allFaces = ListListOps::combine<faceList>
485  (
486  gatheredFaces,
488  );
489  gatheredFaces.clear();
490 
491  labelList allZones = ListListOps::combine<labelList>
492  (
493  gatheredZones,
495  );
496  gatheredZones.clear();
497 
498 
499  // Zones
500  surfZoneIdentifierList surfZones(compactZoneID.size());
501  forAllConstIter(HashTable<label>, compactZoneID, iter)
502  {
503  surfZones[iter()] = surfZoneIdentifier(iter.key(), iter());
504  Info<< "surfZone " << iter() << " : " << surfZones[iter()].name()
505  << endl;
506  }
507 
508  UnsortedMeshedSurface<face> unsortedFace
509  (
510  move(allPoints),
511  move(allFaces),
512  move(allZones),
513  move(surfZones)
514  );
515 
516 
517  MeshedSurface<face> sortedFace(unsortedFace);
518 
519  Info<< "Writing merged surface to "
520  << runTime.globalPath()/outFileName << endl;
521 
522  sortedFace.write(runTime.globalPath()/outFileName);
523  }
524 }
525 
526 
527 // Check writing tolerance before doing any serious work
528 scalar getMergeDistance(const polyMesh& mesh, const scalar mergeTol)
529 {
530  const boundBox& meshBb = mesh.bounds();
531  scalar mergeDist = mergeTol * meshBb.mag();
532 
533  Info<< nl
534  << "Overall mesh bounding box : " << meshBb << nl
535  << "Relative tolerance : " << mergeTol << nl
536  << "Absolute matching distance : " << mergeDist << nl
537  << endl;
538 
539  // check writing tolerance
540  if (mesh.time().writeFormat() == IOstream::ASCII)
541  {
542  const scalar writeTol = std::pow
543  (
544  scalar(10.0),
545  -scalar(IOstream::defaultPrecision())
546  );
547 
548  if (mergeTol < writeTol)
549  {
551  << "Your current settings specify ASCII writing with "
552  << IOstream::defaultPrecision() << " digits precision." << nl
553  << "Your merging tolerance (" << mergeTol
554  << ") is finer than this." << nl
555  << "Change to binary writeFormat, "
556  << "or increase the writePrecision" << endl
557  << "or adjust the merge tolerance (mergeTol)."
558  << exit(FatalError);
559  }
560  }
561 
562  return mergeDist;
563 }
564 
565 
566 void removeZeroSizedPatches(fvMesh& mesh)
567 {
568  // Remove non-constraint zero-sized patches
569 
570  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
571 
572  labelList oldToNew(pbm.size(), -1);
573  label newPatchi = 0;
574  forAll(pbm, patchi)
575  {
576  const polyPatch& pp = pbm[patchi];
577 
578  if
579  (
580  polyPatch::constraintType(pp.type())
581  || returnReduce(pp.size(), sumOp<label>())
582  )
583  {
584  oldToNew[patchi] = newPatchi++;
585  }
586  }
587 
588  const label nKeepPatches = newPatchi;
589 
590  // Shuffle unused ones to end
591  if (nKeepPatches != pbm.size())
592  {
593  Info<< endl
594  << "Removing zero-sized patches:" << endl << incrIndent;
595 
596  forAll(oldToNew, patchi)
597  {
598  if (oldToNew[patchi] == -1)
599  {
600  Info<< indent << pbm[patchi].name()
601  << " type " << pbm[patchi].type()
602  << " at position " << patchi << endl;
603  oldToNew[patchi] = newPatchi++;
604  }
605  }
606  Info<< decrIndent;
607 
608  fvMeshTools::reorderPatches(mesh, oldToNew, nKeepPatches, true);
609  Info<< endl;
610  }
611 }
612 
613 
614 // Write mesh and additional information
615 void writeMesh
616 (
617  const string& msg,
618  const meshRefinement& meshRefiner,
619  const meshRefinement::debugType debugLevel,
620  const meshRefinement::writeType writeLevel
621 )
622 {
623  const fvMesh& mesh = meshRefiner.mesh();
624 
625  meshRefiner.printMeshInfo(debugLevel, msg);
626  Info<< "Writing mesh to time " << meshRefiner.name() << endl;
627 
628  meshRefiner.write
629  (
630  debugLevel,
632  mesh.time().path()/meshRefiner.name()
633  );
634  Info<< "Wrote mesh in = "
635  << mesh.time().cpuTimeIncrement() << " s." << endl;
636 }
637 
638 
639 int main(int argc, char *argv[])
640 {
641  #include "addOverwriteOption.H"
643  (
644  "checkGeometry",
645  "check all surface geometry for quality"
646  );
648  (
649  "surfaceSimplify",
650  "boundBox",
651  "simplify the surface using snappyHexMesh starting from a boundBox"
652  );
654  (
655  "patches",
656  "(patch0 .. patchN)",
657  "only triangulate selected patches (wildcards supported)"
658  );
660  (
661  "outFile",
662  "file",
663  "name of the file to save the simplified surface to"
664  );
665  #include "addDictOption.H"
666 
667  #include "setRootCase.H"
669 
670  const bool overwrite = args.optionFound("overwrite");
671  const bool checkGeometry = args.optionFound("checkGeometry");
672  const bool surfaceSimplify = args.optionFound("surfaceSimplify");
673 
675 
676  {
677  Foam::Info
678  << "Create mesh for time = "
679  << runTime.name() << Foam::nl << Foam::endl;
680 
681  meshPtr.set
682  (
683  new fvMesh
684  (
686  (
688  runTime.name(),
689  runTime,
691  ),
692  false
693  )
694  );
695  }
696 
697  fvMesh& mesh = meshPtr();
698 
699  Info<< "Read mesh in = "
700  << runTime.cpuTimeIncrement() << " s" << endl;
701 
702  // Check that the read mesh is fully 3D
703  // as required for mesh relaxation after snapping
704  if (mesh.nSolutionD() != 3)
705  {
707  << "Mesh provided is not fully 3D"
708  " as required for mesh relaxation after snapping" << nl
709  << "Convert all empty patches to appropriate types for a 3D mesh,"
710  " current patch types are" << nl << mesh.boundaryMesh().types()
711  << exit(FatalError);
712  }
713 
714  // Check patches and faceZones are synchronised
715  mesh.boundaryMesh().checkParallelSync(true);
717 
718 
719  // Read meshing dictionary
720  const IOdictionary meshDict(systemDict("snappyHexMeshDict", args, mesh));
721 
722 
723  // all surface geometry
724  const dictionary& geometryDict = meshDict.subDict("geometry");
725 
726  // refinement parameters
727  const dictionary& refineDict = meshDict.subDict("castellatedMeshControls");
728 
729  // mesh motion and mesh quality parameters
730  const dictionary& motionDict = meshDict.subDict("meshQualityControls");
731 
732  // snap-to-surface parameters
733  const dictionary& snapDict = meshDict.subDict("snapControls");
734 
735  // absolute merge distance
736  const scalar mergeDist = getMergeDistance
737  (
738  mesh,
739  meshDict.lookup<scalar>("mergeTolerance")
740  );
741 
742  const Switch keepPatches(meshDict.lookupOrDefault("keepPatches", false));
743 
744 
745  // Read decomposePar dictionary
746  dictionary decomposeDict;
747  {
748  if (Pstream::parRun())
749  {
750  decomposeDict = decompositionMethod::decomposeParDict(runTime);
751  }
752  else
753  {
754  decomposeDict.add("method", "none");
755  decomposeDict.add("numberOfSubdomains", 1);
756  }
757  }
758 
759 
760  // Debug
761  // ~~~~~
762 
763  // Set debug level
765  (
766  meshDict.lookupOrDefault<label>
767  (
768  "debug",
769  0
770  )
771  );
772  {
773  wordList flags;
774  if (meshDict.readIfPresent("debugFlags", flags))
775  {
776  debugLevel = meshRefinement::debugType
777  (
779  (
781  flags
782  )
783  );
784  }
785  }
786  if (debugLevel > 0)
787  {
788  meshRefinement::debug = debugLevel;
789  snappyRefineDriver::debug = debugLevel;
790  snappySnapDriver::debug = debugLevel;
791  snappyLayerDriver::debug = debugLevel;
792  }
793 
794  // Set file writing level
795  {
796  wordList flags;
797  if (meshDict.readIfPresent("writeFlags", flags))
798  {
800  (
802  (
804  (
806  flags
807  )
808  )
809  );
810  }
811  }
812 
813  // Set output level
814  {
815  wordList flags;
816  if (meshDict.readIfPresent("outputFlags", flags))
817  {
819  (
821  (
823  (
825  flags
826  )
827  )
828  );
829  }
830  }
831 
832 
833  // Read geometry
834  // ~~~~~~~~~~~~~
835 
836  searchableSurfaces allGeometry
837  (
838  IOobject
839  (
840  "abc",
841  mesh.time().constant(),
843  mesh.time(),
846  ),
847  geometryDict,
848  meshDict.lookupOrDefault("singleRegionName", true)
849  );
850 
851 
852  // Read refinement surfaces
853  // ~~~~~~~~~~~~~~~~~~~~~~~~
854 
855  Info<< "Reading refinement surfaces..." << endl;
856 
857  refinementSurfaces surfaces
858  (
859  allGeometry,
860  refineDict.found("refinementSurfaces")
861  ? refineDict.subDict("refinementSurfaces")
863  refineDict.lookupOrDefault("gapLevelIncrement", 0)
864  );
865 
866  Info<< "Read refinement surfaces in = "
867  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
868 
869 
870  // Checking only?
871 
872  if (checkGeometry)
873  {
874  // Extract patchInfo
875  List<wordList> patchTypes(allGeometry.size());
876 
877  const PtrList<dictionary>& patchInfo = surfaces.patchInfo();
878  const labelList& surfaceGeometry = surfaces.surfaces();
879  forAll(surfaceGeometry, surfI)
880  {
881  label geomI = surfaceGeometry[surfI];
882  const wordList& regNames = allGeometry.regionNames()[geomI];
883 
884  patchTypes[geomI].setSize(regNames.size());
885  forAll(regNames, regionI)
886  {
887  label globalRegionI = surfaces.globalRegion(surfI, regionI);
888 
889  if (patchInfo.set(globalRegionI))
890  {
891  patchTypes[geomI][regionI] =
892  word(patchInfo[globalRegionI].lookup("type"));
893  }
894  else
895  {
896  patchTypes[geomI][regionI] = wallPolyPatch::typeName;
897  }
898  }
899  }
900 
901  // Write some stats
902  allGeometry.writeStats(patchTypes, Info);
903  // Check topology
904  allGeometry.checkTopology(true);
905  // Check geometry
906  allGeometry.checkGeometry
907  (
908  100.0, // max size ratio
909  1e-9, // intersection tolerance
910  0.01, // min triangle quality
911  true
912  );
913 
914  return 0;
915  }
916 
917 
918 
919  // Read refinement shells
920  // ~~~~~~~~~~~~~~~~~~~~~~
921 
922  Info<< "Reading refinement regions..." << endl;
923  refinementRegions shells
924  (
925  allGeometry,
926  refineDict.found("refinementRegions")
927  ? refineDict.subDict("refinementRegions")
929  );
930  Info<< "Read refinement regions in = "
931  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
932 
933 
934 
935  // Read feature meshes
936  // ~~~~~~~~~~~~~~~~~~~
937 
938  Info<< "Reading features..." << endl;
939  refinementFeatures features
940  (
941  mesh,
942  refineDict.found("features")
943  ? refineDict.lookup("features")
945  );
946  Info<< "Read features in = "
947  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
948 
949 
950 
951  // Refinement engine
952  // ~~~~~~~~~~~~~~~~~
953 
954  Info<< nl
955  << "Determining initial surface intersections" << nl
956  << "-----------------------------------------" << nl
957  << endl;
958 
959  // Main refinement engine
960  meshRefinement meshRefiner
961  (
962  mesh,
963  refineDict,
964  mergeDist, // tolerance used in sorting coordinates
965  overwrite, // overwrite mesh files?
966  surfaces, // for surface intersection refinement
967  features, // for feature edges/point based refinement
968  shells // for volume (inside/outside) refinement
969  );
970  Info<< "Calculated surface intersections in = "
971  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
972 
973  // Some stats
974  meshRefiner.printMeshInfo(debugLevel, "Initial mesh");
975 
976  meshRefiner.write
977  (
980  mesh.time().path()/meshRefiner.name()
981  );
982 
983 
984  // Add all the surface regions as patches
985  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
986 
987  //- Global surface region to patch (non faceZone surface) or patches
988  // (faceZone surfaces)
989  labelList globalToMasterPatch;
990  labelList globalToSlavePatch;
991  {
992  Info<< nl
993  << "Adding patches for surface regions" << nl
994  << "----------------------------------" << nl
995  << endl;
996 
997  // From global region number to mesh patch.
998  globalToMasterPatch.setSize(surfaces.nRegions(), -1);
999  globalToSlavePatch.setSize(surfaces.nRegions(), -1);
1000 
1001  Info<< setf(ios_base::left)
1002  << setw(6) << "Patch"
1003  << setw(20) << "Type"
1004  << setw(30) << "Region" << nl
1005  << setw(6) << "-----"
1006  << setw(20) << "----"
1007  << setw(30) << "------" << endl;
1008 
1009  const labelList& surfaceGeometry = surfaces.surfaces();
1010  const PtrList<dictionary>& surfacePatchInfo = surfaces.patchInfo();
1011 
1012  forAll(surfaceGeometry, surfI)
1013  {
1014  label geomI = surfaceGeometry[surfI];
1015 
1016  const wordList& regNames = allGeometry.regionNames()[geomI];
1017 
1018  Info<< surfaces.names()[surfI] << ':' << nl << nl;
1019 
1020  if (surfaces.surfZones()[surfI].faceZoneName().empty())
1021  {
1022  // 'Normal' surface
1023  forAll(regNames, i)
1024  {
1025  label globalRegionI = surfaces.globalRegion(surfI, i);
1026 
1027  label patchi;
1028 
1029  if (surfacePatchInfo.set(globalRegionI))
1030  {
1031  patchi = meshRefiner.addMeshedPatch
1032  (
1033  regNames[i],
1034  surfacePatchInfo[globalRegionI]
1035  );
1036  }
1037  else
1038  {
1039  dictionary patchInfo;
1040  patchInfo.set("type", wallPolyPatch::typeName);
1041 
1042  patchi = meshRefiner.addMeshedPatch
1043  (
1044  regNames[i],
1045  patchInfo
1046  );
1047  }
1048 
1049  Info<< setf(ios_base::left)
1050  << setw(6) << patchi
1051  << setw(20) << mesh.boundaryMesh()[patchi].type()
1052  << setw(30) << regNames[i] << nl;
1053 
1054  globalToMasterPatch[globalRegionI] = patchi;
1055  globalToSlavePatch[globalRegionI] = patchi;
1056  }
1057  }
1058  else
1059  {
1060  // Zoned surface
1061  forAll(regNames, i)
1062  {
1063  label globalRegionI = surfaces.globalRegion(surfI, i);
1064 
1065  // Add master side patch
1066  {
1067  label patchi;
1068 
1069  if (surfacePatchInfo.set(globalRegionI))
1070  {
1071  patchi = meshRefiner.addMeshedPatch
1072  (
1073  regNames[i],
1074  surfacePatchInfo[globalRegionI]
1075  );
1076  }
1077  else
1078  {
1079  dictionary patchInfo;
1080  patchInfo.set("type", wallPolyPatch::typeName);
1081 
1082  patchi = meshRefiner.addMeshedPatch
1083  (
1084  regNames[i],
1085  patchInfo
1086  );
1087  }
1088 
1089  Info<< setf(ios_base::left)
1090  << setw(6) << patchi
1091  << setw(20) << mesh.boundaryMesh()[patchi].type()
1092  << setw(30) << regNames[i] << nl;
1093 
1094  globalToMasterPatch[globalRegionI] = patchi;
1095  }
1096  // Add slave side patch
1097  {
1098  const word slaveName = regNames[i] + "_slave";
1099  label patchi;
1100 
1101  if (surfacePatchInfo.set(globalRegionI))
1102  {
1103  patchi = meshRefiner.addMeshedPatch
1104  (
1105  slaveName,
1106  surfacePatchInfo[globalRegionI]
1107  );
1108  }
1109  else
1110  {
1111  dictionary patchInfo;
1112  patchInfo.set("type", wallPolyPatch::typeName);
1113 
1114  patchi = meshRefiner.addMeshedPatch
1115  (
1116  slaveName,
1117  patchInfo
1118  );
1119  }
1120 
1121  Info<< setf(ios_base::left)
1122  << setw(6) << patchi
1123  << setw(20) << mesh.boundaryMesh()[patchi].type()
1124  << setw(30) << slaveName << nl;
1125 
1126  globalToSlavePatch[globalRegionI] = patchi;
1127  }
1128  }
1129  }
1130 
1131  Info<< nl;
1132  }
1133  Info<< "Added patches in = "
1134  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
1135  }
1136 
1137 
1138  // Parallel
1139  // ~~~~~~~~
1140 
1141  // Decomposition
1142  autoPtr<decompositionMethod> decomposerPtr
1143  (
1145  );
1146  decompositionMethod& decomposer = decomposerPtr();
1147 
1148  // Mesh distribution engine (uses tolerance to reconstruct meshes)
1149  fvMeshDistribute distributor(mesh);
1150 
1151 
1152  // Now do the real work -refinement -snapping -layers
1153  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1154 
1155  const Switch wantRefine(meshDict.lookup("castellatedMesh"));
1156  const Switch wantSnap(meshDict.lookup("snap"));
1157  const Switch wantLayers(meshDict.lookup("addLayers"));
1158 
1159  // Refinement parameters
1160  const refinementParameters refineParams(refineDict);
1161 
1162  // Snap parameters
1163  const snapParameters snapParams(snapDict);
1164 
1165  if (wantRefine)
1166  {
1167  cpuTime timer;
1168 
1169  snappyRefineDriver refineDriver
1170  (
1171  meshRefiner,
1172  decomposer,
1173  distributor,
1174  globalToMasterPatch,
1175  globalToSlavePatch
1176  );
1177 
1178  if (!overwrite && !debugLevel)
1179  {
1180  const_cast<Time&>(mesh.time())++;
1181  }
1182 
1183  refineDriver.doRefine
1184  (
1185  refineDict,
1186  refineParams,
1187  snapParams,
1188  refineParams.handleSnapProblems(),
1189  motionDict
1190  );
1191 
1192  if (!keepPatches && !wantSnap && !wantLayers)
1193  {
1194  removeZeroSizedPatches(mesh);
1195  }
1196 
1197  writeMesh
1198  (
1199  "Refined mesh",
1200  meshRefiner,
1201  debugLevel,
1203  );
1204 
1205  Info<< "Mesh refined in = "
1206  << timer.cpuTimeIncrement() << " s." << endl;
1207  }
1208 
1209  if (wantSnap)
1210  {
1211  cpuTime timer;
1212 
1213  snappySnapDriver snapDriver
1214  (
1215  meshRefiner,
1216  globalToMasterPatch,
1217  globalToSlavePatch
1218  );
1219 
1220  if (!overwrite && !debugLevel)
1221  {
1222  const_cast<Time&>(mesh.time())++;
1223  }
1224 
1225  // Use the resolveFeatureAngle from the refinement parameters
1226  scalar curvature = refineParams.curvature();
1227  scalar planarAngle = refineParams.planarAngle();
1228 
1229  snapDriver.doSnap
1230  (
1231  snapDict,
1232  motionDict,
1233  curvature,
1234  planarAngle,
1235  snapParams
1236  );
1237 
1238  if (!keepPatches && !wantLayers)
1239  {
1240  removeZeroSizedPatches(mesh);
1241  }
1242 
1243  writeMesh
1244  (
1245  "Snapped mesh",
1246  meshRefiner,
1247  debugLevel,
1249  );
1250 
1251  Info<< "Mesh snapped in = "
1252  << timer.cpuTimeIncrement() << " s." << endl;
1253  }
1254 
1255  if (wantLayers)
1256  {
1257  cpuTime timer;
1258 
1259  // Layer addition parameters dictionary
1260  const dictionary& layersDict = meshDict.subDict("addLayersControls");
1261 
1262  // Layer addition parameters
1263  const layerParameters layerParams(layersDict, mesh.boundaryMesh());
1264 
1265  snappyLayerDriver layerDriver
1266  (
1267  meshRefiner,
1268  globalToMasterPatch,
1269  globalToSlavePatch
1270  );
1271 
1272  // Use the maxLocalCells from the refinement parameters
1273  bool preBalance = returnReduce
1274  (
1275  (mesh.nCells() >= refineParams.maxLocalCells()),
1276  orOp<bool>()
1277  );
1278 
1279 
1280  if (!overwrite && !debugLevel)
1281  {
1282  const_cast<Time&>(mesh.time())++;
1283  }
1284 
1285  layerDriver.doLayers
1286  (
1287  layersDict,
1288  motionDict,
1289  layerParams,
1290  preBalance,
1291  decomposer,
1292  distributor
1293  );
1294 
1295  if (!keepPatches)
1296  {
1297  removeZeroSizedPatches(mesh);
1298  }
1299 
1300  writeMesh
1301  (
1302  "Layer mesh",
1303  meshRefiner,
1304  debugLevel,
1306  );
1307 
1308  Info<< "Layers added in = "
1309  << timer.cpuTimeIncrement() << " s." << endl;
1310  }
1311 
1312 
1313  {
1314  // Check final mesh
1315  Info<< "Checking final mesh ..." << endl;
1316  faceSet wrongFaces(mesh, "wrongFaces", mesh.nFaces()/100);
1317  motionSmoother::checkMesh(false, mesh, motionDict, wrongFaces);
1318  const label nErrors = returnReduce
1319  (
1320  wrongFaces.size(),
1321  sumOp<label>()
1322  );
1323 
1324  if (nErrors > 0)
1325  {
1326  Info<< "Finished meshing with " << nErrors << " illegal faces"
1327  << " (concave, zero area or negative cell pyramid volume)"
1328  << endl;
1329  wrongFaces.write();
1330  }
1331  else
1332  {
1333  Info<< "Finished meshing without any errors" << endl;
1334  }
1335  }
1336 
1337 
1338  if (surfaceSimplify)
1339  {
1340  const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
1341 
1343 
1344  if (args.optionFound("patches"))
1345  {
1347  (
1348  wordReList(args.optionLookup("patches")())
1349  );
1350  }
1351  else
1352  {
1353  forAll(bMesh, patchi)
1354  {
1355  const polyPatch& patch = bMesh[patchi];
1356 
1357  if (!isA<processorPolyPatch>(patch))
1358  {
1360  }
1361  }
1362  }
1363 
1364  fileName outFileName
1365  (
1367  (
1368  "outFile",
1369  "constant/triSurface/simplifiedSurface.stl"
1370  )
1371  );
1372 
1373  extractSurface
1374  (
1375  mesh,
1376  runTime,
1378  outFileName
1379  );
1380 
1381  pointIOField cellCentres
1382  (
1383  IOobject
1384  (
1385  "internalCellCentres",
1386  runTime.name(),
1387  mesh,
1390  ),
1391  mesh.cellCentres()
1392  );
1393 
1394  cellCentres.write();
1395  }
1396 
1397 
1398  Info<< "Finished meshing in = "
1399  << runTime.elapsedCpuTime() << " s." << endl;
1400 
1401  Info<< "End\n" << endl;
1402 
1403  return 0;
1404 }
1405 
1406 
1407 // ************************************************************************* //
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:111
bool set(const Key &, const T &newElmt)
Assign a new hashedEntry, overwriting existing entries.
Definition: HashTableI.H:91
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:65
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:283
A List with indirect addressing.
Definition: UIndirectList.H:60
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
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 autoPtr< decompositionMethod > NewDistributor(const dictionary &decompositionDict)
Return a reference to the selected decomposition method.
static dictionary decomposeParDict(const Time &time)
Read and return the decomposeParDict.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
autoPtr< dictionary > clone() const
Construct and return clone.
Definition: dictionary.C:548
const entry * lookupEntryPtr(const word &, bool recursive, bool patternMatch) const
Find and return an entry data stream pointer if present.
Definition: dictionary.C:698
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:1307
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:1076
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:998
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1169
wordList toc() const
Return the table of contents.
Definition: dictionary.C:1131
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:659
static const dictionary null
Null dictionary.
Definition: dictionary.H:242
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:101
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:402
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.
static bool checkMesh(const bool report, const polyMesh &mesh, const dictionary &dict, labelHashSet &wrongFaces)
Check mesh with mesh settings in dict. Collects incorrect faces.
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.
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 IDs corresponding to the given names.
label findPatchID(const word &patchName) const
Find patch index given a name.
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
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:268
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1373
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1562
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:1071
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:405
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1360
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:411
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:301
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
int main(int argc, char *argv[])
Definition: financialFoam.C:44
const polyBoundaryMesh & bMesh
label patchi
static fvMesh * meshPtr
Definition: globalFoam.H:52
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:105
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
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:251
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:228
label checkGeometry(const polyMesh &mesh, const bool allGeometry, const autoPtr< surfaceWriter > &, const autoPtr< setWriter > &)
Check the geometry.
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:221
static const char nl
Definition: Ostream.H:260
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)