dynamicKEqn< BasicMomentumTransportModel > Class Template Reference

Dynamic one equation eddy-viscosity model. More...

Inheritance diagram for dynamicKEqn< BasicMomentumTransportModel >:
Collaboration diagram for dynamicKEqn< BasicMomentumTransportModel >:

Public Types

typedef BasicMomentumTransportModel::alphaField alphaField
 
typedef BasicMomentumTransportModel::rhoField rhoField
 
typedef BasicMomentumTransportModel::transportModel transportModel
 
- Public Types inherited from LESeddyViscosity< BasicMomentumTransportModel >
typedef BasicMomentumTransportModel::alphaField alphaField
 
typedef BasicMomentumTransportModel::rhoField rhoField
 
typedef BasicMomentumTransportModel::transportModel transportModel
 
- Public Types inherited from eddyViscosity< LESModel< BasicMomentumTransportModel > >
typedef LESModel< BasicMomentumTransportModel > ::alphaField alphaField
 
typedef LESModel< BasicMomentumTransportModel > ::rhoField rhoField
 
typedef LESModel< BasicMomentumTransportModel > ::transportModel transportModel
 
- Public Types inherited from linearViscousStress< LESModel< BasicMomentumTransportModel > >
typedef LESModel< BasicMomentumTransportModel > ::alphaField alphaField
 
typedef LESModel< BasicMomentumTransportModel > ::rhoField rhoField
 
typedef LESModel< BasicMomentumTransportModel > ::transportModel transportModel
 
- Public Types inherited from LESModel< BasicMomentumTransportModel >
typedef BasicMomentumTransportModel::alphaField alphaField
 
typedef BasicMomentumTransportModel::rhoField rhoField
 
typedef BasicMomentumTransportModel::transportModel transportModel
 

