47#ifndef Foam_meshRefinement_H
48#define Foam_meshRefinement_H
174 const scalar mergeDistance_;
177 const bool overwrite_;
180 const word oldInstance_;
235 template<
class GeoField>
236 static void addPatchFields(
fvMesh&,
const word& patchFieldType);
239 template<
class GeoField>
252 label globalFaceCount(
const labelList& elems)
const;
259 void calcCellCellRays
278 const label exposedPatchi
287 static bool markForRefine
289 const label markValue,
290 const label nAllowRefine,
297 void markFeatureCellLevel
305 label markFeatureRefinement
308 const label nAllowRefine,
315 label markInternalDistanceToFeatureRefinement
317 const label nAllowRefine,
323 label markInternalRefinement
325 const label nAllowRefine,
332 label unmarkInternalRefinement
346 label markSurfaceRefinement
348 const label nAllowRefine,
357 void collectGapCandidates
370 const scalar planarCos,
383 const label nAllowRefine,
394 label markSurfaceGapRefinement
396 const scalar planarCos,
397 const label nAllowRefine,
408 const point& nearPoint,
423 const bool useSurfaceNormal,
425 const point& nearPoint,
444 void selectGapCandidates
468 label markInternalGapRefinement
470 const scalar planarCos,
471 const bool spreadGapSize,
472 const label nAllowRefine,
480 label markSmallFeatureRefinement
482 const scalar planarCos,
483 const label nAllowRefine,
494 label markSurfaceFieldRefinement
496 const label nAllowRefine,
509 const scalar tol = 1
e-6
516 const scalar minCosAngle,
517 const scalar lengthScale,
527 label markSurfaceCurvatureRefinement
529 const scalar curvature,
530 const label nAllowRefine,
540 const scalar planarCos,
541 const label nAllowRefine,
543 const label surfaceLevel,
545 const vector& surfaceNormal,
558 label markProximityRefinement
560 const scalar curvature,
566 const label nAllowRefine,
576 void markProximityCandidateFaces
588 label markFakeGapRefinement
590 const scalar planarCos,
592 const label nAllowRefine,
606 const bool allowBoundary,
613 void getIntersections
626 void getIntersections
637 void getBafflePatches
639 const label nErodeCellZones,
644 const bool exitIfLeakPath,
656 const label nBufferLayers,
670 const label ownPatch,
671 const label neiPatch,
690 void markBoundaryFace
717 const scalar minFaceArea,
718 const scalar maxNonOrtho,
725 const scalar volFraction,
732 void markFacesOnProblemCells
735 const bool removeEdgeConnectedCells,
776 const label defaultRegion
781 void markFacesOnProblemCellsGeometric
801 const scalar freeStandingAngle,
812 void findCellZoneGeometric
823 void findCellZoneInsideWalk
833 void findCellZoneInsideWalk
842 bool calcRegionToZone
844 const label backgroundZoneID,
845 const label surfZoneI,
846 const label ownRegion,
847 const label neiRegion,
856 void findCellZoneTopo
858 const label backgroundZoneID,
870 const label nErodeCellZones,
871 const label backgroundZoneID,
881 const label nGrowCellZones,
882 const label backgroundZoneID,
891 void makeConsistentFaceIndex
901 const bool allowFreeStandingZoneFaces,
902 const label nErodeCellZones,
903 const label backgroundZoneID,
907 const bool exitIfLeakPath,
921 const bitSet& isMasterFace,
925 const bitSet& meshFlipMap,
930 void allocateInterRegionFaceZone
939 void handleSnapProblems
942 const bool useTopologicalSnapDetection,
943 const bool removeEdgeConnectedCells,
965 void calcPatchNumMasterFaces
967 const bitSet& isMasterFace,
982 void consistentOrientation
984 const bitSet& isMasterFace,
994 meshRefinement(
const meshRefinement&) =
delete;
997 void operator=(
const meshRefinement&) =
delete;
1039 return mergeDistance_;
1051 return oldInstance_;
1055 const refinementSurfaces&
surfaces()
const
1061 const refinementFeatures&
features()
const
1067 const shellSurfaces&
shells()
const
1075 return limitShells_;
1099 return faceToCoupledPatch_;
1106 const List<Tuple2<mapType, labelList>>&
userFaceData()
const
1108 return userFaceData_;
1113 return userFaceData_;
1127 autoPtr<mapDistributePolyMesh>
balance
1129 const bool keepZoneFaces,
1130 const bool keepBaffles,
1133 decompositionMethod& decomposer,
1134 fvMeshDistribute& distributor
1144 static autoPtr<indirectPrimitivePatch>
makePatch
1158 const pointMesh& pMesh,
1168 const polyMesh&
mesh,
1169 const bitSet& isMasterEdge,
1178 template<
class Type>
1181 const polyMesh&
mesh,
1182 const bitSet& isMasterEdge,
1186 const Field<Type>& data,
1206 const scalar planarCos,
1219 const scalar curvature,
1220 const scalar planarAngle,
1222 const bool featureRefinement,
1223 const bool featureDistanceRefinement,
1224 const bool internalRefinement,
1225 const bool surfaceRefinement,
1226 const bool curvatureRefinement,
1227 const bool smallFeatureRefinement,
1228 const bool gapRefinement,
1229 const bool bigGapRefinement,
1230 const bool spreadGapSize,
1231 const label maxGlobalCells,
1232 const label maxLocalCells
1251 const bitSet& isOutsideFace,
1259 const bitSet& isOutsideFace,
1267 const scalar planarAngle,
1270 const label growIter
1324 const scalar maxLoadUnbalance,
1325 const label maxCellUnbalance
1335 const scalar maxLoadUnbalance,
1336 const label maxCellUnbalance
1346 const scalar maxLoadUnbalance,
1347 const label maxCellUnbalance
1353 const label maxGlobalCells,
1354 const label maxLocalCells,
1373 const bool handleSnapProblems,
1377 const bool useTopologicalSnapDetection,
1378 const bool removeEdgeConnectedCells,
1380 const label nErodeCellZones,
1388 const bool exitIfLeakPath,
1396 const bool samePatch,
1398 const bool useTopologicalSnapDetection,
1399 const bool removeEdgeConnectedCells,
1401 const scalar planarAngle,
1414 const label nBufferLayers,
1415 const label nErodeCellZones,
1422 const bool exitIfLeakPath,
1430 const label nBufferLayers,
1431 const label nErodeCellZones,
1518 const bool doInternalZones,
1519 const bool doBaffleZones
1529 const bool allowFreeStandingZoneFaces,
1530 const label nErodeCellZones,
1534 const bool exitIfLeakPath,
1547 const label insertPatchi,
1567 const word& masterPatch,
1568 const word& slavePatch,
1577 label& masterPatchID,
1578 label& slavePatchID,
1595 void nearestIntersection
1622 const vector& perturbVec,
1631 const vector& perturbVec,
1640 const vector& perturbVec,
1643 const label nRegions,
1647 const bool exitIfLeakPath,
1660 const bool exitIfLeakPath,
1747 const scalar minCos,
1748 const scalar concaveCos,
1749 const label mergeSize,
1758 const scalar minCos,
1759 const scalar concaveCos,
1800 const scalar minCos,
1863 const bool printCellLevel
1897 const bitSet& isMasterElem,
1911 template<
class EnumContainer>
1914 const EnumContainer& namedEnum,
1919 template<
class Type>
1923 const word& keyword,
1926 const Type& deflt =
Zero
1934 const word& keyword,
1943 const word& keyword,
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Dynamically sized Field. Similar to DynamicList, but inheriting from a Field instead of a List.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
A 1D vector of objects of type <T> with a fixed length <N>.
A HashTable similar to std::unordered_map.
An input stream of tokens.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
A HashTable to objects of type <T> with a label key.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Base class for writing coordSet(s) and tracks with fields.
Abstract base class for domain decomposition.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
A subset of mesh faces organised as a primitive patch.
A face is a list of labels corresponding to mesh vertices.
A class for handling file names.
Sends/receives parts of mesh+fvfields to neighbouring processors. Used in load balancing.
Mesh data needed to do the Finite Volume discretisation.
Refinement of (split) hexes using polyTopoChange.
option
Enumeration for the data type and search/match modes (bitmask).
@ REGEX
Regular expression.
Takes mesh with 'baffles' (= boundary faces sharing points). Determines for selected points on bounda...
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing processor-to-processor mapping information.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Helper class which maintains intersections of (changing) mesh with (static) surfaces.
static void calculateEdgeWeights(const polyMesh &mesh, const bitSet &isMasterEdge, const labelList &meshPoints, const edgeList &edges, scalarField &edgeWeights, scalarField &invSumWeight)
Helper: calculate edge weights (1/length).
label addFaceZone(const word &fzName, const word &masterPatch, const word &slavePatch, const surfaceZonesInfo::faceZoneType &fzType)
Add/lookup faceZone and update information. Return index of.
labelList getZones(const List< surfaceZonesInfo::faceZoneType > &fzTypes) const
Get zones of given type.
bool getFaceZoneInfo(const word &fzName, label &masterPatchID, label &slavePatchID, surfaceZonesInfo::faceZoneType &fzType) const
Lookup faceZone information. Return false if no information.
static T gAverage(const bitSet &isMasterElem, const UList< T > &values)
Helper: calculate average.
static void collectAndPrint(const UList< point > &points, const UList< T > &data)
Print list according to (collected and) sorted coordinate.
labelList intersectedFaces() const
Get faces with intersection.
labelList intersectedPoints() const
Get points on surfaces with intersection and boundary faces.
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.
autoPtr< mapPolyMesh > removeLimitShells(const label nBufferLayers, const label nErodeCellZones, const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const pointField &locationsInMesh, const wordList ®ionsInMesh, const pointField &locationsOutsideMesh)
Remove cells from limitRegions if level -1.
void printMeshInfo(const bool debug, const string &msg, const bool printCellLevel) const
Print some mesh stats.
MeshType
Enumeration for how to operate.
@ CASTELLATEDBUFFERLAYER2
static void removeFiles(const polyMesh &)
Helper: remove all relevant files from mesh instance.
label mergeEdgesUndo(const scalar minCos, const dictionary &motionDict)
Merge edges, maintain mesh quality. Return global number.
static bitSet getMasterPoints(const polyMesh &mesh, const labelList &meshPoints)
Determine master point for subset of points. If coupled.
static void mapBaffles(const polyMesh &mesh, const labelList &faceMap, List< labelPair > &baffles)
Map baffles after layer addition. Gets new-to-old face map.
static const Enum< MeshType > MeshTypeNames
static void testSyncPointList(const string &msg, const polyMesh &mesh, const List< scalar > &fld)
void checkData()
Debugging: check that all faces still obey start()>end().
autoPtr< mapPolyMesh > blockLeakFaces(const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const pointField &locationsInMesh, const wordList &zonesInMesh, const pointField &locationsOutsideMesh, const labelList &selectedSurfaces)
Baffle faces to break any leak from inside to outside.
static List< labelPair > subsetBaffles(const polyMesh &mesh, const labelList &zoneIDs, const List< labelPair > &baffles)
Subset baffles according to zones.
const List< Tuple2< mapType, labelList > > & userFaceData() const
Additional face data that is maintained across.
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.
void updateMesh(const mapPolyMesh &, const labelList &changedFaces)
Update for external change to mesh. changedFaces are in new mesh.
writeType
Enumeration for what to write. Used as a bit-pattern.
void doSplitFaces(const labelList &splitFaces, const labelPairList &splits, const labelPairList &splitPatches, polyTopoChange &meshMod) const
Split faces into two.
void mergeFreeStandingBaffles(const bool samePatch, const snapParameters &snapParams, const bool useTopologicalSnapDetection, const bool removeEdgeConnectedCells, const scalarField &perpendicularAngle, const scalar planarAngle, const dictionary &motionDict, Time &runTime, const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const pointField &locationsInMesh, const pointField &locationsOutsideMesh)
Merge free-standing baffles.
const hexRef8 & meshCutter() const
Reference to meshcutting engine.
autoPtr< mapPolyMesh > mergePoints(const labelList &pointToDuplicate)
Merge duplicate points.
const labelList & surfaceIndex() const
Per start-end edge the index of the surface hit.
autoPtr< mapPolyMesh > doRestorePoints(removePoints &pointRemover, const labelList &facesToRestore)
autoPtr< mapPolyMesh > dupNonManifoldPoints()
Find boundary points that connect to more than one cell.
labelList refineCandidates(const pointField &keepPoints, const scalar curvature, const scalar planarAngle, const bool featureRefinement, const bool featureDistanceRefinement, const bool internalRefinement, const bool surfaceRefinement, const bool curvatureRefinement, const bool smallFeatureRefinement, const bool gapRefinement, const bool bigGapRefinement, const bool spreadGapSize, const label maxGlobalCells, const label maxLocalCells) const
Calculate list of cells to refine.
const shellSurfaces & limitShells() const
Reference to limit shells (regions).
label countFaceDirs(const bitSet &isOutsideFace, const label celli) const
Count number of faces on cell that are in set.
labelList growFaceCellFace(const labelUList &set) const
void dumpRefinementLevel() const
Write refinement level as volScalarFields for postprocessing.
autoPtr< mapPolyMesh > removeGapCells(const scalar planarAngle, const labelList &minSurfaceLevel, const labelList &globalToMasterPatch, const label growIter)
Detect gapRefinement cells and remove them.
autoPtr< mapPolyMesh > dupNonManifoldBoundaryPoints()
Find boundary points that are on faceZones of type boundary.
static const Enum< writeType > writeTypeNames
static int readFlags(const EnumContainer &namedEnum, const wordList &words)
Helper: convert wordList into bit pattern using provided Enum.
void distribute(const mapDistributePolyMesh &)
Update local numbering for mesh redistribution.
static void splitFace(const face &f, const labelPair &split, face &f0, face &f1)
Helper: split face into:
void getZoneFaces(const labelList &zoneIDs, labelList &faceZoneID, labelList &ownPatch, labelList &neiPatch, labelList &nBaffles) const
Get per-face information (faceZone, master/slave patch).
const Map< label > & faceToCoupledPatch() const
For faces originating from processor faces store the original.
autoPtr< mapPolyMesh > mergeBaffles(const List< labelPair > &, const Map< label > &faceToPatch)
Merge baffles. Gets pairs of faces and boundary faces to move.
static label addPatch(fvMesh &, const word &name, const dictionary &)
Helper:add patch to mesh. Update all registered fields.
label countHits() const
Count number of intersections (local).
void markOutsideFaces(const labelList &cellLevel, const labelList &neiLevel, const labelList &refineCell, bitSet &isOutsideFace) const
Mark faces on interface between set and rest.
ClassName("meshRefinement")
Runtime type information.
autoPtr< mapPolyMesh > doRemoveCells(const labelList &cellsToRemove, const labelList &exposedFaces, const labelList &exposedPatchIDs, removeCells &cellRemover)
Remove cells. Put exposedFaces into exposedPatchIDs.
mapType
Enumeration for how userdata is to be mapped upon refinement.
@ KEEPALL
have slaves (upon refinement) from master
@ REMOVE
set value to -1 any face that was refined
@ MASTERONLY
maintain master only
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.
static label findCell(const polyMesh &, const vector &perturbVec, const point &p)
Find cell point is in. Uses optional perturbation to re-test.
const refinementSurfaces & surfaces() const
Reference to surface search engines.
label mergePatchFaces(const scalar minCos, const scalar concaveCos, const label mergeSize, const labelList &patchIDs, const meshRefinement::FaceMergeType mergeType)
Merge coplanar faces if sets are of size mergeSize.
word timeName() const
Replacement for Time::timeName() that returns oldInstance (if overwrite_).
debugType
Enumeration for what to debug. Used as a bit-pattern.
static bitSet getMasterEdges(const polyMesh &mesh, const labelList &meshEdges)
Determine master edge for subset of edges. If coupled.
const shellSurfaces & shells() const
Reference to refinement shells (regions).
labelList meshedPatches() const
Get patchIDs for patches added in addMeshedPatch.
static void checkCoupledFaceZones(const polyMesh &)
Helper function: check that face zones are synced.
labelList countEdgeFaces(const uindirectPrimitivePatch &pp) const
Count number of faces per patch edge. Parallel consistent.
bool isNormalGap(const scalar planarCos, const label level0, const vector &point0, const vector &normal0, const label level1, const vector &point1, const vector &normal1) const
Is local topology a small gap normal to the test vector.
labelList collectFaces(const labelList &candidateFaces, const labelHashSet &set) const
static void weightedSum(const polyMesh &mesh, const bitSet &isMasterEdge, const labelList &meshPoints, const edgeList &edges, const scalarField &edgeWeights, const Field< Type > &data, Field< Type > &sum)
Helper: weighted sum (over all subset of mesh points) by.
autoPtr< mapPolyMesh > splitMeshRegions(const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const pointField &locationsInMesh, const pointField &locationsOutsideMesh, const bool exitIfLeakPath, const refPtr< coordSetWriter > &leakPathFormatter)
Split mesh. Keep part containing point. Return empty map if.
autoPtr< mapDistributePolyMesh > refineAndBalance(const string &msg, decompositionMethod &decomposer, fvMeshDistribute &distributor, const labelList &cellsToRefine, const scalar maxLoadUnbalance, const label maxCellUnbalance)
Refine some cells and rebalance.
label addMeshedPatch(const word &name, const dictionary &)
Add patch originating from meshing. Update meshedPatches_.
List< Tuple2< mapType, labelList > > & userFaceData()
const fvMesh & mesh() const
Reference to mesh.
bool isGap(const scalar, const vector &, const vector &, const vector &, const vector &) const
Is local topology a small gap?
const refinementFeatures & features() const
Reference to feature edge mesh.
void testSyncBoundaryFaceList(const scalar mergeDistance, const string &, const UList< T > &, const UList< T > &) const
Compare two lists over all boundary faces.
void selectSeparatedCoupledFaces(boolList &) const
Select coupled faces that are not collocated.
static tmp< pointVectorField > makeDisplacementField(const pointMesh &pMesh, const labelList &adaptPatchIDs)
Helper function to make a pointVectorField with correct.
const word & oldInstance() const
(points)instance of mesh upon construction
static label findRegion(const polyMesh &, const labelList &cellRegion, const vector &perturbVec, const point &p)
Find region point is in. Uses optional perturbation to re-test.
label addPointZone(const word &name)
Add pointZone if does not exist. Return index of zone.
void selectIntersectedFaces(const labelList &surfaces, boolList &isBlockedFace) const
Faces currently on boundary or intersected by surface.
static label findRegions(const polyMesh &, const vector &perturbVec, const pointField &locationsInMesh, const pointField &locationsOutsideMesh, const label nRegions, labelList &cellRegion, const boolList &blockedFace, const bool exitIfLeakPath, const refPtr< coordSetWriter > &leakPathFormatter)
Find regions points are in.
autoPtr< mapPolyMesh > doRemovePoints(removePoints &pointRemover, const boolList &pointCanBeDeleted)
scalar mergeDistance() const
FaceMergeType
Enumeration for what to do with co-planar patch faces on a single.
autoPtr< mapPolyMesh > directionalRefine(const string &msg, const direction cmpt, const labelList &cellsToRefine)
Directionally refine in direction cmpt.
autoPtr< mapPolyMesh > refine(const labelList &cellsToRefine)
Refine some cells.
static autoPtr< indirectPrimitivePatch > makePatch(const polyMesh &, const labelList &)
Create patch from set of patches.
void storeData(const labelList &pointsToStore, const labelList &facesToStore, const labelList &cellsToStore)
Signal points/face/cells for which to store data.
static void updateList(const labelList &newToOld, const T &nullValue, List< T > &elems)
Helper: reorder list according to map.
bool overwrite() const
Overwrite the mesh?
bool write() const
Write mesh and all data.
autoPtr< mapPolyMesh > mergeZoneBaffles(const bool doInternalZones, const bool doBaffleZones)
Merge all baffles on faceZones.
labelList directionalRefineCandidates(const label maxGlobalCells, const label maxLocalCells, const labelList ¤tLevel, const direction dir) const
Calculate list of cells to directionally refine.
void growSet(const labelList &neiLevel, const bitSet &isOutsideFace, labelList &refineCell, label &nRefine) const
Add one layer of cells to set.
autoPtr< mapDistributePolyMesh > balance(const bool keepZoneFaces, const bool keepBaffles, const labelList &singleProcPoints, const scalarField &cellWeights, decompositionMethod &decomposer, fvMeshDistribute &distributor)
Redecompose according to cell count.
label mergePatchFacesUndo(const scalar minCos, const scalar concaveCos, const labelList &patchIDs, const dictionary &motionDict, const labelList &preserveFaces, const meshRefinement::FaceMergeType mergeType)
Merge coplanar faces. preserveFaces is != -1 for faces.
label splitFacesUndo(const labelList &splitFaces, const labelPairList &splits, const labelPairList &splitPatches, const dictionary &motionDict, labelList &duplicateFace, List< labelPair > &baffles)
Split faces along diagonal. Maintain mesh quality. Return.
autoPtr< mapPolyMesh > removeLeakCells(const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const pointField &locationsInMesh, const pointField &locationsOutsideMesh, const labelList &selectedSurfaces)
Remove minimum amount of cells to break any leak from.
void baffleAndSplitMesh(const bool handleSnapProblems, const snapParameters &snapParams, const bool useTopologicalSnapDetection, const bool removeEdgeConnectedCells, const scalarField &perpendicularAngle, const label nErodeCellZones, const dictionary &motionDict, Time &runTime, const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const pointField &locationsInMesh, const wordList ®ionsInMesh, const pointField &locationsOutsideMesh, const bool exitIfLeakPath, const refPtr< coordSetWriter > &leakPathFormatter, refPtr< surfaceWriter > &surfFormatter)
Split off unreachable areas of mesh.
void setInstance(const fileName &)
Set instance of all local IOobjects.
autoPtr< mapPolyMesh > createBaffles(const labelList &ownPatch, const labelList &neiPatch)
Create baffle for every internal face where ownPatch != -1.
static const Enum< debugType > debugTypeNames
autoPtr< mapDistributePolyMesh > balanceAndRefine(const string &msg, decompositionMethod &decomposer, fvMeshDistribute &distributor, const labelList &cellsToRefine, const scalar maxLoadUnbalance, const label maxCellUnbalance)
Balance before refining some cells.
static label appendPatch(fvMesh &, const label insertPatchi, const word &, const dictionary &)
Helper:append patch to end of mesh.
MeshType meshType() const
Mode of meshing.
static writeType writeLevel()
Get/set write level.
autoPtr< mapPolyMesh > createZoneBaffles(const labelList &zoneIDs, List< labelPair > &baffles, labelList &originatingFaceZone)
Create baffles for faces on faceZones. Return created baffles.
void dumpIntersections(const fileName &prefix) const
Debug: Write intersection information to OBJ format.
labelList detectLeakCells(const boolList &isBlockedFace, const labelList &leakFaces, const labelList &seedCells) const
Return list of cells to block by walking from the seedCells.
Mesh representing a set of points created from polyMesh.
Mesh consisting of general polyhedral cells.
Direct mesh changes based on v1.3 polyTopoChange syntax.
A class for managing references or pointers (no reference counting).
Container with cells to refine. Refinement given as single direction.
Encapsulates queries for features.
Container for data on surfaces used for surface-driven refinement. Contains all the data about the le...
Given list of cells to remove, insert all the topology changes.
Removes selected points from mesh and updates faces using these points.
Encapsulates queries for volume refinement ('refine all cells within shell').
Simple container to keep together snap specific information.
Contains information about location on a triSurface.
faceZoneType
What to do with faceZone faces.
A class for managing temporary objects.
An enumeration wrapper for classification of a location as being inside/outside of a volume.
A class for handling words, derived from Foam::string.
#define ClassName(TypeNameString)
Add typeName information from argument TypeNameString to a class.
const volScalarField & p0
static bool split(const std::string &line, std::string &key, std::string &val)
const labelIOList & zoneIDs
A HashTable to objects of type <T> with a labelPair key. The hashing is based on labelPair (FixedList...
const bitSet isBlockedFace(intersectedFaces())
Namespace for handling debugging switches.
List< edge > edgeList
List of edge.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Pair< label > labelPair
A pair of labels.
List< word > wordList
List of word.
HashTable< word, wordPair, Foam::Hash< wordPair > > wordPairHashTable
HashTable of wordPair.
List< labelPair > labelPairList
List of labelPair.
List< label > labelList
A List of labels.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
IOList< label > labelIOList
IO for a List of label.
HashTable< T, labelPair, Foam::Hash< labelPair > > LabelPairMap
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
List< bool > boolList
A List of bools.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
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.
UList< label > labelUList
A UList of labels.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
DynamicID< faceZoneMesh > faceZoneID
Foam::faceZoneID.
Forwards and collection of common point field types.
#define FOAM_NO_DANGLING_REFERENCE