VoFTurbulenceDamping Class Reference

Free-surface momentumTransport damping function. More...

Inheritance diagram for VoFTurbulenceDamping:
Collaboration diagram for VoFTurbulenceDamping:

Public Member Functions

 TypeName ("compressible::VoFTurbulenceDamping")
 Runtime type information. More...
 
 VoFTurbulenceDamping (const word &sourceName, const word &modelType, const fvMesh &mesh, const dictionary &dict)
 Construct from explicit source name and mesh. More...
 
 VoFTurbulenceDamping (const VoFTurbulenceDamping &)=delete
 Disallow default bitwise copy construction. More...
 
virtual wordList addSupFields () const
 Return the list of fields for which the option adds source term. More...
 
virtual void addSup (const volScalarField &rho, fvMatrix< scalar > &eqn, const word &fieldName) const
 Add explicit contribution to epsilon or omega equation. More...
 
virtual void addSup (const volScalarField &alpha, const volScalarField &rho, fvMatrix< scalar > &eqn, const word &fieldName) const
 Add explicit contribution to phase epsilon or omega equation. More...
 
virtual void topoChange (const polyTopoChangeMap &)
 Update topology using the given map. More...
 
virtual void mapMesh (const polyMeshMap &)
 Update from another mesh using the given map. More...
 
virtual void distribute (const polyDistributionMap &)
 Redistribute or update using the given distribution map. More...
 
virtual bool movePoints ()
 Update for mesh motion. More...
 
void operator= (const VoFTurbulenceDamping &)=delete
 Disallow default bitwise assignment. More...
 
- Public Member Functions inherited from fvModel
 TypeName ("fvModel")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, fvModel, dictionary,(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict),(name, modelType, mesh, dict))
 
 fvModel (const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
 Construct from components. More...
 
autoPtr< fvModelclone () const
 Return clone. More...
 
virtual ~fvModel ()
 Destructor. More...
 
const wordname () const
 Return const access to the source name. More...
 
const fvMeshmesh () const
 Return const access to the mesh database. More...
 
const dictionarycoeffs () const
 Return dictionary. More...
 
virtual bool addsSupToField (const word &fieldName) const
 Return true if the fvModel adds a source term to the given. More...
 
virtual scalar maxDeltaT () const
 Return the maximum time-step for stable operation. More...
 
 FOR_ALL_FIELD_TYPES (DEFINE_FV_MODEL_ADD_SUP)
 Add a source term to an equation. More...
 
 FOR_ALL_FIELD_TYPES (DEFINE_FV_MODEL_ADD_RHO_SUP)
 Add a source term to a compressible equation. More...
 
 FOR_ALL_FIELD_TYPES (DEFINE_FV_MODEL_ADD_ALPHA_RHO_SUP)
 Add a source term to a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const VolField< Type > &field) const
 Return source for an equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const VolField< Type > &field, const word &fieldName) const
 Return source for an equation with a specified name. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const volScalarField &rho, const VolField< Type > &field) const
 Return source for a compressible equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const volScalarField &rho, const VolField< Type > &field, const word &fieldName) const
 Return source for a compressible equation with a specified name. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field) const
 Return source for a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field, const word &fieldName) const
 Return source for a phase equation with a specified name. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const volScalarField &alpha, const geometricOneField &rho, const VolField< Type > &field) const
 Return source for a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const geometricOneField &alpha, const volScalarField &rho, const VolField< Type > &field) const
 Return source for a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > source (const geometricOneField &alpha, const geometricOneField &rho, const VolField< Type > &field) const
 Return source for a phase equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > d2dt2 (const VolField< Type > &field) const
 Return source for an equation with a second time derivative. More...
 
template<class Type >
tmp< fvMatrix< Type > > d2dt2 (const VolField< Type > &field, const word &fieldName) const
 Return source for an equation with a second time derivative. More...
 
virtual void preUpdateMesh ()
 Prepare for mesh update. More...
 
virtual void correct ()
 Correct the fvModel. More...
 
virtual bool read (const dictionary &dict)
 Read source dictionary. More...
 
template<class Type , class AlphaRhoFieldType , class ... AlphaRhoFieldTypes>
Foam::dimensionSet sourceDims (const VolField< Type > &field, const dimensionSet &ds, const AlphaRhoFieldType &alphaRho, const AlphaRhoFieldTypes &... alphaRhos)
 
template<class Type >
Foam::dimensionSet sourceDims (const VolField< Type > &field, const dimensionSet &ds)
 
