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-2019 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 "shellSurfaces.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 "vtkSetWriter.H"
50 #include "faceSet.H"
51 #include "motionSmoother.H"
52 #include "polyTopoChange.H"
53 #include "cellModeller.H"
55 #include "surfZoneIdentifierList.H"
56 #include "UnsortedMeshedSurface.H"
57 #include "MeshedSurface.H"
58 #include "globalIndex.H"
59 #include "IOmanip.H"
60 #include "fvMeshTools.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  readScalar(scsDict.lookup("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  readScalar
232  (
233  scsDict.lookup("surfaceCellSizeCoeff")
234  );
235 
236  const label refLevel = sizeCoeffToRefinement
237  (
238  level0Coeff,
239  surfaceCellSize
240  );
241 
242  regionMinLevel[surfI].insert(regionI, refLevel);
243  regionMaxLevel[surfI].insert(regionI, refLevel);
244  regionLevelIncr[surfI].insert(regionI, 0);
245  }
246  }
247  }
248 
249  surfI++;
250  }
251  }
252 
253  // Calculate local to global region offset
254  label nRegions = 0;
255 
256  forAll(surfaces, surfI)
257  {
258  regionOffset[surfI] = nRegions;
259  nRegions += allGeometry[surfaces[surfI]].regions().size();
260  }
261 
262  // Rework surface specific information into information per global region
263  labelList minLevel(nRegions, 0);
264  labelList maxLevel(nRegions, 0);
265  labelList gapLevel(nRegions, -1);
266  PtrList<dictionary> patchInfo(nRegions);
267 
268  forAll(globalMinLevel, surfI)
269  {
270  label nRegions = allGeometry[surfaces[surfI]].regions().size();
271 
272  // Initialise to global (i.e. per surface)
273  for (label i = 0; i < nRegions; i++)
274  {
275  label globalRegionI = regionOffset[surfI] + i;
276  minLevel[globalRegionI] = globalMinLevel[surfI];
277  maxLevel[globalRegionI] = globalMaxLevel[surfI];
278  gapLevel[globalRegionI] =
279  maxLevel[globalRegionI]
280  + globalLevelIncr[surfI];
281 
282  if (globalPatchInfo.set(surfI))
283  {
284  patchInfo.set
285  (
286  globalRegionI,
287  globalPatchInfo[surfI].clone()
288  );
289  }
290  }
291 
292  // Overwrite with region specific information
293  forAllConstIter(Map<label>, regionMinLevel[surfI], iter)
294  {
295  label globalRegionI = regionOffset[surfI] + iter.key();
296 
297  minLevel[globalRegionI] = iter();
298  maxLevel[globalRegionI] = regionMaxLevel[surfI][iter.key()];
299  gapLevel[globalRegionI] =
300  maxLevel[globalRegionI]
301  + regionLevelIncr[surfI][iter.key()];
302  }
303 
304  const Map<autoPtr<dictionary>>& localInfo = regionPatchInfo[surfI];
305  forAllConstIter(Map<autoPtr<dictionary>>, localInfo, iter)
306  {
307  label globalRegionI = regionOffset[surfI] + iter.key();
308  patchInfo.set(globalRegionI, iter()().clone());
309  }
310  }
311 
312  surfacePtr.set
313  (
315  (
316  allGeometry,
317  surfaces,
318  names,
319  surfZones,
320  regionOffset,
321  minLevel,
322  maxLevel,
323  gapLevel,
324  scalarField(nRegions, -great), // perpendicularAngle,
325  patchInfo
326  )
327  );
328 
329 
330  const refinementSurfaces& rf = surfacePtr();
331 
332  // Determine maximum region name length
333  label maxLen = 0;
334  forAll(rf.surfaces(), surfI)
335  {
336  label geomI = rf.surfaces()[surfI];
337  const wordList& regionNames = allGeometry.regionNames()[geomI];
338  forAll(regionNames, regionI)
339  {
340  maxLen = Foam::max(maxLen, label(regionNames[regionI].size()));
341  }
342  }
343 
344 
345  Info<< setw(maxLen) << "Region"
346  << setw(10) << "Min Level"
347  << setw(10) << "Max Level"
348  << setw(10) << "Gap Level" << nl
349  << setw(maxLen) << "------"
350  << setw(10) << "---------"
351  << setw(10) << "---------"
352  << setw(10) << "---------" << endl;
353 
354  forAll(rf.surfaces(), surfI)
355  {
356  label geomI = rf.surfaces()[surfI];
357 
358  Info<< rf.names()[surfI] << ':' << nl;
359 
360  const wordList& regionNames = allGeometry.regionNames()[geomI];
361 
362  forAll(regionNames, regionI)
363  {
364  label globalI = rf.globalRegion(surfI, regionI);
365 
366  Info<< setw(maxLen) << regionNames[regionI]
367  << setw(10) << rf.minLevel()[globalI]
368  << setw(10) << rf.maxLevel()[globalI]
369  << setw(10) << rf.gapLevel()[globalI] << endl;
370  }
371  }
372 
373 
374  return surfacePtr;
375 }
376 
377 
378 void extractSurface
379 (
380  const polyMesh& mesh,
381  const Time& runTime,
382  const labelHashSet& includePatches,
383  const fileName& outFileName
384 )
385 {
386  const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
387 
388  // Collect sizes. Hash on names to handle local-only patches (e.g.
389  // processor patches)
390  HashTable<label> patchSize(1000);
391  label nFaces = 0;
392  forAllConstIter(labelHashSet, includePatches, iter)
393  {
394  const polyPatch& pp = bMesh[iter.key()];
395  patchSize.insert(pp.name(), pp.size());
396  nFaces += pp.size();
397  }
399 
400 
401  // Allocate zone/patch for all patches
402  HashTable<label> compactZoneID(1000);
403  forAllConstIter(HashTable<label>, patchSize, iter)
404  {
405  label sz = compactZoneID.size();
406  compactZoneID.insert(iter.key(), sz);
407  }
408  Pstream::mapCombineScatter(compactZoneID);
409 
410 
411  // Rework HashTable into labelList just for speed of conversion
412  labelList patchToCompactZone(bMesh.size(), -1);
413  forAllConstIter(HashTable<label>, compactZoneID, iter)
414  {
415  label patchi = bMesh.findPatchID(iter.key());
416  if (patchi != -1)
417  {
418  patchToCompactZone[patchi] = iter();
419  }
420  }
421 
422  // Collect faces on zones
423  DynamicList<label> faceLabels(nFaces);
424  DynamicList<label> compactZones(nFaces);
425  forAllConstIter(labelHashSet, includePatches, iter)
426  {
427  const polyPatch& pp = bMesh[iter.key()];
428  forAll(pp, i)
429  {
430  faceLabels.append(pp.start()+i);
431  compactZones.append(patchToCompactZone[pp.index()]);
432  }
433  }
434 
435  // Addressing engine for all faces
436  uindirectPrimitivePatch allBoundary
437  (
438  UIndirectList<face>(mesh.faces(), faceLabels),
439  mesh.points()
440  );
441 
442 
443  // Find correspondence to master points
444  labelList pointToGlobal;
445  labelList uniqueMeshPoints;
446  autoPtr<globalIndex> globalNumbers = mesh.globalData().mergePoints
447  (
448  allBoundary.meshPoints(),
449  allBoundary.meshPointMap(),
450  pointToGlobal,
451  uniqueMeshPoints
452  );
453 
454  // Gather all unique points on master
455  List<pointField> gatheredPoints(Pstream::nProcs());
456  gatheredPoints[Pstream::myProcNo()] = pointField
457  (
458  mesh.points(),
459  uniqueMeshPoints
460  );
461  Pstream::gatherList(gatheredPoints);
462 
463  // Gather all faces
464  List<faceList> gatheredFaces(Pstream::nProcs());
465  gatheredFaces[Pstream::myProcNo()] = allBoundary.localFaces();
466  forAll(gatheredFaces[Pstream::myProcNo()], i)
467  {
468  inplaceRenumber(pointToGlobal, gatheredFaces[Pstream::myProcNo()][i]);
469  }
470  Pstream::gatherList(gatheredFaces);
471 
472  // Gather all ZoneIDs
473  List<labelList> gatheredZones(Pstream::nProcs());
474  gatheredZones[Pstream::myProcNo()] = move(compactZones);
475  Pstream::gatherList(gatheredZones);
476 
477  // On master combine all points, faces, zones
478  if (Pstream::master())
479  {
480  pointField allPoints = ListListOps::combine<pointField>
481  (
482  gatheredPoints,
484  );
485  gatheredPoints.clear();
486 
487  faceList allFaces = ListListOps::combine<faceList>
488  (
489  gatheredFaces,
491  );
492  gatheredFaces.clear();
493 
494  labelList allZones = ListListOps::combine<labelList>
495  (
496  gatheredZones,
498  );
499  gatheredZones.clear();
500 
501 
502  // Zones
503  surfZoneIdentifierList surfZones(compactZoneID.size());
504  forAllConstIter(HashTable<label>, compactZoneID, iter)
505  {
506  surfZones[iter()] = surfZoneIdentifier(iter.key(), iter());
507  Info<< "surfZone " << iter() << " : " << surfZones[iter()].name()
508  << endl;
509  }
510 
511  UnsortedMeshedSurface<face> unsortedFace
512  (
513  move(allPoints),
514  move(allFaces),
515  move(allZones),
516  move(surfZones)
517  );
518 
519 
520  MeshedSurface<face> sortedFace(unsortedFace);
521 
522  fileName globalCasePath
523  (
524  runTime.processorCase()
525  ? runTime.path()/".."/outFileName
526  : runTime.path()/outFileName
527  );
528  globalCasePath.clean();
529 
530  Info<< "Writing merged surface to " << globalCasePath << endl;
531 
532  sortedFace.write(globalCasePath);
533  }
534 }
535 
536 
537 // Check writing tolerance before doing any serious work
538 scalar getMergeDistance(const polyMesh& mesh, const scalar mergeTol)
539 {
540  const boundBox& meshBb = mesh.bounds();
541  scalar mergeDist = mergeTol * meshBb.mag();
542 
543  Info<< nl
544  << "Overall mesh bounding box : " << meshBb << nl
545  << "Relative tolerance : " << mergeTol << nl
546  << "Absolute matching distance : " << mergeDist << nl
547  << endl;
548 
549  // check writing tolerance
550  if (mesh.time().writeFormat() == IOstream::ASCII)
551  {
552  const scalar writeTol = std::pow
553  (
554  scalar(10.0),
555  -scalar(IOstream::defaultPrecision())
556  );
557 
558  if (mergeTol < writeTol)
559  {
561  << "Your current settings specify ASCII writing with "
562  << IOstream::defaultPrecision() << " digits precision." << nl
563  << "Your merging tolerance (" << mergeTol
564  << ") is finer than this." << nl
565  << "Change to binary writeFormat, "
566  << "or increase the writePrecision" << endl
567  << "or adjust the merge tolerance (mergeTol)."
568  << exit(FatalError);
569  }
570  }
571 
572  return mergeDist;
573 }
574 
575 
576 void removeZeroSizedPatches(fvMesh& mesh)
577 {
578  // Remove any zero-sized ones. Assumes
579  // - processor patches are already only there if needed
580  // - all other patches are available on all processors
581  // - but coupled ones might still be needed, even if zero-size
582  // (e.g. processorCyclic)
583  // See also logic in createPatch.
584  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
585 
586  labelList oldToNew(pbm.size(), -1);
587  label newPatchi = 0;
588  forAll(pbm, patchi)
589  {
590  const polyPatch& pp = pbm[patchi];
591 
592  if (!isA<processorPolyPatch>(pp))
593  {
594  if
595  (
596  isA<coupledPolyPatch>(pp)
597  || returnReduce(pp.size(), sumOp<label>())
598  )
599  {
600  // Coupled (and unknown size) or uncoupled and used
601  oldToNew[patchi] = newPatchi++;
602  }
603  }
604  }
605 
606  forAll(pbm, patchi)
607  {
608  const polyPatch& pp = pbm[patchi];
609 
610  if (isA<processorPolyPatch>(pp))
611  {
612  oldToNew[patchi] = newPatchi++;
613  }
614  }
615 
616 
617  const label nKeepPatches = newPatchi;
618 
619  // Shuffle unused ones to end
620  if (nKeepPatches != pbm.size())
621  {
622  Info<< endl
623  << "Removing zero-sized patches:" << endl << incrIndent;
624 
625  forAll(oldToNew, patchi)
626  {
627  if (oldToNew[patchi] == -1)
628  {
629  Info<< indent << pbm[patchi].name()
630  << " type " << pbm[patchi].type()
631  << " at position " << patchi << endl;
632  oldToNew[patchi] = newPatchi++;
633  }
634  }
635  Info<< decrIndent;
636 
637  fvMeshTools::reorderPatches(mesh, oldToNew, nKeepPatches, true);
638  Info<< endl;
639  }
640 }
641 
642 
643 // Write mesh and additional information
644 void writeMesh
645 (
646  const string& msg,
647  const meshRefinement& meshRefiner,
648  const meshRefinement::debugType debugLevel,
649  const meshRefinement::writeType writeLevel
650 )
651 {
652  const fvMesh& mesh = meshRefiner.mesh();
653 
654  meshRefiner.printMeshInfo(debugLevel, msg);
655  Info<< "Writing mesh to time " << meshRefiner.timeName() << endl;
656 
657  meshRefiner.write
658  (
659  debugLevel,
661  mesh.time().path()/meshRefiner.timeName()
662  );
663  Info<< "Wrote mesh in = "
664  << mesh.time().cpuTimeIncrement() << " s." << endl;
665 }
666 
667 
668 int main(int argc, char *argv[])
669 {
670  #include "addOverwriteOption.H"
672  (
673  "checkGeometry",
674  "check all surface geometry for quality"
675  );
677  (
678  "surfaceSimplify",
679  "boundBox",
680  "simplify the surface using snappyHexMesh starting from a boundBox"
681  );
683  (
684  "patches",
685  "(patch0 .. patchN)",
686  "only triangulate selected patches (wildcards supported)"
687  );
689  (
690  "outFile",
691  "file",
692  "name of the file to save the simplified surface to"
693  );
694  #include "addDictOption.H"
695 
696  #include "setRootCase.H"
697  #include "createTime.H"
698  runTime.functionObjects().off();
699 
700  const bool overwrite = args.optionFound("overwrite");
701  const bool checkGeometry = args.optionFound("checkGeometry");
702  const bool surfaceSimplify = args.optionFound("surfaceSimplify");
703 
705 
706  {
707  Foam::Info
708  << "Create mesh for time = "
709  << runTime.timeName() << Foam::nl << Foam::endl;
710 
711  meshPtr.set
712  (
713  new fvMesh
714  (
716  (
718  runTime.timeName(),
719  runTime,
721  )
722  )
723  );
724  }
725 
726  fvMesh& mesh = meshPtr();
727 
728  Info<< "Read mesh in = "
729  << runTime.cpuTimeIncrement() << " s" << endl;
730 
731  // Check patches and faceZones are synchronised
732  mesh.boundaryMesh().checkParallelSync(true);
734 
735 
736  // Read meshing dictionary
737  const word dictName("snappyHexMeshDict");
738  #include "setSystemMeshDictionaryIO.H"
739  const IOdictionary meshDict(dictIO);
740 
741 
742  // all surface geometry
743  const dictionary& geometryDict = meshDict.subDict("geometry");
744 
745  // refinement parameters
746  const dictionary& refineDict = meshDict.subDict("castellatedMeshControls");
747 
748  // mesh motion and mesh quality parameters
749  const dictionary& motionDict = meshDict.subDict("meshQualityControls");
750 
751  // snap-to-surface parameters
752  const dictionary& snapDict = meshDict.subDict("snapControls");
753 
754  // layer addition parameters
755  const dictionary& layerDict = meshDict.subDict("addLayersControls");
756 
757  // absolute merge distance
758  const scalar mergeDist = getMergeDistance
759  (
760  mesh,
761  readScalar(meshDict.lookup("mergeTolerance"))
762  );
763 
764  const Switch keepPatches(meshDict.lookupOrDefault("keepPatches", false));
765 
766 
767 
768  // Read decomposePar dictionary
769  dictionary decomposeDict;
770  {
771  if (Pstream::parRun())
772  {
773  decomposeDict = IOdictionary
774  (
775  IOobject
776  (
777  "decomposeParDict",
778  runTime.system(),
779  mesh,
782  )
783  );
784  }
785  else
786  {
787  decomposeDict.add("method", "none");
788  decomposeDict.add("numberOfSubdomains", 1);
789  }
790  }
791 
792 
793  // Debug
794  // ~~~~~
795 
796  // Set debug level
798  (
799  meshDict.lookupOrDefault<label>
800  (
801  "debug",
802  0
803  )
804  );
805  {
806  wordList flags;
807  if (meshDict.readIfPresent("debugFlags", flags))
808  {
809  debugLevel = meshRefinement::debugType
810  (
812  (
814  flags
815  )
816  );
817  }
818  }
819  if (debugLevel > 0)
820  {
821  meshRefinement::debug = debugLevel;
822  snappyRefineDriver::debug = debugLevel;
823  snappySnapDriver::debug = debugLevel;
824  snappyLayerDriver::debug = debugLevel;
825  }
826 
827  // Set file writing level
828  {
829  wordList flags;
830  if (meshDict.readIfPresent("writeFlags", flags))
831  {
833  (
835  (
837  (
839  flags
840  )
841  )
842  );
843  }
844  }
845 
846  // Set output level
847  {
848  wordList flags;
849  if (meshDict.readIfPresent("outputFlags", flags))
850  {
852  (
854  (
856  (
858  flags
859  )
860  )
861  );
862  }
863  }
864 
865 
866  // Read geometry
867  // ~~~~~~~~~~~~~
868 
869  searchableSurfaces allGeometry
870  (
871  IOobject
872  (
873  "abc", // dummy name
874  mesh.time().constant(), // instance
875  // mesh.time().findInstance("triSurface", word::null),// instance
876  "triSurface", // local
877  mesh.time(), // registry
880  ),
881  geometryDict,
882  meshDict.lookupOrDefault("singleRegionName", true)
883  );
884 
885 
886  // Read refinement surfaces
887  // ~~~~~~~~~~~~~~~~~~~~~~~~
888 
889  autoPtr<refinementSurfaces> surfacesPtr;
890 
891  Info<< "Reading refinement surfaces." << endl;
892 
893  if (surfaceSimplify)
894  {
895  IOdictionary foamyHexMeshDict
896  (
897  IOobject
898  (
899  "foamyHexMeshDict",
900  runTime.system(),
901  runTime,
904  )
905  );
906 
907  const dictionary& conformationDict =
908  foamyHexMeshDict.subDict("surfaceConformation").subDict
909  (
910  "geometryToConformTo"
911  );
912 
913  const dictionary& motionDict =
914  foamyHexMeshDict.subDict("motionControl");
915 
916  const dictionary& shapeControlDict =
917  motionDict.subDict("shapeControlFunctions");
918 
919  // Calculate current ratio of hex cells v.s. wanted cell size
920  const scalar defaultCellSize =
921  readScalar(motionDict.lookup("defaultCellSize"));
922 
923  const scalar initialCellSize = ::pow(meshPtr().V()[0], 1.0/3.0);
924 
925  // Info<< "Wanted cell size = " << defaultCellSize << endl;
926  // Info<< "Current cell size = " << initialCellSize << endl;
927  // Info<< "Fraction = " << initialCellSize/defaultCellSize
928  // << endl;
929 
930  surfacesPtr =
931  createRefinementSurfaces
932  (
933  allGeometry,
934  conformationDict,
935  shapeControlDict,
936  refineDict.lookupOrDefault("gapLevelIncrement", 0),
937  initialCellSize/defaultCellSize
938  );
939  }
940  else
941  {
942  surfacesPtr.set
943  (
945  (
946  allGeometry,
947  refineDict.subDict("refinementSurfaces"),
948  refineDict.lookupOrDefault("gapLevelIncrement", 0)
949  )
950  );
951 
952  Info<< "Read refinement surfaces in = "
953  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
954  }
955 
956  refinementSurfaces& surfaces = surfacesPtr();
957 
958 
959  // Checking only?
960 
961  if (checkGeometry)
962  {
963  // Extract patchInfo
964  List<wordList> patchTypes(allGeometry.size());
965 
966  const PtrList<dictionary>& patchInfo = surfaces.patchInfo();
967  const labelList& surfaceGeometry = surfaces.surfaces();
968  forAll(surfaceGeometry, surfI)
969  {
970  label geomI = surfaceGeometry[surfI];
971  const wordList& regNames = allGeometry.regionNames()[geomI];
972 
973  patchTypes[geomI].setSize(regNames.size());
974  forAll(regNames, regionI)
975  {
976  label globalRegionI = surfaces.globalRegion(surfI, regionI);
977 
978  if (patchInfo.set(globalRegionI))
979  {
980  patchTypes[geomI][regionI] =
981  word(patchInfo[globalRegionI].lookup("type"));
982  }
983  else
984  {
985  patchTypes[geomI][regionI] = wallPolyPatch::typeName;
986  }
987  }
988  }
989 
990  // Write some stats
991  allGeometry.writeStats(patchTypes, Info);
992  // Check topology
993  allGeometry.checkTopology(true);
994  // Check geometry
995  allGeometry.checkGeometry
996  (
997  100.0, // max size ratio
998  1e-9, // intersection tolerance
1000  0.01, // min triangle quality
1001  true
1002  );
1003 
1004  return 0;
1005  }
1006 
1007 
1008 
1009  // Read refinement shells
1010  // ~~~~~~~~~~~~~~~~~~~~~~
1011 
1012  Info<< "Reading refinement shells." << endl;
1013  shellSurfaces shells
1014  (
1015  allGeometry,
1016  refineDict.subDict("refinementRegions")
1017  );
1018  Info<< "Read refinement shells in = "
1019  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
1020 
1021 
1022  Info<< "Setting refinement level of surface to be consistent"
1023  << " with shells." << endl;
1024  surfaces.setMinLevelFields(shells);
1025  Info<< "Checked shell refinement in = "
1026  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
1027 
1028 
1029  // Read feature meshes
1030  // ~~~~~~~~~~~~~~~~~~~
1031 
1032  Info<< "Reading features." << endl;
1033  refinementFeatures features
1034  (
1035  mesh,
1036  refineDict.lookup("features")
1037  );
1038  Info<< "Read features in = "
1039  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
1040 
1041 
1042 
1043  // Refinement engine
1044  // ~~~~~~~~~~~~~~~~~
1045 
1046  Info<< nl
1047  << "Determining initial surface intersections" << nl
1048  << "-----------------------------------------" << nl
1049  << endl;
1050 
1051  // Main refinement engine
1052  meshRefinement meshRefiner
1053  (
1054  mesh,
1055  mergeDist, // tolerance used in sorting coordinates
1056  overwrite, // overwrite mesh files?
1057  surfaces, // for surface intersection refinement
1058  features, // for feature edges/point based refinement
1059  shells // for volume (inside/outside) refinement
1060  );
1061  Info<< "Calculated surface intersections in = "
1062  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
1063 
1064  // Some stats
1065  meshRefiner.printMeshInfo(debugLevel, "Initial mesh");
1066 
1067  meshRefiner.write
1068  (
1069  meshRefinement::debugType(debugLevel&meshRefinement::OBJINTERSECTIONS),
1071  mesh.time().path()/meshRefiner.timeName()
1072  );
1073 
1074 
1075  // Add all the surface regions as patches
1076  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1077 
1078  //- Global surface region to patch (non faceZone surface) or patches
1079  // (faceZone surfaces)
1080  labelList globalToMasterPatch;
1081  labelList globalToSlavePatch;
1082  {
1083  Info<< nl
1084  << "Adding patches for surface regions" << nl
1085  << "----------------------------------" << nl
1086  << endl;
1087 
1088  // From global region number to mesh patch.
1089  globalToMasterPatch.setSize(surfaces.nRegions(), -1);
1090  globalToSlavePatch.setSize(surfaces.nRegions(), -1);
1091 
1092  Info<< setf(ios_base::left)
1093  << setw(6) << "Patch"
1094  << setw(20) << "Type"
1095  << setw(30) << "Region" << nl
1096  << setw(6) << "-----"
1097  << setw(20) << "----"
1098  << setw(30) << "------" << endl;
1099 
1100  const labelList& surfaceGeometry = surfaces.surfaces();
1101  const PtrList<dictionary>& surfacePatchInfo = surfaces.patchInfo();
1102 
1103  forAll(surfaceGeometry, surfI)
1104  {
1105  label geomI = surfaceGeometry[surfI];
1106 
1107  const wordList& regNames = allGeometry.regionNames()[geomI];
1108 
1109  Info<< surfaces.names()[surfI] << ':' << nl << nl;
1110 
1111  if (surfaces.surfZones()[surfI].faceZoneName().empty())
1112  {
1113  // 'Normal' surface
1114  forAll(regNames, i)
1115  {
1116  label globalRegionI = surfaces.globalRegion(surfI, i);
1117 
1118  label patchi;
1119 
1120  if (surfacePatchInfo.set(globalRegionI))
1121  {
1122  patchi = meshRefiner.addMeshedPatch
1123  (
1124  regNames[i],
1125  surfacePatchInfo[globalRegionI]
1126  );
1127  }
1128  else
1129  {
1130  dictionary patchInfo;
1131  patchInfo.set("type", wallPolyPatch::typeName);
1132 
1133  patchi = meshRefiner.addMeshedPatch
1134  (
1135  regNames[i],
1136  patchInfo
1137  );
1138  }
1139 
1140  Info<< setf(ios_base::left)
1141  << setw(6) << patchi
1142  << setw(20) << mesh.boundaryMesh()[patchi].type()
1143  << setw(30) << regNames[i] << nl;
1144 
1145  globalToMasterPatch[globalRegionI] = patchi;
1146  globalToSlavePatch[globalRegionI] = patchi;
1147  }
1148  }
1149  else
1150  {
1151  // Zoned surface
1152  forAll(regNames, i)
1153  {
1154  label globalRegionI = surfaces.globalRegion(surfI, i);
1155 
1156  // Add master side patch
1157  {
1158  label patchi;
1159 
1160  if (surfacePatchInfo.set(globalRegionI))
1161  {
1162  patchi = meshRefiner.addMeshedPatch
1163  (
1164  regNames[i],
1165  surfacePatchInfo[globalRegionI]
1166  );
1167  }
1168  else
1169  {
1170  dictionary patchInfo;
1171  patchInfo.set("type", wallPolyPatch::typeName);
1172 
1173  patchi = meshRefiner.addMeshedPatch
1174  (
1175  regNames[i],
1176  patchInfo
1177  );
1178  }
1179 
1180  Info<< setf(ios_base::left)
1181  << setw(6) << patchi
1182  << setw(20) << mesh.boundaryMesh()[patchi].type()
1183  << setw(30) << regNames[i] << nl;
1184 
1185  globalToMasterPatch[globalRegionI] = patchi;
1186  }
1187  // Add slave side patch
1188  {
1189  const word slaveName = regNames[i] + "_slave";
1190  label patchi;
1191 
1192  if (surfacePatchInfo.set(globalRegionI))
1193  {
1194  patchi = meshRefiner.addMeshedPatch
1195  (
1196  slaveName,
1197  surfacePatchInfo[globalRegionI]
1198  );
1199  }
1200  else
1201  {
1202  dictionary patchInfo;
1203  patchInfo.set("type", wallPolyPatch::typeName);
1204 
1205  patchi = meshRefiner.addMeshedPatch
1206  (
1207  slaveName,
1208  patchInfo
1209  );
1210  }
1211 
1212  Info<< setf(ios_base::left)
1213  << setw(6) << patchi
1214  << setw(20) << mesh.boundaryMesh()[patchi].type()
1215  << setw(30) << slaveName << nl;
1216 
1217  globalToSlavePatch[globalRegionI] = patchi;
1218  }
1219  }
1220  }
1221 
1222  Info<< nl;
1223  }
1224  Info<< "Added patches in = "
1225  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
1226  }
1227 
1228 
1229  // Parallel
1230  // ~~~~~~~~
1231 
1232  // Decomposition
1233  autoPtr<decompositionMethod> decomposerPtr
1234  (
1236  (
1237  decomposeDict
1238  )
1239  );
1240  decompositionMethod& decomposer = decomposerPtr();
1241 
1242  if (Pstream::parRun() && !decomposer.parallelAware())
1243  {
1245  << "You have selected decomposition method "
1246  << decomposer.typeName
1247  << " which is not parallel aware." << endl
1248  << "Please select one that is (hierarchical, ptscotch)"
1249  << exit(FatalError);
1250  }
1251 
1252  // Mesh distribution engine (uses tolerance to reconstruct meshes)
1253  fvMeshDistribute distributor(mesh, mergeDist);
1254 
1255 
1256 
1257 
1258 
1259  // Now do the real work -refinement -snapping -layers
1260  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1261 
1262  const Switch wantRefine(meshDict.lookup("castellatedMesh"));
1263  const Switch wantSnap(meshDict.lookup("snap"));
1264  const Switch wantLayers(meshDict.lookup("addLayers"));
1265 
1266  // Refinement parameters
1267  const refinementParameters refineParams(refineDict);
1268 
1269  // Snap parameters
1270  const snapParameters snapParams(snapDict);
1271 
1272  // Layer addition parameters
1273  const layerParameters layerParams(layerDict, mesh.boundaryMesh());
1274 
1275 
1276  if (wantRefine)
1277  {
1278  cpuTime timer;
1279 
1280  snappyRefineDriver refineDriver
1281  (
1282  meshRefiner,
1283  decomposer,
1284  distributor,
1285  globalToMasterPatch,
1286  globalToSlavePatch
1287  );
1288 
1289 
1290  if (!overwrite && !debugLevel)
1291  {
1292  const_cast<Time&>(mesh.time())++;
1293  }
1294 
1295 
1296  refineDriver.doRefine
1297  (
1298  refineDict,
1299  refineParams,
1300  snapParams,
1301  refineParams.handleSnapProblems(),
1302  motionDict
1303  );
1304 
1305 
1306  if (!keepPatches && !wantSnap && !wantLayers)
1307  {
1308  removeZeroSizedPatches(mesh);
1309  }
1310 
1311  writeMesh
1312  (
1313  "Refined mesh",
1314  meshRefiner,
1315  debugLevel,
1317  );
1318 
1319  Info<< "Mesh refined in = "
1320  << timer.cpuTimeIncrement() << " s." << endl;
1321  }
1322 
1323  if (wantSnap)
1324  {
1325  cpuTime timer;
1326 
1327  snappySnapDriver snapDriver
1328  (
1329  meshRefiner,
1330  globalToMasterPatch,
1331  globalToSlavePatch
1332  );
1333 
1334  if (!overwrite && !debugLevel)
1335  {
1336  const_cast<Time&>(mesh.time())++;
1337  }
1338 
1339  // Use the resolveFeatureAngle from the refinement parameters
1340  scalar curvature = refineParams.curvature();
1341  scalar planarAngle = refineParams.planarAngle();
1342 
1343  snapDriver.doSnap
1344  (
1345  snapDict,
1346  motionDict,
1347  curvature,
1348  planarAngle,
1349  snapParams
1350  );
1351 
1352  if (!keepPatches && !wantLayers)
1353  {
1354  removeZeroSizedPatches(mesh);
1355  }
1356 
1357  writeMesh
1358  (
1359  "Snapped mesh",
1360  meshRefiner,
1361  debugLevel,
1363  );
1364 
1365  Info<< "Mesh snapped in = "
1366  << timer.cpuTimeIncrement() << " s." << endl;
1367  }
1368 
1369  if (wantLayers)
1370  {
1371  cpuTime timer;
1372 
1373  snappyLayerDriver layerDriver
1374  (
1375  meshRefiner,
1376  globalToMasterPatch,
1377  globalToSlavePatch
1378  );
1379 
1380  // Use the maxLocalCells from the refinement parameters
1381  bool preBalance = returnReduce
1382  (
1383  (mesh.nCells() >= refineParams.maxLocalCells()),
1384  orOp<bool>()
1385  );
1386 
1387 
1388  if (!overwrite && !debugLevel)
1389  {
1390  const_cast<Time&>(mesh.time())++;
1391  }
1392 
1393  layerDriver.doLayers
1394  (
1395  layerDict,
1396  motionDict,
1397  layerParams,
1398  preBalance,
1399  decomposer,
1400  distributor
1401  );
1402 
1403  if (!keepPatches)
1404  {
1405  removeZeroSizedPatches(mesh);
1406  }
1407 
1408  writeMesh
1409  (
1410  "Layer mesh",
1411  meshRefiner,
1412  debugLevel,
1414  );
1415 
1416  Info<< "Layers added in = "
1417  << timer.cpuTimeIncrement() << " s." << endl;
1418  }
1419 
1420 
1421  {
1422  // Check final mesh
1423  Info<< "Checking final mesh ..." << endl;
1424  faceSet wrongFaces(mesh, "wrongFaces", mesh.nFaces()/100);
1425  motionSmoother::checkMesh(false, mesh, motionDict, wrongFaces);
1426  const label nErrors = returnReduce
1427  (
1428  wrongFaces.size(),
1429  sumOp<label>()
1430  );
1431 
1432  if (nErrors > 0)
1433  {
1434  Info<< "Finished meshing with " << nErrors << " illegal faces"
1435  << " (concave, zero area or negative cell pyramid volume)"
1436  << endl;
1437  wrongFaces.write();
1438  }
1439  else
1440  {
1441  Info<< "Finished meshing without any errors" << endl;
1442  }
1443  }
1444 
1445 
1446  if (surfaceSimplify)
1447  {
1448  const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
1449 
1450  labelHashSet includePatches(bMesh.size());
1451 
1452  if (args.optionFound("patches"))
1453  {
1454  includePatches = bMesh.patchSet
1455  (
1456  wordReList(args.optionLookup("patches")())
1457  );
1458  }
1459  else
1460  {
1461  forAll(bMesh, patchi)
1462  {
1463  const polyPatch& patch = bMesh[patchi];
1464 
1465  if (!isA<processorPolyPatch>(patch))
1466  {
1467  includePatches.insert(patchi);
1468  }
1469  }
1470  }
1471 
1472  fileName outFileName
1473  (
1475  (
1476  "outFile",
1477  "constant/triSurface/simplifiedSurface.stl"
1478  )
1479  );
1480 
1481  extractSurface
1482  (
1483  mesh,
1484  runTime,
1485  includePatches,
1486  outFileName
1487  );
1488 
1489  pointIOField cellCentres
1490  (
1491  IOobject
1492  (
1493  "internalCellCentres",
1494  runTime.timeName(),
1495  mesh,
1498  ),
1499  mesh.cellCentres()
1500  );
1501 
1502  cellCentres.write();
1503  }
1504 
1505 
1506  Info<< "Finished meshing in = "
1507  << runTime.elapsedCpuTime() << " s." << endl;
1508 
1509  Info<< "End\n" << endl;
1510 
1511  return 0;
1512 }
1513 
1514 
1515 // ************************************************************************* //
static bool checkMesh(const bool report, const polyMesh &mesh, const dictionary &dict, labelHashSet &wrongFaces)
Check mesh with mesh settings in dict. Collects incorrect faces.
A surface geometry mesh, in which the surface zone information is conveyed by the &#39;zoneId&#39; associated...
Definition: MeshedSurface.H:74
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
void printMeshInfo(const bool, const string &) const
Print some mesh stats.
A HashTable with keys but without contents.
Definition: HashSet.H:59
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:140
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:438
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
bool processorCase() const
Return true if this is a processor case.
Definition: TimePaths.H:90
fileName path() const
Return path.
Definition: Time.H:266
Simple container to keep together layer specific information.
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
const word & name() const
Return name.
const word & name() const
Return name.
Definition: IOobject.H:295
A class for handling file names.
Definition: fileName.H:79
const labelList & surfaces() const
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:226
A list of face labels.
Definition: faceSet.H:48
virtual bool parallelAware() const =0
Is method parallel aware (i.e. does it synchronize domains across.
const entry * lookupEntryPtr(const word &, bool recursive, bool patternMatch) const
Find and return an entry data stream pointer if present.
Definition: dictionary.C:477
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
An identifier for a surface zone on a meshed surface.
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
const keyType & keyword() const
Return keyword.
Definition: entry.H:123
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:124
void off()
Switch the function objects off.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
bool clean()
Cleanup file name.
Definition: fileName.C:81
label nFaces() const
static fvMesh * meshPtr
Definition: globalFoam.H:52
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:461
static const NamedEnum< IOoutputType, 1 > IOoutputTypeNames
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
wordList patchTypes(nPatches)
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
label nCells() const
engineTime & runTime
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:309
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static writeType writeLevel()
Get/set write level.
virtual const dictionary & dict() const =0
Return dictionary if this entry is a dictionary.
label findPatchID(const word &patchName) const
Find patch index given a name.
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
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.
const labelList & minLevel() const
From global region number to refinement level.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
Simple container to keep together refinement specific information.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:108
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.
bool erase(T *p)
Remove the specified element from the list and delete it.
Definition: ILList.C:101
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
autoPtr< dictionary > clone() const
Construct and return clone.
Definition: dictionary.C:347
wordList toc() const
Return the table of contents.
Definition: dictionary.C:783
const fvMesh & mesh() const
Reference to mesh.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
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:1131
static void mapCombineScatter(const List< commsStruct > &comms, Container &Values, const int tag, const label comm)
Scatter data. Reverse of combineGather.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:699
Encapsulates queries for features.
T optionLookupOrDefault(const word &opt, const T &deflt) const
Read a value from the named option if present.
Definition: argListI.H:237
A list of faces which address into the list of points.
Sends/receives parts of mesh+fvfields to neighbouring processors. Used in load balancing.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
stressControl lookup("compactNormalStress") >> compactNormalStress
T clone(const T &t)
Definition: List.H:54
dynamicFvMesh & mesh
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:766
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Simple container to keep together snap specific information.
double cpuTimeIncrement() const
Return CPU time (in seconds) since last call to cpuTimeIncrement()
Definition: cpuTime.C:74
const labelList & gapLevel() const
From global region number to small gap refinement level.
const List< wordList > & regionNames() const
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:90
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:127
static const NamedEnum< IOdebugType, 5 > IOdebugTypeNames
static void checkCoupledFaceZones(const polyMesh &)
Helper function: check that face zones are synced.
Container for searchableSurfaces.
Smanip< ios_base::fmtflags > setf(const ios_base::fmtflags flags)
Definition: IOmanip.H:164
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1394
const word dictName("particleTrackDict")
Abstract base class for decomposition.
Encapsulates queries for volume refinement (&#39;refine all cells within shell&#39;).
Definition: shellSurfaces.H:52
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
static void mapCombineGather(const List< commsStruct > &comms, Container &Values, const CombineOp &cop, const int tag, const label comm)
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1156
const PtrList< surfaceZonesInfo > & surfZones() const
const vectorField & cellCentres() const
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
word timeName() const
Replacement for Time::timeName() : return oldInstance (if.
const word & system() const
Return system name.
Definition: TimePaths.H:114
static autoPtr< decompositionMethod > New(const dictionary &decompositionDict)
Return a reference to the selected decomposition method.
Foam::polyBoundaryMesh.
Istream and Ostream manipulators taking arguments.
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order. Return.
static const char nl
Definition: Ostream.H:265
IOstream::streamFormat writeFormat() const
Default write format.
Definition: Time.H:288
const Time & time() const
Return time.
static int readFlags(const Enum &namedEnum, const wordList &)
Helper: convert wordList into bit pattern using provided.
label addMeshedPatch(const word &name, const dictionary &)
Add patch originating from meshing. Update meshedPatches_.
void set(T *)
Set pointer to that given.
Definition: autoPtrI.H:99
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:240
labelHashSet includePatches
All to do with adding layers.
static outputType outputLevel()
Get/set output level.
IOobject dictIO(dictName, runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE)
Helper class which maintains intersections of (changing) mesh with (static) surfaces.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
label checkGeometry(const scalar maxRatio, const scalar tolerance, const autoPtr< writer< scalar >> &setWriter, const scalar minQuality, const bool report) const
All geometric checks. Return number of failed checks.
All to do with snapping to surface.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void setSize(const label)
Reset size of List.
Definition: List.C:281
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:440
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
const wordList & names() const
Names of surfaces.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
label patchi
double elapsedCpuTime() const
Return CPU time (in seconds) from the start.
Definition: cpuTime.C:67
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const functionObjectList & functionObjects() const
Return the list of function objects.
Definition: Time.H:410
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
const wordList & names() const
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A List with indirect addressing.
Definition: fvMatrix.H:106
label index() const
Return the index of this patch in the boundaryMesh.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:303
label globalRegion(const label surfI, const label regionI) const
From surface and region on surface to global region.
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:948
messageStream Info
const PtrList< dictionary > & patchInfo() const
From global region number to patch type.
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
virtual Ostream & write(const token &)=0
Write next token to stream.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:117
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Starts timing CPU usage and return elapsed time from start.
Definition: cpuTime.H:52
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
bool write() const
Write mesh and all data.
virtual bool write(const bool write=true) const
Write using setting from DB.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
void setMinLevelFields(const shellSurfaces &shells)
Calculate minLevelFields.
Foam::argList args(argc, argv)
void writeStats(const List< wordList > &, Ostream &) const
Write some stats.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:233
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
Namespace for OpenFOAM.
A keyword and a list of tokens is an &#39;entry&#39;.
Definition: entry.H:65
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583
static const NamedEnum< IOwriteType, 5 > IOwriteTypeNames
const labelList & maxLevel() const
From global region number to refinement level.
IStringStream optionLookup(const word &opt) const
Return an IStringStream from the named option.
Definition: argListI.H:114