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);
269 virtual void mapFields(
const mapPolyMesh& mpm);
277 IOstreamOption streamOpt,
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
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
TypeName("dynamicRefineFvMesh")
virtual bool init(const bool doInit)
Abstract base class for geometry and/or topology changing fvMesh.
virtual void mapFields(const mapPolyMesh &mpm)
void mapNewInternalFaces(const labelList &faceMap, GeometricField< T, fvsPatchField, surfaceMesh > &)
virtual labelList selectRefineCells(const label maxCells, const label maxRefinement, const bitSet &candidateCell) const
HashTable< word > correctFluxes_
scalar getRefineLevel(const label maxCells, const label maxRefinement, const scalar refineLevel, const scalarField &) const
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
void checkEightAnchorPoints(bitSet &protectedCell) const
Generic templated field type.
const surfaceScalarField & magSf() const
virtual void selectRefineCandidates(const scalar lowerRefineLevel, const scalar upperRefineLevel, const scalarField &vFld, bitSet &candidateCell) const
The IOstreamOption is a simple container for options an IOstream can normally have.
scalarField maxCellField(const volScalarField &) const
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 &)
const bitSet & protectedCell() const
virtual labelList selectUnrefinePoints(const scalar unrefineLevel, const bitSet &markedCell, const scalarField &pFld) const
virtual ~dynamicRefineFvMesh()=default
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
Refinement of (split) hexes using polyTopoChange.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
A fvMesh with built-in refinement.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
virtual autoPtr< mapPolyMesh > unrefine(const labelList &)
scalarField cellToPoint(const scalarField &vFld) const
const hexRef8 & meshCutter() const
Generic GeometricField class.
scalarField maxPointField(const scalarField &) const
label nRefinementIterations_
const surfaceVectorField & Sf() const