refiner Class Reference

Dynamic mesh refinement/unrefinement based on volScalarField values. More...

Inheritance diagram for refiner:
Collaboration diagram for refiner:

Public Member Functions

 TypeName ("refiner")
 Runtime type information. More...
 
 refiner (fvMesh &mesh, const dictionary &dict)
 Construct from fvMesh and dictionary. More...
 
 refiner (const refiner &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~refiner ()
 Destructor. More...
 
const hexRef8meshCutter () const
 Direct access to the refinement engine. More...
 
const PackedBoolListprotectedCell () const
 Cells which should not be refined/unrefined. More...
 
PackedBoolListprotectedCell ()
 Cells which should not be refined/unrefined. More...
 
virtual bool update ()
 Update the mesh for both mesh motion and topology change. More...
 
virtual void topoChange (const polyTopoChangeMap &)
 Update corresponding to the given map. More...
 
virtual void mapMesh (const polyMeshMap &)
 Update from another mesh using the given map. More...
 
virtual void distribute (const polyDistributionMap &)
 Update corresponding to the given distribution map. More...
 
virtual bool write (const bool write=true) const
 Write using given format, version and compression. More...
 
void operator= (const refiner &)=delete
 Disallow default bitwise assignment. More...
 
- Public Member Functions inherited from fvMeshTopoChanger
 TypeName ("fvMeshTopoChanger")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, fvMeshTopoChanger, fvMesh,(fvMesh &mesh, const dictionary &dict),(mesh, dict))
 
 fvMeshTopoChanger (fvMesh &)
 Construct from fvMesh. More...
 
 fvMeshTopoChanger (const fvMeshTopoChanger &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~fvMeshTopoChanger ()
 Destructor. More...
 
fvMeshmesh ()
 Return the fvMesh. More...
 
const fvMeshmesh () const
 Return the fvMesh. More...
 
virtual bool dynamic () const
 Is mesh dynamic, i.e. might it change? More...
 
void operator= (const fvMeshTopoChanger &)=delete
 Disallow default bitwise assignment. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from fvMeshTopoChanger
static autoPtr< fvMeshTopoChangerNew (fvMesh &, const dictionary &dict)
 Select, construct and return the fvMeshTopoChanger. More...
 
static autoPtr< fvMeshTopoChangerNew (fvMesh &)
 Select, construct and return the fvMeshTopoChanger. More...
 

Detailed Description

Dynamic mesh refinement/unrefinement based on volScalarField values.

Refinement can optionally be specified in a cellZone or in multiple regions, each controlled by a different volScalarField.

Usage
Example of single field based refinement in all cells:
topoChanger
{
    type            refiner;

    libs            ("libfvMeshTopoChangers.so");

    // How often to refine
    refineInterval  1;

    // Field to be refinement on
    field           alpha.water;

    // Refine field in between lower..upper
    lowerRefineLevel 0.001;
    upperRefineLevel 0.999;

    // Have slower than 2:1 refinement
    nBufferLayers   1;

    // Refine cells only up to maxRefinement levels
    maxRefinement   2;

    // Stop refinement if maxCells reached
    maxCells        200000;

    // Flux field and corresponding velocity field. Fluxes on changed
    // faces get recalculated by interpolating the velocity. Use 'none'
    // on surfaceScalarFields that do not need to be reinterpolated.
    correctFluxes
    (
        (phi none)
        (nHatf none)
        (rhoPhi none)
        (alphaPhi.water none)
        (ghf none)
    );

    // Write the refinement level as a volScalarField
    dumpLevel       true;
}

Example of single field based refinement in two regions:

topoChanger
{
    type            refiner;

    // How often to refine
    refineInterval  1;

    refinementRegions
    {
        region1
        {
            cellZone        refinementRegion1;

            // Field to be refinement on
            field           alpha.water;

            // Refine field in between lower..upper
            lowerRefineLevel 0.001;
            upperRefineLevel 0.999;

            // Refine cells only up to maxRefinement levels
            maxRefinement   1;

            // If value < unrefineLevel unrefine
            unrefineLevel   10;
        }

        region2
        {
            cellZone        refinementRegion2;

            // Field to be refinement on
            field           alpha.water;

            // Refine field in between lower..upper
            lowerRefineLevel 0.001;
            upperRefineLevel 0.999;

            // Refine cells only up to maxRefinement levels
            maxRefinement   2;

            // If value < unrefineLevel unrefine
            unrefineLevel   10;
        }
    }

    // If value < unrefineLevel (default=great) unrefine
    // unrefineLevel   10;
    // Have slower than 2:1 refinement
    nBufferLayers   1;

    // Stop refinement if maxCells reached
    maxCells        200000;

    // Flux field and corresponding velocity field. Fluxes on changed
    // faces get recalculated by interpolating the velocity. Use 'none'
    // on surfaceScalarFields that do not need to be reinterpolated.
    correctFluxes
    (
        (phi none)
        (nHatf none)
        (rhoPhi none)
        (alphaPhi.water none)
        (ghf none)
    );

    // Write the refinement level as a volScalarField
    dumpLevel       true;
}
Source files

Definition at line 175 of file fvMeshTopoChangersRefiner.H.

Constructor & Destructor Documentation

◆ refiner() [1/2]

◆ refiner() [2/2]

refiner ( const refiner )
delete

Disallow default bitwise copy construction.

◆ ~refiner()

~refiner ( )
virtual

Destructor.

Definition at line 1475 of file fvMeshTopoChangersRefiner.C.

Member Function Documentation

◆ TypeName()

TypeName ( "refiner"  )

Runtime type information.

◆ meshCutter()

const hexRef8& meshCutter ( ) const
inline

Direct access to the refinement engine.

Definition at line 358 of file fvMeshTopoChangersRefiner.H.

Referenced by refiner::update().

Here is the caller graph for this function:

◆ protectedCell() [1/2]

const PackedBoolList& protectedCell ( ) const
inline

Cells which should not be refined/unrefined.

Definition at line 364 of file fvMeshTopoChangersRefiner.H.

◆ protectedCell() [2/2]

PackedBoolList& protectedCell ( )
inline

Cells which should not be refined/unrefined.

Definition at line 370 of file fvMeshTopoChangersRefiner.H.

References refiner::distribute(), refiner::mapMesh(), refiner::operator=(), refiner::topoChange(), refiner::update(), and refiner::write().

Here is the call graph for this function:

◆ update()

bool update ( )
virtual

◆ topoChange()

void topoChange ( const polyTopoChangeMap map)
virtual

Update corresponding to the given map.

Implements fvMeshTopoChanger.

Definition at line 1659 of file fvMeshTopoChangersRefiner.C.

References hexRef8::topoChange().

Referenced by refiner::protectedCell().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ mapMesh()

void mapMesh ( const polyMeshMap map)
virtual

Update from another mesh using the given map.

Implements fvMeshTopoChanger.

Definition at line 1666 of file fvMeshTopoChangersRefiner.C.

References refiner::distribute(), and NotImplemented.

Referenced by refiner::protectedCell().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ distribute()

void distribute ( const polyDistributionMap map)
virtual

Update corresponding to the given distribution map.

Implements fvMeshTopoChanger.

Definition at line 1676 of file fvMeshTopoChangersRefiner.C.

References hexRef8::distribute().

Referenced by refiner::mapMesh(), and refiner::protectedCell().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ write()

bool write ( const bool  write = true) const
virtual

Write using given format, version and compression.

Reimplemented from fvMeshTopoChanger.

Definition at line 1685 of file fvMeshTopoChangersRefiner.C.

References IOobject::AUTO_WRITE, hexRef8::cellLevel(), Foam::dimless, forAll, fvMeshTopoChanger::mesh(), IOobject::NO_READ, fvMesh::time(), timeName, Time::timeName(), regIOobject::write(), and hexRef8::write().

Referenced by refiner::protectedCell().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=()

void operator= ( const refiner )
delete

Disallow default bitwise assignment.

Referenced by refiner::protectedCell().

Here is the caller graph for this function:

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