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