Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions
kOmegaSSTIDDES< BasicTurbulenceModel > Class Template Reference

k-omega-SST IDDES turbulence model for incompressible and compressible flows More...

Inheritance diagram for kOmegaSSTIDDES< BasicTurbulenceModel >:
Inheritance graph
[legend]
Collaboration diagram for kOmegaSSTIDDES< BasicTurbulenceModel >:
Collaboration graph
[legend]

Public Types

typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 
- Public Types inherited from kOmegaSSTDES< BasicTurbulenceModel >
typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 
- Public Types inherited from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >
typedef DESModel< BasicTurbulenceModel > ::alphaField alphaField
 
typedef DESModel< BasicTurbulenceModel > ::rhoField rhoField
 
typedef DESModel< BasicTurbulenceModel > ::transportModel transportModel
 

Public Member Functions

 TypeName ("kOmegaSSTIDDES")
 Runtime type information. More...
 
 kOmegaSSTIDDES (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 ~kOmegaSSTIDDES ()
 Destructor. More...
 
virtual bool read ()
 Re-read model coefficients if they have changed. More...
 
- Public Member Functions inherited from kOmegaSSTDES< BasicTurbulenceModel >
 TypeName ("kOmegaSSTDES")
 Runtime type information. More...
 
 kOmegaSSTDES (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 ~kOmegaSSTDES ()
 Destructor. More...
 
virtual void correct ()
 Solve the turbulence equations and correct the turbulence viscosity. More...
 
virtual tmp< volScalarFieldLESRegion () const
 Return the LES field indicator. More...
 
- Public Member Functions inherited from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >
 kOmegaSSTBase (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 ~kOmegaSSTBase ()
 Destructor. More...
 
tmp< volScalarFieldDkEff (const volScalarField &F1) const
 Return the effective diffusivity for k. More...
 
tmp< volScalarFieldDomegaEff (const volScalarField &F1) const
 Return the effective diffusivity for omega. More...
 
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy. More...
 
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate. More...
 
virtual tmp< volScalarFieldomega () const
 Return the turbulence kinetic energy dissipation rate. More...
 

Protected Member Functions

virtual tmp< volScalarFielddTilda (const volScalarField &magGradU, const volScalarField &CDES) const
 Length scale. More...
 
- Protected Member Functions inherited from kOmegaSSTDES< BasicTurbulenceModel >
virtual tmp< volScalarFieldCDES (const volScalarField &F1) const
 Blending for CDES parameter. More...
 
virtual void correctNut (const volScalarField &S2)
 
virtual void correctNut ()
 
- Protected Member Functions inherited from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >
tmp< volScalarFieldF1 (const volScalarField &CDkOmega) const
 
tmp< volScalarFieldF2 () const
 
tmp< volScalarFieldF3 () const
 
tmp< volScalarFieldF23 () const
 
tmp< volScalarFieldblend (const volScalarField &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
 
tmp< volScalarFieldalphaK (const volScalarField &F1) const
 
tmp< volScalarFieldalphaOmega (const volScalarField &F1) const
 
tmp< volScalarFieldbeta (const volScalarField &F1) const
 
tmp< volScalarFieldgamma (const volScalarField &F1) const
 
virtual tmp< fvScalarMatrixkSource () const
 
virtual tmp< fvScalarMatrixomegaSource () const
 
virtual tmp< fvScalarMatrixQsas (const volScalarField &S2, const volScalarField &gamma, const volScalarField &beta) const
 

Protected Attributes

dimensionedScalar Cdt1_
 
dimensionedScalar Cdt2_
 
dimensionedScalar Cl_
 
dimensionedScalar Ct_
 
const IDDESDeltaIDDESDelta_
 
- Protected Attributes inherited from kOmegaSSTDES< BasicTurbulenceModel >
dimensionedScalar kappa_
 
dimensionedScalar CDESkom_
 
dimensionedScalar CDESkeps_
 
- Protected Attributes inherited from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >
dimensionedScalar alphaK1_
 
dimensionedScalar alphaK2_
 
dimensionedScalar alphaOmega1_
 
dimensionedScalar alphaOmega2_
 
dimensionedScalar gamma1_
 
dimensionedScalar gamma2_
 
dimensionedScalar beta1_
 
dimensionedScalar beta2_
 
dimensionedScalar betaStar_
 
dimensionedScalar a1_
 
dimensionedScalar b1_
 
dimensionedScalar c1_
 
Switch F3_
 
const volScalarFieldy_
 Wall distance. More...
 
volScalarField k_
 
volScalarField omega_
 

Private Member Functions

const IDDESDeltasetDelta () const
 Check that the supplied delta is an IDDESDelta. More...
 
tmp< volScalarFieldalpha () const
 
tmp< volScalarFieldft (const volScalarField &magGradU) const
 
tmp< volScalarFieldfl (const volScalarField &magGradU) const
 
tmp< volScalarFieldrd (const volScalarField &nur, const volScalarField &magGradU) const
 
tmp< volScalarFieldfdt (const volScalarField &magGradU) const
 Delay function. More...
 
 kOmegaSSTIDDES (const kOmegaSSTIDDES &)
 
kOmegaSSTIDDESoperator= (const kOmegaSSTIDDES &)
 

Detailed Description

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

k-omega-SST IDDES turbulence model for incompressible and compressible flows

Reference:

    Gritskevich, M.S., Garbaruk, A.V., Schuetze, J., Menter, F.R. (2011)
    Development of DDES and IDDES Formulations for the k-omega
    Shear Stress Transport Model, Flow, Turbulence and Combustion,
    pp. 1-19
Source files

Definition at line 63 of file kOmegaSSTIDDES.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 117 of file kOmegaSSTIDDES.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 118 of file kOmegaSSTIDDES.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 119 of file kOmegaSSTIDDES.H.

Constructor & Destructor Documentation

◆ kOmegaSSTIDDES() [1/2]

kOmegaSSTIDDES ( const kOmegaSSTIDDES< BasicTurbulenceModel > &  )
private

◆ kOmegaSSTIDDES() [2/2]

kOmegaSSTIDDES ( 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 163 of file kOmegaSSTIDDES.C.

References Foam::type().

Here is the call graph for this function:

◆ ~kOmegaSSTIDDES()

virtual ~kOmegaSSTIDDES ( )
inlinevirtual

Destructor.

Definition at line 143 of file kOmegaSSTIDDES.H.

Member Function Documentation

◆ setDelta()

const IDDESDelta & setDelta
private

Check that the supplied delta is an IDDESDelta.

Definition at line 38 of file kOmegaSSTIDDES.C.

References Foam::exit(), Foam::FatalError, and FatalErrorInFunction.

Here is the call graph for this function:

◆ alpha()

tmp< volScalarField > alpha
private

Definition at line 52 of file kOmegaSSTIDDES.C.

References Foam::max().

Here is the call graph for this function:

◆ ft()

tmp< volScalarField > ft ( const volScalarField magGradU) const
private

Definition at line 60 of file kOmegaSSTIDDES.C.

References Foam::pow3(), Foam::sqr(), and Foam::tanh().

Here is the call graph for this function:

◆ fl()

tmp< volScalarField > fl ( const volScalarField magGradU) const
private

Definition at line 70 of file kOmegaSSTIDDES.C.

References nu, Foam::pow(), Foam::sqr(), and Foam::tanh().

Here is the call graph for this function:

◆ rd()

tmp< volScalarField > rd ( const volScalarField nur,
const volScalarField magGradU 
) const
private

Definition at line 80 of file kOmegaSSTIDDES.C.

References Foam::max(), Foam::min(), Foam::sqr(), and Foam::tr().

Here is the call graph for this function:

◆ fdt()

tmp< volScalarField > fdt ( const volScalarField magGradU) const
private

Delay function.

Definition at line 111 of file kOmegaSSTIDDES.C.

References Foam::pow(), and Foam::tanh().

Here is the call graph for this function:

◆ operator=()

kOmegaSSTIDDES& operator= ( const kOmegaSSTIDDES< BasicTurbulenceModel > &  )
private

◆ dTilda()

tmp< volScalarField > dTilda ( const volScalarField magGradU,
const volScalarField CDES 
) const
protectedvirtual

Length scale.

Reimplemented from kOmegaSSTDES< BasicTurbulenceModel >.

Definition at line 121 of file kOmegaSSTIDDES.C.

References Foam::constant::atomic::alpha, delta, Foam::dimLength, Foam::exp(), k, Foam::max(), Foam::min(), Foam::neg(), Foam::pos(), Foam::pow(), Foam::sqr(), and Foam::sqrt().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "kOmegaSSTIDDES< BasicTurbulenceModel >"  )

Runtime type information.

◆ read()

bool read
virtual

Re-read model coefficients if they have changed.

Reimplemented from kOmegaSSTDES< BasicTurbulenceModel >.

Definition at line 234 of file kOmegaSSTIDDES.C.

Field Documentation

◆ Cdt1_

dimensionedScalar Cdt1_
protected

Definition at line 96 of file kOmegaSSTIDDES.H.

◆ Cdt2_

dimensionedScalar Cdt2_
protected

Definition at line 97 of file kOmegaSSTIDDES.H.

◆ Cl_

dimensionedScalar Cl_
protected

Definition at line 98 of file kOmegaSSTIDDES.H.

◆ Ct_

dimensionedScalar Ct_
protected

Definition at line 99 of file kOmegaSSTIDDES.H.

◆ IDDESDelta_

const IDDESDelta& IDDESDelta_
protected

Definition at line 103 of file kOmegaSSTIDDES.H.


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