Go to the documentation of this file.
75 #ifndef dynamicRefineFvMesh_H
76 #define dynamicRefineFvMesh_H
136 const label maxCells,
137 const label maxRefinement,
138 const scalar refineLevel,
153 const scalar minLevel,
154 const scalar maxLevel
160 const scalar lowerRefineLevel,
161 const scalar upperRefineLevel,
169 const label maxCells,
170 const label maxRefinement,
171 const bitSet& candidateCell
177 const scalar unrefineLevel,
234 const bool doInit=
true
245 virtual bool init(
const bool doInit);
Defines the attributes of an object for which implicit objectRegistry management is supported,...
scalarField error(const scalarField &fld, const scalar minLevel, const scalar maxLevel) const
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
void extendMarkedCells(bitSet &markedCell) const
Extend markedCell with cell-face-cell.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
bitSet protectedCell_
Protected cells (usually since not hexes)
TypeName("dynamicRefineFvMesh")
Runtime type information.
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Abstract base class for geometry and/or topology changing fvMesh.
hexRef8 meshCutter_
Mesh cutting engine.
virtual void mapFields(const mapPolyMesh &mpm)
Map all fields in time using given map.
virtual bool update()
Update the mesh for both mesh motion and topology change.
void mapNewInternalFaces(const labelList &faceMap, GeometricField< T, fvsPatchField, surfaceMesh > &)
Map single non-flux surface<Type>Field.
virtual labelList selectRefineCells(const label maxCells, const label maxRefinement, const bitSet &candidateCell) const
Subset candidate cells for refinement.
HashTable< word > correctFluxes_
Fluxes to map.
scalar getRefineLevel(const label maxCells, const label maxRefinement, const scalar refineLevel, const scalarField &) const
Calculates approximate value for refinement level so.
void readDict()
Read the projection parameters from dictionary.
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write using stream options.
void checkEightAnchorPoints(bitSet &protectedCell) const
Check all cells have 8 anchor points.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
virtual void selectRefineCandidates(const scalar lowerRefineLevel, const scalar upperRefineLevel, const scalarField &vFld, bitSet &candidateCell) const
Select candidate cells for refinement.
The IOstreamOption is a simple container for options an IOstream can normally have.
scalarField maxCellField(const volScalarField &) const
Get point max of connected cell.
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;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
virtual autoPtr< mapPolyMesh > refine(const labelList &)
Refine cells. Update mesh and fields.
const bitSet & protectedCell() const
Cells which should not be refined/unrefined.
bitSet & protectedCell()
Cells which should not be refined/unrefined.
virtual labelList selectUnrefinePoints(const scalar unrefineLevel, const bitSet &markedCell, const scalarField &pFld) const
Select points that can be unrefined.
virtual ~dynamicRefineFvMesh()=default
Destructor.
A HashTable similar to std::unordered_map.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
void calculateProtectedCells(bitSet &unrefineableCell) const
Calculate cells that cannot be refined since would trigger.
Refinement of (split) hexes using polyTopoChange.
A fvMesh with built-in refinement.
bool dumpLevel_
Dump cellLevel for post-processing.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
virtual autoPtr< mapPolyMesh > unrefine(const labelList &)
Unrefine cells. Gets passed in centre points of cells to combine.
scalarField cellToPoint(const scalarField &vFld) const
const hexRef8 & meshCutter() const
Direct access to the refinement engine.
scalarField maxPointField(const scalarField &) const
Get per cell max of connected point.
label nRefinementIterations_
Number of refinement/unrefinement steps done so far.
const surfaceVectorField & Sf() const
Return cell face area vectors.