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-2022 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 non-constraint zero-sized patches
569 
570  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
571 
572  labelList oldToNew(pbm.size(), -1);
573  label newPatchi = 0;
574  forAll(pbm, patchi)
575  {
576  const polyPatch& pp = pbm[patchi];
577 
578  if
579  (
580  polyPatch::constraintType(pp.type())
581  || returnReduce(pp.size(), sumOp<label>())
582  )
583  {
584  oldToNew[patchi] = newPatchi++;
585  }
586  }
587 
588  const label nKeepPatches = newPatchi;
589 
590  // Shuffle unused ones to end
591  if (nKeepPatches != pbm.size())
592  {
593  Info<< endl
594  << "Removing zero-sized patches:" << endl << incrIndent;
595 
596  forAll(oldToNew, patchi)
597  {
598  if (oldToNew[patchi] == -1)
599  {
600  Info<< indent << pbm[patchi].name()
601  << " type " << pbm[patchi].type()
602  << " at position " << patchi << endl;
603  oldToNew[patchi] = newPatchi++;
604  }
605  }
606  Info<< decrIndent;
607 
608  fvMeshTools::reorderPatches(mesh, oldToNew, nKeepPatches, true);
609  Info<< endl;
610  }
611 }
612 
613 
614 // Write mesh and additional information
615 void writeMesh
616 (
617  const string& msg,
618  const meshRefinement& meshRefiner,
619  const meshRefinement::debugType debugLevel,
620  const meshRefinement::writeType writeLevel
621 )
622 {
623  const fvMesh& mesh = meshRefiner.mesh();
624 
625  meshRefiner.printMeshInfo(debugLevel, msg);
626  Info<< "Writing mesh to time " << meshRefiner.timeName() << endl;
627 
628  meshRefiner.write
629  (
630  debugLevel,
632  mesh.time().path()/meshRefiner.timeName()
633  );
634  Info<< "Wrote mesh in = "
635  << mesh.time().cpuTimeIncrement() << " s." << endl;
636 }
637 
638 
639 int main(int argc, char *argv[])
640 {
641  #include "addOverwriteOption.H"
643  (
644  "checkGeometry",
645  "check all surface geometry for quality"
646  );
648  (
649  "surfaceSimplify",
650  "boundBox",
651  "simplify the surface using snappyHexMesh starting from a boundBox"
652  );
654  (
655  "patches",
656  "(patch0 .. patchN)",
657  "only triangulate selected patches (wildcards supported)"
658  );
660  (
661  "outFile",
662  "file",
663  "name of the file to save the simplified surface to"
664  );
665  #include "addDictOption.H"
666 
667  #include "setRootCase.H"
668  #include "createTime.H"
669  runTime.functionObjects().off();
670 
671  const bool overwrite = args.optionFound("overwrite");
672  const bool checkGeometry = args.optionFound("checkGeometry");
673  const bool surfaceSimplify = args.optionFound("surfaceSimplify");
674 
676 
677  {
678  Foam::Info
679  << "Create mesh for time = "
680  << runTime.timeName() << Foam::nl << Foam::endl;
681 
682  meshPtr.set
683  (
684  new fvMesh
685  (
687  (
689  runTime.timeName(),
690  runTime,
692  ),
693  false
694  )
695  );
696  }
697 
698  fvMesh& mesh = meshPtr();
699 
700  Info<< "Read mesh in = "
701  << runTime.cpuTimeIncrement() << " s" << endl;
702 
703  // Check that the read mesh is fully 3D
704  // as required for mesh relaxation after snapping
705  if (mesh.nSolutionD() != 3)
706  {
708  << "Mesh provided is not fully 3D"
709  " as required for mesh relaxation after snapping" << nl
710  << "Convert all empty patches to appropriate types for a 3D mesh,"
711  " current patch types are" << nl << mesh.boundaryMesh().types()
712  << exit(FatalError);
713  }
714 
715  // Check patches and faceZones are synchronised
716  mesh.boundaryMesh().checkParallelSync(true);
718 
719 
720  // Read meshing dictionary
721  const IOdictionary meshDict(systemDict("snappyHexMeshDict", args, mesh));
722 
723 
724  // all surface geometry
725  const dictionary& geometryDict = meshDict.subDict("geometry");
726 
727  // refinement parameters
728  const dictionary& refineDict = meshDict.subDict("castellatedMeshControls");
729 
730  // mesh motion and mesh quality parameters
731  const dictionary& motionDict = meshDict.subDict("meshQualityControls");
732 
733  // snap-to-surface parameters
734  const dictionary& snapDict = meshDict.subDict("snapControls");
735 
736  // absolute merge distance
737  const scalar mergeDist = getMergeDistance
738  (
739  mesh,
740  meshDict.lookup<scalar>("mergeTolerance")
741  );
742 
743  const Switch keepPatches(meshDict.lookupOrDefault("keepPatches", false));
744 
745 
746  // Read decomposePar dictionary
747  dictionary decomposeDict;
748  {
749  if (Pstream::parRun())
750  {
751  decomposeDict = decompositionMethod::decomposeParDict(runTime);
752  }
753  else
754  {
755  decomposeDict.add("method", "none");
756  decomposeDict.add("numberOfSubdomains", 1);
757  }
758  }
759 
760 
761  // Debug
762  // ~~~~~
763 
764  // Set debug level
766  (
767  meshDict.lookupOrDefault<label>
768  (
769  "debug",
770  0
771  )
772  );
773  {
774  wordList flags;
775  if (meshDict.readIfPresent("debugFlags", flags))
776  {
777  debugLevel = meshRefinement::debugType
778  (
780  (
782  flags
783  )
784  );
785  }
786  }
787  if (debugLevel > 0)
788  {
789  meshRefinement::debug = debugLevel;
790  snappyRefineDriver::debug = debugLevel;
791  snappySnapDriver::debug = debugLevel;
792  snappyLayerDriver::debug = debugLevel;
793  }
794 
795  // Set file writing level
796  {
797  wordList flags;
798  if (meshDict.readIfPresent("writeFlags", flags))
799  {
801  (
803  (
805  (
807  flags
808  )
809  )
810  );
811  }
812  }
813 
814  // Set output level
815  {
816  wordList flags;
817  if (meshDict.readIfPresent("outputFlags", flags))
818  {
820  (
822  (
824  (
826  flags
827  )
828  )
829  );
830  }
831  }
832 
833 
834  // Read geometry
835  // ~~~~~~~~~~~~~
836 
837  searchableSurfaces allGeometry
838  (
839  IOobject
840  (
841  "abc",
842  mesh.time().constant(),
844  mesh.time(),
847  ),
848  geometryDict,
849  meshDict.lookupOrDefault("singleRegionName", true)
850  );
851 
852 
853  // Read refinement surfaces
854  // ~~~~~~~~~~~~~~~~~~~~~~~~
855 
856  autoPtr<refinementSurfaces> surfacesPtr;
857 
858  Info<< "Reading refinement surfaces..." << endl;
859 
860  if (surfaceSimplify)
861  {
862  IOdictionary foamyHexMeshDict
863  (
864  IOobject
865  (
866  "foamyHexMeshDict",
867  runTime.system(),
868  runTime,
871  )
872  );
873 
874  const dictionary& conformationDict =
875  foamyHexMeshDict.subDict("surfaceConformation").subDict
876  (
877  "geometryToConformTo"
878  );
879 
880  const dictionary& motionDict =
881  foamyHexMeshDict.subDict("motionControl");
882 
883  const dictionary& shapeControlDict =
884  motionDict.subDict("shapeControlFunctions");
885 
886  // Calculate current ratio of hex cells v.s. wanted cell size
887  const scalar defaultCellSize =
888  motionDict.lookup<scalar>("defaultCellSize");
889 
890  const scalar initialCellSize = ::pow(meshPtr().V()[0], 1.0/3.0);
891 
892  // Info<< "Wanted cell size = " << defaultCellSize << endl;
893  // Info<< "Current cell size = " << initialCellSize << endl;
894  // Info<< "Fraction = " << initialCellSize/defaultCellSize
895  // << endl;
896 
897  surfacesPtr =
898  createRefinementSurfaces
899  (
900  allGeometry,
901  conformationDict,
902  shapeControlDict,
903  refineDict.lookupOrDefault("gapLevelIncrement", 0),
904  initialCellSize/defaultCellSize
905  );
906  }
907  else
908  {
909  surfacesPtr.set
910  (
912  (
913  allGeometry,
914  refineDict.found("refinementSurfaces")
915  ? refineDict.subDict("refinementSurfaces")
917  refineDict.lookupOrDefault("gapLevelIncrement", 0)
918  )
919  );
920 
921  Info<< "Read refinement surfaces in = "
922  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
923  }
924 
925  refinementSurfaces& surfaces = surfacesPtr();
926 
927 
928  // Checking only?
929 
930  if (checkGeometry)
931  {
932  // Extract patchInfo
933  List<wordList> patchTypes(allGeometry.size());
934 
935  const PtrList<dictionary>& patchInfo = surfaces.patchInfo();
936  const labelList& surfaceGeometry = surfaces.surfaces();
937  forAll(surfaceGeometry, surfI)
938  {
939  label geomI = surfaceGeometry[surfI];
940  const wordList& regNames = allGeometry.regionNames()[geomI];
941 
942  patchTypes[geomI].setSize(regNames.size());
943  forAll(regNames, regionI)
944  {
945  label globalRegionI = surfaces.globalRegion(surfI, regionI);
946 
947  if (patchInfo.set(globalRegionI))
948  {
949  patchTypes[geomI][regionI] =
950  word(patchInfo[globalRegionI].lookup("type"));
951  }
952  else
953  {
954  patchTypes[geomI][regionI] = wallPolyPatch::typeName;
955  }
956  }
957  }
958 
959  // Write some stats
960  allGeometry.writeStats(patchTypes, Info);
961  // Check topology
962  allGeometry.checkTopology(true);
963  // Check geometry
964  allGeometry.checkGeometry
965  (
966  100.0, // max size ratio
967  1e-9, // intersection tolerance
968  0.01, // min triangle quality
969  true
970  );
971 
972  return 0;
973  }
974 
975 
976 
977  // Read refinement shells
978  // ~~~~~~~~~~~~~~~~~~~~~~
979 
980  Info<< "Reading refinement regions..." << endl;
981  refinementRegions shells
982  (
983  allGeometry,
984  refineDict.found("refinementRegions")
985  ? refineDict.subDict("refinementRegions")
987  );
988  Info<< "Read refinement regions in = "
989  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
990 
991 
992 
993  // Read feature meshes
994  // ~~~~~~~~~~~~~~~~~~~
995 
996  Info<< "Reading features..." << endl;
997  refinementFeatures features
998  (
999  mesh,
1000  refineDict.found("features")
1001  ? refineDict.lookup("features")
1003  );
1004  Info<< "Read features in = "
1005  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
1006 
1007 
1008 
1009  // Refinement engine
1010  // ~~~~~~~~~~~~~~~~~
1011 
1012  Info<< nl
1013  << "Determining initial surface intersections" << nl
1014  << "-----------------------------------------" << nl
1015  << endl;
1016 
1017  // Main refinement engine
1018  meshRefinement meshRefiner
1019  (
1020  mesh,
1021  refineDict,
1022  mergeDist, // tolerance used in sorting coordinates
1023  overwrite, // overwrite mesh files?
1024  surfaces, // for surface intersection refinement
1025  features, // for feature edges/point based refinement
1026  shells // for volume (inside/outside) refinement
1027  );
1028  Info<< "Calculated surface intersections in = "
1029  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
1030 
1031  // Some stats
1032  meshRefiner.printMeshInfo(debugLevel, "Initial mesh");
1033 
1034  meshRefiner.write
1035  (
1036  meshRefinement::debugType(debugLevel&meshRefinement::OBJINTERSECTIONS),
1038  mesh.time().path()/meshRefiner.timeName()
1039  );
1040 
1041 
1042  // Add all the surface regions as patches
1043  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1044 
1045  //- Global surface region to patch (non faceZone surface) or patches
1046  // (faceZone surfaces)
1047  labelList globalToMasterPatch;
1048  labelList globalToSlavePatch;
1049  {
1050  Info<< nl
1051  << "Adding patches for surface regions" << nl
1052  << "----------------------------------" << nl
1053  << endl;
1054 
1055  // From global region number to mesh patch.
1056  globalToMasterPatch.setSize(surfaces.nRegions(), -1);
1057  globalToSlavePatch.setSize(surfaces.nRegions(), -1);
1058 
1059  Info<< setf(ios_base::left)
1060  << setw(6) << "Patch"
1061  << setw(20) << "Type"
1062  << setw(30) << "Region" << nl
1063  << setw(6) << "-----"
1064  << setw(20) << "----"
1065  << setw(30) << "------" << endl;
1066 
1067  const labelList& surfaceGeometry = surfaces.surfaces();
1068  const PtrList<dictionary>& surfacePatchInfo = surfaces.patchInfo();
1069 
1070  forAll(surfaceGeometry, surfI)
1071  {
1072  label geomI = surfaceGeometry[surfI];
1073 
1074  const wordList& regNames = allGeometry.regionNames()[geomI];
1075 
1076  Info<< surfaces.names()[surfI] << ':' << nl << nl;
1077 
1078  if (surfaces.surfZones()[surfI].faceZoneName().empty())
1079  {
1080  // 'Normal' surface
1081  forAll(regNames, i)
1082  {
1083  label globalRegionI = surfaces.globalRegion(surfI, i);
1084 
1085  label patchi;
1086 
1087  if (surfacePatchInfo.set(globalRegionI))
1088  {
1089  patchi = meshRefiner.addMeshedPatch
1090  (
1091  regNames[i],
1092  surfacePatchInfo[globalRegionI]
1093  );
1094  }
1095  else
1096  {
1097  dictionary patchInfo;
1098  patchInfo.set("type", wallPolyPatch::typeName);
1099 
1100  patchi = meshRefiner.addMeshedPatch
1101  (
1102  regNames[i],
1103  patchInfo
1104  );
1105  }
1106 
1107  Info<< setf(ios_base::left)
1108  << setw(6) << patchi
1109  << setw(20) << mesh.boundaryMesh()[patchi].type()
1110  << setw(30) << regNames[i] << nl;
1111 
1112  globalToMasterPatch[globalRegionI] = patchi;
1113  globalToSlavePatch[globalRegionI] = patchi;
1114  }
1115  }
1116  else
1117  {
1118  // Zoned surface
1119  forAll(regNames, i)
1120  {
1121  label globalRegionI = surfaces.globalRegion(surfI, i);
1122 
1123  // Add master side patch
1124  {
1125  label patchi;
1126 
1127  if (surfacePatchInfo.set(globalRegionI))
1128  {
1129  patchi = meshRefiner.addMeshedPatch
1130  (
1131  regNames[i],
1132  surfacePatchInfo[globalRegionI]
1133  );
1134  }
1135  else
1136  {
1137  dictionary patchInfo;
1138  patchInfo.set("type", wallPolyPatch::typeName);
1139 
1140  patchi = meshRefiner.addMeshedPatch
1141  (
1142  regNames[i],
1143  patchInfo
1144  );
1145  }
1146 
1147  Info<< setf(ios_base::left)
1148  << setw(6) << patchi
1149  << setw(20) << mesh.boundaryMesh()[patchi].type()
1150  << setw(30) << regNames[i] << nl;
1151 
1152  globalToMasterPatch[globalRegionI] = patchi;
1153  }
1154  // Add slave side patch
1155  {
1156  const word slaveName = regNames[i] + "_slave";
1157  label patchi;
1158 
1159  if (surfacePatchInfo.set(globalRegionI))
1160  {
1161  patchi = meshRefiner.addMeshedPatch
1162  (
1163  slaveName,
1164  surfacePatchInfo[globalRegionI]
1165  );
1166  }
1167  else
1168  {
1169  dictionary patchInfo;
1170  patchInfo.set("type", wallPolyPatch::typeName);
1171 
1172  patchi = meshRefiner.addMeshedPatch
1173  (
1174  slaveName,
1175  patchInfo
1176  );
1177  }
1178 
1179  Info<< setf(ios_base::left)
1180  << setw(6) << patchi
1181  << setw(20) << mesh.boundaryMesh()[patchi].type()
1182  << setw(30) << slaveName << nl;
1183 
1184  globalToSlavePatch[globalRegionI] = patchi;
1185  }
1186  }
1187  }
1188 
1189  Info<< nl;
1190  }
1191  Info<< "Added patches in = "
1192  << mesh.time().cpuTimeIncrement() << " s" << nl << endl;
1193  }
1194 
1195 
1196  // Parallel
1197  // ~~~~~~~~
1198 
1199  // Decomposition
1200  autoPtr<decompositionMethod> decomposerPtr
1201  (
1203  );
1204  decompositionMethod& decomposer = decomposerPtr();
1205 
1206  // Mesh distribution engine (uses tolerance to reconstruct meshes)
1207  fvMeshDistribute distributor(mesh);
1208 
1209 
1210  // Now do the real work -refinement -snapping -layers
1211  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1212 
1213  const Switch wantRefine(meshDict.lookup("castellatedMesh"));
1214  const Switch wantSnap(meshDict.lookup("snap"));
1215  const Switch wantLayers(meshDict.lookup("addLayers"));
1216 
1217  // Refinement parameters
1218  const refinementParameters refineParams(refineDict);
1219 
1220  // Snap parameters
1221  const snapParameters snapParams(snapDict);
1222 
1223  if (wantRefine)
1224  {
1225  cpuTime timer;
1226 
1227  snappyRefineDriver refineDriver
1228  (
1229  meshRefiner,
1230  decomposer,
1231  distributor,
1232  globalToMasterPatch,
1233  globalToSlavePatch
1234  );
1235 
1236  if (!overwrite && !debugLevel)
1237  {
1238  const_cast<Time&>(mesh.time())++;
1239  }
1240 
1241  refineDriver.doRefine
1242  (
1243  refineDict,
1244  refineParams,
1245  snapParams,
1246  refineParams.handleSnapProblems(),
1247  motionDict
1248  );
1249 
1250  if (!keepPatches && !wantSnap && !wantLayers)
1251  {
1252  removeZeroSizedPatches(mesh);
1253  }
1254 
1255  writeMesh
1256  (
1257  "Refined mesh",
1258  meshRefiner,
1259  debugLevel,
1261  );
1262 
1263  Info<< "Mesh refined in = "
1264  << timer.cpuTimeIncrement() << " s." << endl;
1265  }
1266 
1267  if (wantSnap)
1268  {
1269  cpuTime timer;
1270 
1271  snappySnapDriver snapDriver
1272  (
1273  meshRefiner,
1274  globalToMasterPatch,
1275  globalToSlavePatch
1276  );
1277 
1278  if (!overwrite && !debugLevel)
1279  {
1280  const_cast<Time&>(mesh.time())++;
1281  }
1282 
1283  // Use the resolveFeatureAngle from the refinement parameters
1284  scalar curvature = refineParams.curvature();
1285  scalar planarAngle = refineParams.planarAngle();
1286 
1287  snapDriver.doSnap
1288  (
1289  snapDict,
1290  motionDict,
1291  curvature,
1292  planarAngle,
1293  snapParams
1294  );
1295 
1296  if (!keepPatches && !wantLayers)
1297  {
1298  removeZeroSizedPatches(mesh);
1299  }
1300 
1301  writeMesh
1302  (
1303  "Snapped mesh",
1304  meshRefiner,
1305  debugLevel,
1307  );
1308 
1309  Info<< "Mesh snapped in = "
1310  << timer.cpuTimeIncrement() << " s." << endl;
1311  }
1312 
1313  if (wantLayers)
1314  {
1315  cpuTime timer;
1316 
1317  // Layer addition parameters dictionary
1318  const dictionary& layersDict = meshDict.subDict("addLayersControls");
1319 
1320  // Layer addition parameters
1321  const layerParameters layerParams(layersDict, mesh.boundaryMesh());
1322 
1323  snappyLayerDriver layerDriver
1324  (
1325  meshRefiner,
1326  globalToMasterPatch,
1327  globalToSlavePatch
1328  );
1329 
1330  // Use the maxLocalCells from the refinement parameters
1331  bool preBalance = returnReduce
1332  (
1333  (mesh.nCells() >= refineParams.maxLocalCells()),
1334  orOp<bool>()
1335  );
1336 
1337 
1338  if (!overwrite && !debugLevel)
1339  {
1340  const_cast<Time&>(mesh.time())++;
1341  }
1342 
1343  layerDriver.doLayers
1344  (
1345  layersDict,
1346  motionDict,
1347  layerParams,
1348  preBalance,
1349  decomposer,
1350  distributor
1351  );
1352 
1353  if (!keepPatches)
1354  {
1355  removeZeroSizedPatches(mesh);
1356  }
1357 
1358  writeMesh
1359  (
1360  "Layer mesh",
1361  meshRefiner,
1362  debugLevel,
1364  );
1365 
1366  Info<< "Layers added in = "
1367  << timer.cpuTimeIncrement() << " s." << endl;
1368  }
1369 
1370 
1371  {
1372  // Check final mesh
1373  Info<< "Checking final mesh ..." << endl;
1374  faceSet wrongFaces(mesh, "wrongFaces", mesh.nFaces()/100);
1375  motionSmoother::checkMesh(false, mesh, motionDict, wrongFaces);
1376  const label nErrors = returnReduce
1377  (
1378  wrongFaces.size(),
1379  sumOp<label>()
1380  );
1381 
1382  if (nErrors > 0)
1383  {
1384  Info<< "Finished meshing with " << nErrors << " illegal faces"
1385  << " (concave, zero area or negative cell pyramid volume)"
1386  << endl;
1387  wrongFaces.write();
1388  }
1389  else
1390  {
1391  Info<< "Finished meshing without any errors" << endl;
1392  }
1393  }
1394 
1395 
1396  if (surfaceSimplify)
1397  {
1398  const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
1399 
1400  labelHashSet includePatches(bMesh.size());
1401 
1402  if (args.optionFound("patches"))
1403  {
1404  includePatches = bMesh.patchSet
1405  (
1406  wordReList(args.optionLookup("patches")())
1407  );
1408  }
1409  else
1410  {
1411  forAll(bMesh, patchi)
1412  {
1413  const polyPatch& patch = bMesh[patchi];
1414 
1415  if (!isA<processorPolyPatch>(patch))
1416  {
1417  includePatches.insert(patchi);
1418  }
1419  }
1420  }
1421 
1422  fileName outFileName
1423  (
1425  (
1426  "outFile",
1427  "constant/triSurface/simplifiedSurface.stl"
1428  )
1429  );
1430 
1431  extractSurface
1432  (
1433  mesh,
1434  runTime,
1435  includePatches,
1436  outFileName
1437  );
1438 
1439  pointIOField cellCentres
1440  (
1441  IOobject
1442  (
1443  "internalCellCentres",
1444  runTime.timeName(),
1445  mesh,
1448  ),
1449  mesh.cellCentres()
1450  );
1451 
1452  cellCentres.write();
1453  }
1454 
1455 
1456  Info<< "Finished meshing in = "
1457  << runTime.elapsedCpuTime() << " s." << endl;
1458 
1459  Info<< "End\n" << endl;
1460 
1461  return 0;
1462 }
1463 
1464 
1465 // ************************************************************************* //
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:453
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:141
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
const word & executable() const
Name of executable without the path.
Definition: argListI.H:36
#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.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
virtual Ostream & write(const char)=0
Write character.
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const word & name() const
Return name.
const word & name() const
Return name.
Definition: IOobject.H:315
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
const entry * lookupEntryPtr(const word &, bool recursive, bool patternMatch) const
Find and return an entry data stream pointer if present.
Definition: dictionary.C:702
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
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
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
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:325
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
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.
wordList types() const
Return a list of patch types.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
Simple container to keep together refinement specific information.
fvMesh & mesh
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:69
autoPtr< dictionary > clone() const
Construct and return clone.
Definition: dictionary.C:552
wordList toc() const
Return the table of contents.
Definition: dictionary.C:1115
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:1211
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:1002
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: polyPatch.C:265
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
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:1060
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
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:60
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
static autoPtr< decompositionMethod > NewDistributor(const dictionary &decompositionDict)
Return a reference to the selected decomposition method.
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:1435
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:1224
An STL-conforming hash table.
Definition: HashTable.H:61
const PtrList< surfaceZonesInfo > & surfZones() const
const vectorField & cellCentres() const
static dictionary decomposeParDict(const Time &time)
Read and return the decomposeParDict.
word timeName() const
Replacement for Time::timeName() : return oldInstance (if.
const word & system() const
Return system name.
Definition: TimePaths.H:113
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:293
label globalRegion(const label surfi, const label regioni) const
From surface and region on surface to global region.
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
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:922
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:459
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:53
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:433
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:95
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:306
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:1291
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:76
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:54
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:98
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:301
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:864
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