Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
dynamicKEqn< BasicTurbulenceModel > Class Template Reference

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

Inheritance diagram for dynamicKEqn< BasicTurbulenceModel >:
Inheritance graph
[legend]
Collaboration diagram for dynamicKEqn< BasicTurbulenceModel >:
Collaboration graph
[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 LESModel< BasicTurbulenceModel > ::alphaField alphaField
 
typedef LESModel< BasicTurbulenceModel > ::rhoField rhoField
 
typedef LESModel< BasicTurbulenceModel > ::transportModel transportModel
 
- Public Types inherited from linearViscousStress< LESModel< BasicTurbulenceModel > >
typedef LESModel< BasicTurbulenceModel > ::alphaField alphaField
 
typedef LESModel< BasicTurbulenceModel > ::rhoField rhoField
 
typedef LESModel< BasicTurbulenceModel > ::transportModel transportModel
 
- Public Types inherited from LESModel< BasicTurbulenceModel >
typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 

Public Member Functions

 TypeName ("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)
 
virtual ~dynamicKEqn ()=default
 
virtual bool read ()
 
virtual tmp< volScalarFieldk () const
 
tmp< volScalarFieldDkEff () const
 
virtual void correct ()
 
- 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)
 
virtual ~LESeddyViscosity ()=default
 
- 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)
 
virtual ~eddyViscosity ()=default
 
virtual tmp< volScalarFieldnut () const
 
virtual tmp< scalarFieldnut (const label patchi) const
 
virtual tmp< volScalarFieldk () const=0
 
virtual tmp< volSymmTensorFieldR () const
 
virtual void validate ()
 
- Public Member Functions inherited from linearViscousStress< LESModel< 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)
 
virtual ~linearViscousStress ()=default
 
virtual tmp< volSymmTensorFielddevRhoReff () const
 
virtual tmp< volSymmTensorFielddevRhoReff (const volVectorField &U) const
 
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 
virtual tmp< fvVectorMatrixdivDevRhoReff (const volScalarField &rho, volVectorField &U) const
 
- Public Member Functions inherited from LESModel< BasicTurbulenceModel >
 TypeName ("LES")
 
 declareRunTimeSelectionTable (autoPtr, LESModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName),(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName))
 
 LESModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
 
virtual ~LESModel ()=default
 
virtual const dictionarycoeffDict () const
 
const dimensionedScalarCe () const noexcept
 
const dimensionedScalarkMin () const
 
dimensionedScalarkMin ()
 
const volScalarFielddelta () const
 
virtual tmp< volScalarFieldnuEff () const
 
virtual tmp< scalarFieldnuEff (const label patchi) const
 
virtual tmp< volScalarFieldepsilon () const
 
virtual tmp< volScalarFieldomega () const
 

Protected Member Functions

volScalarField Ck (const volSymmTensorField &D, const volScalarField &KK) const
 
volScalarField Ce (const volSymmTensorField &D, const volScalarField &KK) const
 
volScalarField Ce () const
 
void correctNut (const volSymmTensorField &D, const volScalarField &KK)
 
virtual void correctNut ()
 
virtual tmp< fvScalarMatrixkSource () const
 
- Protected Member Functions inherited from LESModel< BasicTurbulenceModel >
virtual void printCoeffs (const word &type)
 
 LESModel (const LESModel &)=delete
 
void operator= (const LESModel &)=delete
 

Protected Attributes

volScalarField k_
 
simpleFilter simpleFilter_
 
autoPtr< LESfilterfilterPtr_
 
LESfilterfilter_
 
- Protected Attributes inherited from eddyViscosity< LESModel< BasicTurbulenceModel > >
volScalarField nut_
 
- Protected Attributes inherited from LESModel< BasicTurbulenceModel >
dictionary LESDict_
 
Switch turbulence_
 
Switch printCoeffs_
 
dictionary coeffDict_
 
dimensionedScalar Ce_
 
dimensionedScalar kMin_
 
dimensionedScalar epsilonMin_
 
dimensionedScalar omegaMin_
 
autoPtr< Foam::LESdeltadelta_
 

Additional Inherited Members

- Static Public Member Functions inherited from LESModel< BasicTurbulenceModel >
static autoPtr< LESModelNew (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName)
 

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

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 136 of file dynamicKEqn.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 137 of file dynamicKEqn.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 138 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 
)

Definition at line 150 of file dynamicKEqn.C.

References Foam::bound(), and Foam::type().

Here is the call graph for this function:

◆ ~dynamicKEqn()

virtual ~dynamicKEqn ( )
virtualdefault

Member Function Documentation

◆ Ck()

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

Definition at line 36 of file dynamicKEqn.C.

References D, delta, Foam::dev(), 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

Definition at line 73 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 90 of file dynamicKEqn.C.

References D, Foam::dev(), Foam::fvc::grad(), Foam::magSqr(), GeometricField::max(), and Foam::symm().

Here is the call graph for this function:

◆ correctNut() [1/2]

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

Definition at line 106 of file dynamicKEqn.C.

References optionList::correct(), D, delta, options::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 120 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 132 of file dynamicKEqn.C.

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

◆ TypeName()

TypeName ( "dynamicKEqn< BasicTurbulenceModel >"  )

◆ read()

bool read
virtual

Reimplemented from LESeddyViscosity< BasicTurbulenceModel >.

Definition at line 202 of file dynamicKEqn.C.

◆ k()

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Definition at line 171 of file dynamicKEqn.H.

References dynamicKEqn< BasicTurbulenceModel >::k_.

◆ DkEff()

tmp<volScalarField> DkEff ( ) const
inline

Definition at line 177 of file dynamicKEqn.H.

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

◆ correct()

void correct
virtual

Member Data Documentation

◆ k_

volScalarField k_
protected

Definition at line 94 of file dynamicKEqn.H.

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

◆ simpleFilter_

simpleFilter simpleFilter_
protected

Definition at line 99 of file dynamicKEqn.H.

◆ filterPtr_

autoPtr<LESfilter> filterPtr_
protected

Definition at line 100 of file dynamicKEqn.H.

◆ filter_

LESfilter& filter_
protected

Definition at line 101 of file dynamicKEqn.H.


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