Go to the documentation of this file.
68 #ifndef dynamicRefineFvMesh_H
69 #define dynamicRefineFvMesh_H
133 const label maxCells,
134 const label maxRefinement,
135 const scalar refineLevel,
150 const scalar minLevel,
151 const scalar maxLevel
157 const scalar lowerRefineLevel,
158 const scalar upperRefineLevel,
166 const label maxCells,
167 const label maxRefinement,
174 const scalar unrefineLevel,
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
scalarField error(const scalarField &fld, const scalar minLevel, const scalar maxLevel) const
const PackedBoolList & protectedCell() const
Cells which should not be refined/unrefined.
virtual labelList selectRefineCells(const label maxCells, const label maxRefinement, const PackedBoolList &candidateCell) const
Subset candidate cells for refinement.
static label count(const PackedBoolList &, const unsigned int)
Count set/unset elements in packedlist.
virtual ~dynamicRefineFvMesh()
Destructor.
virtual labelList selectUnrefinePoints(const scalar unrefineLevel, const PackedBoolList &markedCell, const scalarField &pFld) const
Select points that can be unrefined.
TypeName("dynamicRefineFvMesh")
Runtime type information.
compressionType
Enumeration for the format of data in the stream.
Abstract base class for geometry and/or topology changing fvMesh.
hexRef8 meshCutter_
Mesh cutting engine.
virtual bool update()
Update the mesh for both mesh motion and topology change.
virtual void selectRefineCandidates(const scalar lowerRefineLevel, const scalar upperRefineLevel, const scalarField &vFld, PackedBoolList &candidateCell) const
Select 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.
PackedBoolList protectedCell_
Protected cells (usually since not hexes)
Switch dumpLevel_
Dump cellLevel for postprocessing.
void readDict()
Read the projection parameters from dictionary.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp) const
Write using given format, version and compression.
void calculateProtectedCells(PackedBoolList &unrefineableCell) const
Calculate cells that cannot be refined since would trigger.
scalarField maxCellField(const volScalarField &) const
Get point max of connected cell.
dynamicRefineFvMesh(const dynamicRefineFvMesh &)
Disallow default bitwise copy construct.
autoPtr< mapPolyMesh > refine(const labelList &)
Refine cells. Update mesh and fields.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
void checkEightAnchorPoints(PackedBoolList &protectedCell, label &nProtected) const
Check all cells have 8 anchor points.
void extendMarkedCells(PackedBoolList &markedCell) const
Extend markedCell with cell-face-cell.
An STL-conforming hash table.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
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.
void operator=(const dynamicRefineFvMesh &)
Disallow default bitwise assignment.
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.
Generic GeometricField class.
scalarField maxPointField(const scalarField &) const
Get per cell max of connected point.
label nRefinementIterations_
Number of refinement/unrefinement steps done so far.
PackedBoolList & protectedCell()
Cells which should not be refined/unrefined.
streamFormat
Enumeration for the format of data in the stream.