FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

A fvMesh with built-in refinement. More...

#include <dynamicFvMesh/dynamicRefineFvMesh.H>


Detailed Description

A fvMesh with built-in refinement.

Determines which cells to refine/unrefine and does all in update().

Source files

Definition at line 54 of file dynamicRefineFvMesh.H.

Inheritance diagram for dynamicRefineFvMesh:
Collaboration diagram for dynamicRefineFvMesh:

List of all members.

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)

Constructor & Destructor Documentation

~dynamicRefineFvMesh (  ) [virtual]

Definition at line 1009 of file dynamicRefineFvMesh.C.


Member Function Documentation

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]
autoPtr< mapPolyMesh > unrefine ( const labelList &   splitPoints  ) [protected]
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]
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]
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]
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 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().


Member Data Documentation

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().

Protected cells (usually since not hexes)

Definition at line 73 of file dynamicRefineFvMesh.H.

Referenced by dynamicRefineFvMesh::dynamicRefineFvMesh(), and dynamicRefineFvMesh::protectedCell().


The documentation for this class was generated from the following files: