dynamicKEqn< BasicTurbulenceModel > Class Template Reference

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

Inheritance diagram for dynamicKEqn< BasicTurbulenceModel >:
[legend]
Collaboration diagram for dynamicKEqn< BasicTurbulenceModel >:
[legend]

Public Types

typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 
- Public Types inherited from LESeddyViscosity< BasicTurbulenceModel >
typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 
- Public Types inherited from eddyViscosity< LESModel< BasicTurbulenceModel > >
typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 
- Public Types inherited from linearViscousStress< BasicTurbulenceModel >
typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::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 &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
 Construct from components. More...
 
virtual ~dynamicKEqn ()=default
 Destructor. More...
 
virtual bool read ()
 Read model coefficients if they have changed. More...
 
virtual tmp< volScalarFieldk () const
 Return SGS kinetic energy. More...
 
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k. More...
 
virtual void correct ()
 Correct Eddy-Viscosity and related properties. More...
 
- Public Member Functions inherited from LESeddyViscosity< BasicTurbulenceModel >
 LESeddyViscosity (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName)
 Construct from components. More...
 
virtual ~LESeddyViscosity ()=default
 Destructor. More...
 
virtual bool read ()
 Read model coefficients if they have changed. More...
 
- Public Member Functions inherited from eddyViscosity< LESModel< BasicTurbulenceModel > >
 eddyViscosity (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
 Construct from components. More...
 
virtual ~eddyViscosity ()=default
 Destructor. More...
 
virtual bool read ()=0
 Re-read model coefficients if they have changed. 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< volScalarFieldk () const=0
 Return the turbulence kinetic energy. More...
 
virtual tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor. More...
 
virtual void validate ()
 Validate the turbulence fields after construction. More...
 
virtual void correct ()=0
 Solve the turbulence equations and correct the turbulence viscosity. More...
 
- Public Member Functions inherited from linearViscousStress< BasicTurbulenceModel >
 linearViscousStress (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
 Construct from components. More...
 
virtual ~linearViscousStress ()=default
 Destructor. More...
 
virtual bool read ()=0
 Re-read model coefficients if they have changed. More...
 
virtual tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor. More...
 
virtual tmp< volSymmTensorFielddevRhoReff (const volVectorField &U) const
 Return the effective stress tensor based on a given velocity field. More...
 
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 Return the source term for the momentum equation. More...
 
virtual tmp< fvVectorMatrixdivDevRhoReff (const volScalarField &rho, volVectorField &U) const
 Return the source term for the momentum equation. More...
 
virtual void correct ()=0
 Solve the turbulence equations and correct the turbulence viscosity. More...
 

Protected Member Functions

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
 
virtual void correctNut ()=0
 

Protected Attributes

volScalarField k_
 
simpleFilter simpleFilter_
 
autoPtr< LESfilterfilterPtr_
 
LESfilterfilter_
 
- Protected Attributes inherited from eddyViscosity< LESModel< BasicTurbulenceModel > >
volScalarField nut_
 

Detailed Description

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

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 79 of file dynamicKEqn.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 140 of file dynamicKEqn.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 141 of file dynamicKEqn.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 142 of file dynamicKEqn.H.

Constructor & Destructor Documentation

◆ dynamicKEqn()

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

Construct from components.

Definition at line 156 of file dynamicKEqn.C.

References Foam::bound(), dynamicKEqn< BasicTurbulenceModel >::k_, and Foam::type().

Here is the call graph for this function:

◆ ~dynamicKEqn()

virtual ~dynamicKEqn ( )
virtualdefault

Destructor.

Member Function Documentation

◆ Ck()

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

Calculate Ck by filtering the velocity field U.

Definition at line 42 of file dynamicKEqn.C.

References D, delta, Foam::dev(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::mag(), Foam::magSqr(), Foam::max(), Foam::sqr(), Foam::sqrt(), and Foam::Zero.

Here is the call 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 79 of file dynamicKEqn.C.

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

Here is the call graph for this function:

◆ Ce() [2/2]

volScalarField Ce
protected

Definition at line 97 of file dynamicKEqn.C.

References D, Foam::dev(), DimensionedField< Type, GeoMesh >::dimensions(), Foam::fvc::grad(), Foam::magSqr(), GeometricField< Type, PatchField, GeoMesh >::max(), and Foam::symm().

Here is the call graph for this function:

◆ correctNut() [1/2]

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

Update sub-grid eddy-viscosity.

Definition at line 112 of file dynamicKEqn.C.

References D, delta, Time::New(), and Foam::sqrt().

Here is the call graph for this function:

◆ correctNut() [2/2]

void correctNut
protectedvirtual

Implements eddyViscosity< LESModel< BasicTurbulenceModel > >.

Definition at line 127 of file dynamicKEqn.C.

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

Here is the call graph for this function:

◆ kSource()

tmp< fvScalarMatrix > kSource
protectedvirtual

Definition at line 139 of file dynamicKEqn.C.

References Foam::dimTime, and Foam::dimVolume.

◆ TypeName()

TypeName ( "dynamicKEqn< BasicTurbulenceModel >"  )

Runtime type information.

◆ read()

bool read
virtual

Read model coefficients if they have changed.

Reimplemented from LESeddyViscosity< BasicTurbulenceModel >.

Definition at line 209 of file dynamicKEqn.C.

◆ k()

virtual tmp< volScalarField > k ( ) const
inlinevirtual

Return SGS kinetic energy.

Implements eddyViscosity< LESModel< BasicTurbulenceModel > >.

Definition at line 175 of file dynamicKEqn.H.

References dynamicKEqn< BasicTurbulenceModel >::k_.

◆ DkEff()

tmp< volScalarField > DkEff ( ) const
inline

Return the effective diffusivity for k.

Definition at line 181 of file dynamicKEqn.H.

References nu, and eddyViscosity< LESModel< BasicTurbulenceModel > >::nut_.

◆ correct()

Member Data Documentation

◆ k_

◆ simpleFilter_

simpleFilter simpleFilter_
protected

Definition at line 103 of file dynamicKEqn.H.

◆ filterPtr_

autoPtr<LESfilter> filterPtr_
protected

Definition at line 104 of file dynamicKEqn.H.

◆ filter_

LESfilter& filter_
protected

Definition at line 105 of file dynamicKEqn.H.


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