A fvMesh with built-in refinement. More...
#include <dynamicFvMesh/dynamicRefineFvMesh.H>
A fvMesh with built-in refinement.
Determines which cells to refine/unrefine and does all in update().
Definition at line 54 of file dynamicRefineFvMesh.H.
Public Member Functions | |
TypeName ("dynamicRefineFvMesh") | |
Runtime type information.
| |
dynamicRefineFvMesh (const IOobject &io) | |
Construct from IOobject.
| |
virtual | ~dynamicRefineFvMesh () |
const hexRef8 & | meshCutter () const |
Direct access to the refinement engine.
| |
const PackedBoolList & | protectedCell () const |
Cells which should not be refined/unrefined.
| |
PackedBoolList & | protectedCell () |
Cells which should not be refined/unrefined.
| |
virtual bool | update () |
Update the mesh for both mesh motion and topology change.
| |
virtual bool | writeObject (IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp) const |
Write using given format, version and compression.
| |
Protected Member Functions | |
void | calculateProtectedCells (PackedBoolList &unrefineableCell) const |
Calculate cells that cannot be refined since would trigger.
| |
void | readDict () |
Read the projection parameters from dictionary.
| |
autoPtr< mapPolyMesh > | refine (const labelList &) |
Refine cells. Update mesh and fields.
| |
autoPtr< mapPolyMesh > | unrefine (const labelList &) |
Unrefine cells. Gets passed in centre points of cells to combine.
| |
scalar | getRefineLevel (const label maxCells, const label maxRefinement, const scalar refineLevel, const scalarField &) const |
Calculates approximate value for refinement level so.
| |
scalarField | maxPointField (const scalarField &) const |
Get per cell max of connected point.
| |
scalarField | minCellField (const volScalarField &) const |
Get point min of connected cell.
| |
scalarField | cellToPoint (const scalarField &vFld) const |
scalarField | error (const scalarField &fld, const scalar minLevel, const scalar maxLevel) const |
virtual void | selectRefineCandidates (const scalar lowerRefineLevel, const scalar upperRefineLevel, const scalarField &vFld, PackedBoolList &candidateCell) const |
Select candidate cells for refinement.
| |
virtual labelList | selectRefineCells (const label maxCells, const label maxRefinement, const PackedBoolList &candidateCell) const |
Subset candidate cells for refinement.
| |
virtual labelList | selectUnrefinePoints (const scalar unrefineLevel, const PackedBoolList &markedCell, const scalarField &pFld) const |
Select points that can be unrefined.
| |
void | extendMarkedCells (PackedBoolList &markedCell) const |
Extend markedCell with cell-face-cell.
| |
Static Protected Member Functions | |
static label | count (const PackedBoolList &, const unsigned int) |
Count set/unset elements in packedlist.
| |
Protected Attributes | |
hexRef8 | meshCutter_ |
Mesh cutting engine.
| |
Switch | dumpLevel_ |
Dump cellLevel for postprocessing.
| |
List< Pair< word > > | correctFluxes_ |
Fluxes to map.
| |
label | nRefinementIterations_ |
Number of refinement/unrefinement steps done so far.
| |
PackedBoolList | protectedCell_ |
Protected cells (usually since not hexes)
|
dynamicRefineFvMesh | ( | const IOobject & | io ) | [explicit]
|
Construct from IOobject.
Definition at line 865 of file dynamicRefineFvMesh.C.
References hexRef8::cellLevel(), primitiveMesh::cellPoints(), PackedList< nBits >::clear(), f(), polyMesh::faceNeighbour(), polyMesh::faceOwner(), polyMesh::faces(), forAll, Foam::max(), dynamicRefineFvMesh::meshCutter_, primitiveMesh::nFaces(), primitiveMesh::nInternalFaces(), hexRef8::pointLevel(), dynamicRefineFvMesh::protectedCell_, dynamicRefineFvMesh::readDict(), Foam::reduce(), PackedList< nBits >::set(), syncTools::swapFaceList(), and syncTools::syncFaceList().
~dynamicRefineFvMesh | ( | ) | [virtual]
|
Definition at line 1009 of file dynamicRefineFvMesh.C.
label count | ( | const PackedBoolList & | l, |
const unsigned int | val | ||
) | [static, protected]
|
Count set/unset elements in packedlist.
Definition at line 50 of file dynamicRefineFvMesh.C.
References forAll, and PackedList< nBits >::get().
void calculateProtectedCells | ( | PackedBoolList & | unrefineableCell ) | const [protected]
|
Calculate cells that cannot be refined since would trigger.
refinement of protectedCell_ (since 2:1 refinement cascade)
Definition at line 68 of file dynamicRefineFvMesh.C.
References PackedList< nBits >::clear(), forAll, PackedList< nBits >::get(), Foam::returnReduce(), PackedList< nBits >::set(), syncTools::swapBoundaryFaceList(), and syncTools::syncFaceList().
void readDict | ( | ) | [protected]
|
Read the projection parameters from dictionary.
Definition at line 173 of file dynamicRefineFvMesh.C.
References dynamicRefineFvMesh::correctFluxes_, dynamicRefineFvMesh::dumpLevel_, dictionary::lookup(), IOobject::MUST_READ, IOobject::NO_WRITE, dictionary::subDict(), and fvMesh::time().
Referenced by dynamicRefineFvMesh::dynamicRefineFvMesh().
autoPtr< mapPolyMesh > refine | ( | const labelList & | cellsToRefine ) | [protected]
|
Refine cells. Update mesh and fields.
Definition at line 199 of file dynamicRefineFvMesh.C.
References Foam::abort(), GeometricField< Type, PatchField, GeoMesh >::boundaryField(), polyTopoChange::changeMesh(), Foam::endl(), Foam::FatalError, FatalErrorIn, forAll, forAllConstIter, Foam::Info, Foam::fvc::interpolate(), Foam::nl, fvPatch::patch(), fvsPatchField< Type >::patch(), phi, phiU(), Foam::returnReduce(), PackedList< nBits >::set(), List< T >::size(), and polyPatch::start().
Referenced by dynamicRefineFvMesh::update().
autoPtr< mapPolyMesh > unrefine | ( | const labelList & | splitPoints ) | [protected]
|
Unrefine cells. Gets passed in centre points of cells to combine.
Definition at line 411 of file dynamicRefineFvMesh.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), polyTopoChange::changeMesh(), Foam::endl(), forAll, forAllConstIter, Foam::Info, Foam::fvc::interpolate(), pFaces, phi, phiU(), Foam::returnReduce(), PackedList< nBits >::set(), and List< T >::size().
Referenced by dynamicRefineFvMesh::update().
scalar getRefineLevel | ( | const label | maxCells, |
const label | maxRefinement, | ||
const scalar | refineLevel, | ||
const scalarField & | |||
) | const [protected]
|
Calculates approximate value for refinement level so.
we don't go above maxCell
scalarField maxPointField | ( | const scalarField & | pFld ) | const [protected]
|
Get per cell max of connected point.
Definition at line 565 of file dynamicRefineFvMesh.C.
References forAll, Foam::max(), primitiveMesh::nCells(), and primitiveMesh::pointCells().
scalarField minCellField | ( | const volScalarField & | vFld ) | const [protected]
|
Get point min of connected cell.
Definition at line 583 of file dynamicRefineFvMesh.C.
References forAll, Foam::min(), primitiveMesh::nPoints(), and primitiveMesh::pointCells().
Referenced by dynamicRefineFvMesh::update().
scalarField cellToPoint | ( | const scalarField & | vFld ) | const [protected]
|
Definition at line 601 of file dynamicRefineFvMesh.C.
References forAll, primitiveMesh::nPoints(), primitiveMesh::pointCells(), List< T >::size(), and Foam::sum().
scalarField error | ( | const scalarField & | fld, |
const scalar | minLevel, | ||
const scalar | maxLevel | ||
) | const [protected]
|
Definition at line 622 of file dynamicRefineFvMesh.C.
References forAll, Foam::min(), and List< T >::size().
void selectRefineCandidates | ( | const scalar | lowerRefineLevel, |
const scalar | upperRefineLevel, | ||
const scalarField & | vFld, | ||
PackedBoolList & | candidateCell | ||
) | const [protected, virtual]
|
Select candidate cells for refinement.
Definition at line 644 of file dynamicRefineFvMesh.C.
References forAll, and PackedList< nBits >::set().
Referenced by dynamicRefineFvMesh::update().
labelList selectRefineCells | ( | const label | maxCells, |
const label | maxRefinement, | ||
const PackedBoolList & | candidateCell | ||
) | const [protected, virtual]
|
Subset candidate cells for refinement.
Definition at line 678 of file dynamicRefineFvMesh.C.
References DynamicList< T, SizeInc, SizeMult, SizeDiv >::append(), PackedList< nBits >::empty(), Foam::endl(), forAll, PackedList< nBits >::get(), Foam::Info, Foam::returnReduce(), DynamicList< T, SizeInc, SizeMult, SizeDiv >::shrink(), and List< T >::size().
Referenced by dynamicRefineFvMesh::update().
labelList selectUnrefinePoints | ( | const scalar | unrefineLevel, |
const PackedBoolList & | markedCell, | ||
const scalarField & | pFld | ||
) | const [protected, virtual]
|
Select points that can be unrefined.
Definition at line 765 of file dynamicRefineFvMesh.C.
References Foam::endl(), forAll, PackedList< nBits >::get(), Foam::Info, Foam::returnReduce(), and List< T >::size().
Referenced by dynamicRefineFvMesh::update().
void extendMarkedCells | ( | PackedBoolList & | markedCell ) | const [protected]
|
Extend markedCell with cell-face-cell.
Definition at line 824 of file dynamicRefineFvMesh.C.
References primitiveMesh::cells(), polyMesh::faceNeighbour(), polyMesh::faceOwner(), forAll, PackedList< nBits >::get(), primitiveMesh::nFaces(), primitiveMesh::nInternalFaces(), PackedList< nBits >::set(), and syncTools::syncFaceList().
Referenced by dynamicRefineFvMesh::update().
TypeName | ( | "dynamicRefineFvMesh" | ) |
Runtime type information.
const hexRef8& meshCutter | ( | ) | const [inline]
|
Direct access to the refinement engine.
Definition at line 180 of file dynamicRefineFvMesh.H.
References dynamicRefineFvMesh::meshCutter_.
Referenced by dynamicRefineFvMesh::update().
const PackedBoolList& protectedCell | ( | ) | const [inline]
|
Cells which should not be refined/unrefined.
Definition at line 186 of file dynamicRefineFvMesh.H.
References dynamicRefineFvMesh::protectedCell_.
PackedBoolList& protectedCell | ( | ) | [inline]
|
Cells which should not be refined/unrefined.
Definition at line 192 of file dynamicRefineFvMesh.H.
References dynamicRefineFvMesh::protectedCell_.
bool update | ( | ) | [virtual]
|
Update the mesh for both mesh motion and topology change.
Implements dynamicFvMesh.
Definition at line 1015 of file dynamicRefineFvMesh.C.
References polyMesh::changing(), Foam::exit(), dynamicRefineFvMesh::extendMarkedCells(), Foam::FatalError, FatalErrorIn, forAll, polyMesh::globalData(), hexRef8::history(), dictionary::lookup(), dynamicRefineFvMesh::meshCutter(), dynamicRefineFvMesh::minCellField(), IOobject::MUST_READ, primitiveMesh::nCells(), Foam::nl, IOobject::NO_WRITE, dynamicRefineFvMesh::nRefinementIterations_, Foam::readScalar(), dynamicRefineFvMesh::refine(), Foam::returnReduce(), dynamicRefineFvMesh::selectRefineCandidates(), dynamicRefineFvMesh::selectRefineCells(), dynamicRefineFvMesh::selectUnrefinePoints(), List< T >::size(), dictionary::subDict(), fvMesh::time(), timeIndex, List< T >::transfer(), and dynamicRefineFvMesh::unrefine().
bool writeObject | ( | IOstream::streamFormat | fmt, |
IOstream::versionNumber | ver, | ||
IOstream::compressionType | cmp | ||
) | const [virtual]
|
Write using given format, version and compression.
Reimplemented from objectRegistry.
Definition at line 1218 of file dynamicRefineFvMesh.C.
References IOobject::AUTO_WRITE, Foam::dimless, forAll, IOobject::NO_READ, timeName, regIOobject::write(), and fvMesh::writeObjects().
hexRef8 meshCutter_ [protected]
|
Mesh cutting engine.
Definition at line 61 of file dynamicRefineFvMesh.H.
Referenced by dynamicRefineFvMesh::dynamicRefineFvMesh(), and dynamicRefineFvMesh::meshCutter().
Switch dumpLevel_ [protected]
|
Dump cellLevel for postprocessing.
Definition at line 64 of file dynamicRefineFvMesh.H.
Referenced by dynamicRefineFvMesh::readDict().
List<Pair<word> > correctFluxes_ [protected]
|
Fluxes to map.
Definition at line 67 of file dynamicRefineFvMesh.H.
Referenced by dynamicRefineFvMesh::readDict().
label nRefinementIterations_ [protected]
|
Number of refinement/unrefinement steps done so far.
Definition at line 70 of file dynamicRefineFvMesh.H.
Referenced by dynamicRefineFvMesh::update().
PackedBoolList protectedCell_ [protected]
|
Protected cells (usually since not hexes)
Definition at line 73 of file dynamicRefineFvMesh.H.
Referenced by dynamicRefineFvMesh::dynamicRefineFvMesh(), and dynamicRefineFvMesh::protectedCell().