interfaceTurbulenceDamping Class Reference

Free-surface phase turbulence damping function. More...

Inheritance diagram for interfaceTurbulenceDamping:
Collaboration diagram for interfaceTurbulenceDamping:

Public Member Functions

 TypeName ("interfaceTurbulenceDamping")
 Runtime type information. More...
 
 interfaceTurbulenceDamping (const word &sourceName, const word &modelType, const fvMesh &mesh, const dictionary &dict)
 Construct from explicit source name and mesh. More...
 
 interfaceTurbulenceDamping (const interfaceTurbulenceDamping &)=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 &field, fvMatrix< scalar > &eqn) const
 Add source to mixture epsilon or omega equation. More...
 
virtual void addSup (const volScalarField &rho, const volScalarField &field, fvMatrix< scalar > &eqn) const
 Add source to compressible mixture epsilon or omega equation. More...
 
virtual void addSup (const volScalarField &alpha, const volScalarField &rho, const volScalarField &field, fvMatrix< scalar > &eqn) const
 Add source 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 interfaceTurbulenceDamping &)=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...
 
virtual void addSup (fvMatrix< scalar > &eqn) const
 Add a source term to a field-less proxy equation. More...
 
template<class Type >
tmp< fvMatrix< Type > > sourceProxy (const VolField< Type > &eqnField) const
 Add a source term to an 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 > > sourceProxy (const VolField< Type > &field, const VolField< Type > &eqnField) const
 Return source for an equation. 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 > > sourceProxy (const volScalarField &rho, const VolField< Type > &field, const VolField< Type > &eqnField) const
 Return source for a compressible equation. 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 > > sourceProxy (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field, const VolField< Type > &eqnField) const
 Return source for a phase equation. 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...
 
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...
 
virtual bool write (const bool write=true) const
 Write fvModel data. More...
 
template<class AlphaRhoFieldType , class ... AlphaRhoFieldTypes>
Foam::dimensionSet sourceDims (const dimensionSet &ds, const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields)
 
template<class AlphaRhoFieldType , class ... AlphaRhoFieldTypes>
const Foam::wordfieldName (const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields)
 
template<class AlphaRhoFieldType >
const Foam::wordfieldName (const AlphaRhoFieldType &alphaRhoField)
 
template<class Type , class ... AlphaRhoFieldTypes>
Foam::tmp< Foam::fvMatrix< Type > > sourceTerm (const VolField< Type > &eqnField, const dimensionSet &ds, const AlphaRhoFieldTypes &... alphaRhoFields) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > sourceProxy (const VolField< Type > &eqnField) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > source (const VolField< Type > &field) const
 
template<class Type >
Foam::tmp< Foam::fvMatrix< Type > > sourceProxy (const VolField< Type > &field, const VolField< Type > &eqnField) 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 > > sourceProxy (const volScalarField &rho, const VolField< Type > &field, const VolField< Type > &eqnField) 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 > > sourceProxy (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field, const VolField< Type > &eqnField) 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
 

Additional Inherited Members

- Static Public Member Functions inherited from fvModel
template<class AlphaRhoFieldType , class ... AlphaRhoFieldTypes>
static dimensionSet sourceDims (const dimensionSet &ds, const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields)
 Return the dimensions of the matrix of a source term. More...
 
static const dimensionSetsourceDims (const dimensionSet &ds)
 Return the dimensions of the matrix of a source term (base. More...
 
template<class AlphaRhoFieldType , class ... AlphaRhoFieldTypes>
static const wordfieldName (const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields)
 Return the name of the field associated with a source term. More...
 
template<class AlphaRhoFieldType >
static const wordfieldName (const AlphaRhoFieldType &alphaRhoField)
 Return the name of the field associated with a source term (base. More...
 
static const wordfieldName ()
 Return the name of the field associated with a source term. Special. 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 (const VolField< Type > &field, fvMatrix< Type > &eqn) const
 Add a source term to an equation. More...
 
template<class Type >
void addSupType (const volScalarField &rho, const VolField< Type > &field, fvMatrix< Type > &eqn) const
 Add a source term to a compressible equation. More...
 
template<class Type >
void addSupType (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field, fvMatrix< Type > &eqn) const
 Add a source term to a phase equation. More...
 
template<class Type , class ... AlphaRhoFieldTypes>
tmp< fvMatrix< Type > > sourceTerm (const VolField< Type > &eqnField, const dimensionSet &ds, const AlphaRhoFieldTypes &... alphaRhoFields) const
 Return a source for an equation. More...
 

Detailed Description

Free-surface phase turbulence damping function.

Adds an extra source term to the mixture or phase epsilon or omega equation to reduce turbulence 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 turbulence 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 turbulence conditions.

Usage
Example usage:
interfaceTurbulenceDamping
{
    type    interfaceTurbulenceDamping;

    libs    ("libmultiphaseEulerFvModels.so");

    phase   water;

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

Definition at line 90 of file interfaceTurbulenceDamping.H.

Constructor & Destructor Documentation

◆ interfaceTurbulenceDamping() [1/2]

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

◆ interfaceTurbulenceDamping() [2/2]

Disallow default bitwise copy construction.

Member Function Documentation

◆ TypeName()

TypeName ( "interfaceTurbulenceDamping"  )

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 237 of file interfaceTurbulenceDamping.C.

◆ addSup() [1/3]

void addSup ( const volScalarField field,
fvMatrix< scalar > &  eqn 
) const
virtual

Add source to mixture epsilon or omega equation.

Definition at line 243 of file interfaceTurbulenceDamping.C.

◆ addSup() [2/3]

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

Add source to compressible mixture epsilon or omega equation.

Definition at line 253 of file interfaceTurbulenceDamping.C.

References rho.

◆ addSup() [3/3]

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

Add source to phase epsilon or omega equation.

Definition at line 264 of file interfaceTurbulenceDamping.C.

References alpha(), Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, IOobject::groupName(), Foam::Info, IOobject::name(), Foam::pow4(), fvMatrix< Type >::psi(), rho, 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 301 of file interfaceTurbulenceDamping.C.

◆ mapMesh()

void mapMesh ( const polyMeshMap map)
virtual

Update from another mesh using the given map.

Implements fvModel.

Definition at line 305 of file interfaceTurbulenceDamping.C.

◆ distribute()

void distribute ( const polyDistributionMap )
virtual

Redistribute or update using the given distribution map.

Implements fvModel.

Definition at line 309 of file interfaceTurbulenceDamping.C.

◆ movePoints()

bool movePoints ( )
virtual

Update for mesh motion.

Implements fvModel.

Definition at line 316 of file interfaceTurbulenceDamping.C.

◆ operator=()

void operator= ( const interfaceTurbulenceDamping )
delete

Disallow default bitwise assignment.


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