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