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