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