template<class Type , class ... AlphaRhoFieldTypes>
Foam::tmp< Foam::fvMatrix< Type > > source (const VolField< Type > &field, const word &fieldName, const dimensionSet &ds, const AlphaRhoFieldTypes &... alphaRhos) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const VolField< Type > &field, const word &fieldName) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const volScalarField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const volScalarField &rho, const VolField< Type > &field, const word &fieldName) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field, const word &fieldName) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const geometricOneField &alpha, const geometricOneField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const volScalarField &alpha, const geometricOneField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const geometricOneField &alpha, const volScalarField &rho, const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > d2dt2 (const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > d2dt2 (const VolField< Type > &field, const word &fieldName) const
 

Additional Inherited Members

- Static Public Member Functions inherited from fvModel
template<class Type , class AlphaRhoFieldType , class ... AlphaRhoFieldTypes>
static dimensionSet sourceDims (const VolField< Type > &field, const dimensionSet &ds, const AlphaRhoFieldType &alphaRho, const AlphaRhoFieldTypes &... alphaRhos)
 Return the dimensions of the matrix of a source term. More...
 
template<class Type >
static dimensionSet sourceDims (const VolField< Type > &field, const dimensionSet &ds)
 Return the dimensions of the matrix of a source term (base. More...
 
static autoPtr< fvModelNew (const word &name, const fvMesh &mesh, const dictionary &dict)
 Return a reference to the selected fvModel. More...
 
- Protected Member Functions inherited from fvModel
template<class Type >
void addSupType (fvMatrix< Type > &eqn, const word &fieldName) const
 Add a source term to an equation. More...
 
template<class Type >
void addSupType (const volScalarField &rho, fvMatrix< Type > &eqn, const word &fieldName) const
 Add a source term to a compressible equation. More...
 
template<class Type >
void addSupType (const volScalarField &alpha, const volScalarField &rho, fvMatrix< Type > &eqn, const word &fieldName) const
 Add a source term to a phase equation. More...
 
template<class Type , class ... AlphaRhoFieldTypes>
tmp< fvMatrix< Type > > source (const VolField< Type > &field, const word &fieldName, const dimensionSet &ds, const AlphaRhoFieldTypes &... alphaRhos) const
 Return source for equation with specified name and dimensions. More...
 

Detailed Description

Free-surface momentumTransport damping function.

Adds an extra source term to the mixture or phase epsilon or omega equation to reduce momentumTransport generated near a free-surface. The implementation is based on

Reference:

    Frederix, E. M. A., Mathur, A., Dovizio, D., Geurts, B. J.,
    & Komen, E. M. J. (2018).
    Reynolds-averaged modeling of momentumTransport damping
    near a large-scale interface in two-phase flow.
    Nuclear engineering and design, 333, 122-130.

but with an improved formulation for the coefficient A appropriate for unstructured meshes including those with split-cell refinement patterns. However the dimensioned length-scale coefficient delta remains and must be set appropriately for the case by performing test runs and comparing with known results. Clearly this model is far from general and more research is needed in order that delta can be obtained directly from the interface flow and momentumTransport conditions.

Usage
Example usage:
VoFTurbulenceDamping
{
    type    VoFTurbulenceDamping;

    // Interface momentumTransport damping length scale
    // This is a required input as described in section 3.3 of the paper
    delta   1e-4;

    // phase   water; // Optional phase name
}
Source files

Definition at line 93 of file VoFTurbulenceDamping.H.

Constructor & Destructor Documentation

◆ VoFTurbulenceDamping() [1/2]

VoFTurbulenceDamping ( const word sourceName,
const word modelType,
const fvMesh mesh,
const dictionary dict 
)

◆ VoFTurbulenceDamping() [2/2]

Disallow default bitwise copy construction.

Member Function Documentation

◆ TypeName()

TypeName ( "compressible::VoFTurbulenceDamping"  )

Runtime type information.

◆ addSupFields()

Foam::wordList addSupFields ( ) const
virtual

Return the list of fields for which the option adds source term.

to the transport equation

Reimplemented from fvModel.

Definition at line 121 of file VoFTurbulenceDamping.C.

◆ addSup() [1/2]

void addSup ( const volScalarField rho,
fvMatrix< scalar > &  eqn,
const word fieldName 
) const
virtual

Add explicit contribution to epsilon or omega equation.

Definition at line 127 of file VoFTurbulenceDamping.C.

References Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, Foam::Info, Foam::pow4(), fvMatrix< Type >::psi(), Foam::sqr(), and Foam::type().

Here is the call graph for this function:

◆ addSup() [2/2]

void addSup ( const volScalarField alpha,
const volScalarField rho,
fvMatrix< scalar > &  eqn,
const word fieldName 
) const
virtual

Add explicit contribution to phase epsilon or omega equation.

Definition at line 164 of file VoFTurbulenceDamping.C.

References alpha(), Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, IOobject::groupName(), Foam::Info, Foam::pow4(), fvMatrix< Type >::psi(), Foam::sqr(), and Foam::type().

Here is the call graph for this function:

◆ topoChange()

void topoChange ( const polyTopoChangeMap )
virtual

Update topology using the given map.

Implements fvModel.

Definition at line 215 of file VoFTurbulenceDamping.C.

◆ mapMesh()

void mapMesh ( const polyMeshMap map)
virtual

Update from another mesh using the given map.

Implements fvModel.

Definition at line 222 of file VoFTurbulenceDamping.C.

◆ distribute()

void distribute ( const polyDistributionMap )
virtual

Redistribute or update using the given distribution map.

Implements fvModel.

Definition at line 229 of file VoFTurbulenceDamping.C.

◆ movePoints()

bool movePoints ( )
virtual

Update for mesh motion.

Implements fvModel.

Definition at line 236 of file VoFTurbulenceDamping.C.

◆ operator=()

void operator= ( const VoFTurbulenceDamping )
delete

Disallow default bitwise assignment.


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