Public Member Functions

 TypeName ("dynamicKEqn")
 Runtime type information. More...
 
 dynamicKEqn (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
 Construct from components. More...
 
 dynamicKEqn (const dynamicKEqn &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~dynamicKEqn ()
 Destructor. More...
 
virtual bool read ()
 Read model coefficients if they have changed. More...
 
virtual tmp< volScalarFieldk () const
 Return SGS kinetic energy. More...
 
virtual tmp< volScalarFieldepsilon () const
 Return sub-grid disipation rate. More...
 
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k. More...
 
virtual void correct ()
 Correct Eddy-Viscosity and related properties. More...
 
void operator= (const dynamicKEqn &)=delete
 Disallow default bitwise assignment. More...
 
- Public Member Functions inherited from LESeddyViscosity< BasicMomentumTransportModel >
 LESeddyViscosity (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Construct from components. More...
 
 LESeddyViscosity (const LESeddyViscosity &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~LESeddyViscosity ()
 Destructor. More...
 
void operator= (const LESeddyViscosity &)=delete
 Disallow default bitwise assignment. More...
 
- Public Member Functions inherited from eddyViscosity< LESModel< BasicMomentumTransportModel > >
 eddyViscosity (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Construct from components. More...
 
virtual ~eddyViscosity ()
 Destructor. More...
 
virtual tmp< volScalarFieldnut () const
 Return the turbulence viscosity. More...
 
virtual tmp< scalarFieldnut (const label patchi) const
 Return the turbulence viscosity on patch. More...
 
virtual tmp< volSymmTensorFieldsigma () const
 Return the Reynolds stress tensor [m^2/s^2]. More...
 
virtual void validate ()
 Validate the turbulence fields after construction. More...
 
- Public Member Functions inherited from linearViscousStress< LESModel< BasicMomentumTransportModel > >
 linearViscousStress (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Construct from components. More...
 
virtual ~linearViscousStress ()
 Destructor. More...
 
virtual tmp< volSymmTensorFielddevTau () const
 Return the effective stress tensor. More...
 
virtual tmp< fvVectorMatrixdivDevTau (volVectorField &U) const
 Return the source term for the momentum equation. More...
 
virtual tmp< fvVectorMatrixdivDevTau (const volScalarField &rho, volVectorField &U) const
 Return the source term for the momentum equation. More...
 
- Public Member Functions inherited from LESModel< BasicMomentumTransportModel >
 TypeName ("LES")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, LESModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport),(alpha, rho, U, alphaRhoPhi, phi, transport))
 
 LESModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Construct from components. More...
 
 LESModel (const LESModel &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~LESModel ()
 Destructor. More...
 
virtual const dictionarycoeffDict () const
 Const access to the coefficients dictionary. More...
 
const dimensionedScalarkMin () const
 Return the lower allowable limit for k (default: small) More...
 
dimensionedScalarkMin ()
 Allow kMin to be changed. More...
 
const volScalarFielddelta () const
 Access function to filter width. More...
 
virtual tmp< volScalarFieldnuEff () const
 Return the effective viscosity. More...
 
virtual tmp< scalarFieldnuEff (const label patchi) const
 Return the effective viscosity on patch. More...
 
void operator= (const LESModel &)=delete
 Disallow default bitwise assignment. More...
 

Protected Member Functions

tmp< volScalarFieldKK () const
 Return the KK field obtained by filtering the velocity field U. More...
 
volScalarField Ck (const volSymmTensorField &D, const volScalarField &KK) const
 Calculate Ck by filtering the velocity field U. More...
 
volScalarField Ce (const volSymmTensorField &D, const volScalarField &KK) const
 Calculate Ce by filtering the velocity field U. More...
 
volScalarField Ce () const
 
void correctNut (const volSymmTensorField &D, const volScalarField &KK)
 Update sub-grid eddy-viscosity. More...
 
virtual void correctNut ()
 
virtual tmp< fvScalarMatrixkSource () const
 
- Protected Member Functions inherited from LESModel< BasicMomentumTransportModel >
virtual void printCoeffs (const word &type)
 Print model coefficients. More...
 

Protected Attributes

volScalarField k_
 
simpleFilter simpleFilter_
 
autoPtr< LESfilterfilterPtr_
 
LESfilterfilter_
 
- Protected Attributes inherited from LESeddyViscosity< BasicMomentumTransportModel >
dimensionedScalar Ce_
 
- Protected Attributes inherited from eddyViscosity< LESModel< BasicMomentumTransportModel > >
volScalarField nut_
 
- Protected Attributes inherited from LESModel< BasicMomentumTransportModel >
dictionary LESDict_
 LES coefficients dictionary. More...
 
Switch turbulence_
 Turbulence on/off flag. More...
 
Switch printCoeffs_
 Flag to print the model coeffs at run-time. More...
 
dictionary coeffDict_
 Model coefficients dictionary. More...
 
dimensionedScalar kMin_
 Lower limit of k. More...
 
dimensionedScalar epsilonMin_
 Lower limit of epsilon. More...
 
dimensionedScalar omegaMin_
 Lower limit for omega. More...
 
autoPtr< Foam::LESdeltadelta_
 Run-time selectable delta model. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from LESModel< BasicMomentumTransportModel >
static autoPtr< LESModelNew (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Return a reference to the selected LES model. More...
 

Detailed Description

template<class BasicMomentumTransportModel>
class Foam::LESModels::dynamicKEqn< BasicMomentumTransportModel >

Dynamic one equation eddy-viscosity model.

Eddy viscosity SGS model using a modeled balance equation to simulate the behaviour of k in which a dynamic procedure is applied to evaluate the coefficients.

Reference:

    Kim, W and Menon, S. (1995).
    A new dynamic one-equation subgrid-scale model for
    large eddy simulation.
    In 33rd Aerospace Sciences Meeting and Exhibit, Reno, NV, 1995.

There are no default model coefficients but the filter used for KK must be supplied, e.g.

    dynamicKEqnCoeffs
    {
        filter simple;
    }
Source files

Definition at line 73 of file dynamicKEqn.H.

Member Typedef Documentation

◆ alphaField

typedef BasicMomentumTransportModel::alphaField alphaField

Definition at line 128 of file dynamicKEqn.H.

◆ rhoField

typedef BasicMomentumTransportModel::rhoField rhoField

Definition at line 129 of file dynamicKEqn.H.

◆ transportModel

typedef BasicMomentumTransportModel::transportModel transportModel

Definition at line 130 of file dynamicKEqn.H.

Constructor & Destructor Documentation

◆ dynamicKEqn() [1/2]

dynamicKEqn ( const alphaField alpha,
const rhoField rho,
const volVectorField U,
const surfaceScalarField alphaRhoPhi,
const surfaceScalarField phi,
const transportModel transport,
const word type = typeName 
)

Construct from components.

Definition at line 147 of file dynamicKEqn.C.

References Foam::bound().

Referenced by dynamicKEqn< BasicMomentumTransportModel >::kSource().

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

◆ dynamicKEqn() [2/2]

dynamicKEqn ( const dynamicKEqn< BasicMomentumTransportModel > &  )
delete

Disallow default bitwise copy construction.

◆ ~dynamicKEqn()

virtual ~dynamicKEqn ( )
inlinevirtual

Destructor.

Definition at line 156 of file dynamicKEqn.H.

References dynamicKEqn< BasicMomentumTransportModel >::read().

Here is the call graph for this function:

Member Function Documentation

◆ KK()

Foam::tmp< Foam::volScalarField > KK ( ) const
protected

Return the KK field obtained by filtering the velocity field U.

Definition at line 41 of file dynamicKEqn.C.

References dynamicKEqn< BasicMomentumTransportModel >::Ck(), Foam::dimVelocity, Foam::magSqr(), Foam::max(), and Foam::sqr().

Here is the call graph for this function:

◆ Ck()

volScalarField Ck ( const volSymmTensorField D,
const volScalarField KK 
) const
protected

Calculate Ck by filtering the velocity field U.

Definition at line 53 of file dynamicKEqn.C.

References dynamicKEqn< BasicMomentumTransportModel >::Ce(), delta, Foam::dev(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::mag(), Foam::magSqr(), Foam::sqr(), and Foam::sqrt().

Referenced by dynamicKEqn< BasicMomentumTransportModel >::KK().

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

◆ Ce() [1/2]

volScalarField Ce ( const volSymmTensorField D,
const volScalarField KK 
) const
protected

Calculate Ce by filtering the velocity field U.

Definition at line 84 of file dynamicKEqn.C.

References delta, Foam::mag(), Foam::magSqr(), and Foam::pow().

Here is the call graph for this function:

◆ Ce() [2/2]

volScalarField Ce ( ) const
protected

Definition at line 101 of file dynamicKEqn.C.

References dynamicKEqn< BasicMomentumTransportModel >::correctNut(), Foam::dev(), Foam::fvc::grad(), and Foam::symm().

Referenced by dynamicKEqn< BasicMomentumTransportModel >::Ck().

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

◆ correctNut() [1/2]

void correctNut ( const volSymmTensorField D,
const volScalarField KK 
)
protected

Update sub-grid eddy-viscosity.

Definition at line 110 of file dynamicKEqn.C.

References delta, dictionary::New(), and Foam::sqrt().

Here is the call graph for this function:

◆ correctNut() [2/2]

void correctNut ( )
protectedvirtual

Implements eddyViscosity< LESModel< BasicMomentumTransportModel > >.

Definition at line 122 of file dynamicKEqn.C.

References Foam::fvc::grad(), and Foam::symm().

Referenced by dynamicKEqn< BasicMomentumTransportModel >::Ce().

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

◆ kSource()

tmp< fvScalarMatrix > kSource ( ) const
protectedvirtual

Definition at line 129 of file dynamicKEqn.C.

References Foam::dimTime, Foam::dimVolume, and dynamicKEqn< BasicMomentumTransportModel >::dynamicKEqn().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "dynamicKEqn< BasicMomentumTransportModel >"  )

Runtime type information.

◆ read()

bool read ( )
virtual

Read model coefficients if they have changed.

Reimplemented from LESeddyViscosity< BasicMomentumTransportModel >.

Definition at line 197 of file dynamicKEqn.C.

Referenced by dynamicKEqn< BasicMomentumTransportModel >::~dynamicKEqn().

Here is the caller graph for this function:

◆ k()

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Return SGS kinetic energy.

Implements eddyViscosity< LESModel< BasicMomentumTransportModel > >.

Definition at line 166 of file dynamicKEqn.H.

References dynamicKEqn< BasicMomentumTransportModel >::epsilon(), and dynamicKEqn< BasicMomentumTransportModel >::k_.

Here is the call graph for this function:

◆ epsilon()

tmp< volScalarField > epsilon ( ) const
virtual

Return sub-grid disipation rate.

Reimplemented from LESeddyViscosity< BasicMomentumTransportModel >.

Definition at line 213 of file dynamicKEqn.C.

References delta, IOobject::groupName(), k, GeometricField< scalar, fvPatchField, volMesh >::New(), and Foam::sqrt().

Referenced by dynamicKEqn< BasicMomentumTransportModel >::k().

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

◆ DkEff()

◆ correct()

◆ operator=()

void operator= ( const dynamicKEqn< BasicMomentumTransportModel > &  )
delete

Disallow default bitwise assignment.

Referenced by dynamicKEqn< BasicMomentumTransportModel >::DkEff().

Here is the caller graph for this function:

Member Data Documentation

◆ k_

volScalarField k_
protected

Definition at line 83 of file dynamicKEqn.H.

Referenced by dynamicKEqn< BasicMomentumTransportModel >::k().

◆ simpleFilter_

simpleFilter simpleFilter_
protected

Definition at line 88 of file dynamicKEqn.H.

◆ filterPtr_

autoPtr<LESfilter> filterPtr_
protected

Definition at line 89 of file dynamicKEqn.H.

◆ filter_

LESfilter& filter_
protected

Definition at line 90 of file dynamicKEqn.H.


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