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

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

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

Public Types

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 ("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 ()=default
 Destructor. More...
 
virtual bool read ()
 Re-read model coefficients if they have changed. More...
 
virtual tmp< volScalarFieldlengthScaleRAS () const
 RAS length scale. More...
 
virtual tmp< volScalarFieldlengthScaleLES (const volScalarField &CDES) const
 LES length scale. More...
 
virtual tmp< volScalarFieldLESRegion () const
 Return the LES field indicator. More...
 
- Public Member Functions inherited from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >
 kOmegaSSTBase (const kOmegaSSTBase &)=delete
 No copy construct. More...
 
 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...
 
void operator= (const kOmegaSSTBase &)=delete
 No copy assignment. More...
 
virtual ~kOmegaSSTBase ()=default
 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< volScalarFieldomega () const
 Return the turbulence kinetic energy dissipation rate. More...
 
virtual void correct ()
 Solve the turbulence equations and correct the turbulence viscosity. More...
 

Protected Member Functions

virtual tmp< volScalarFieldCDES (const volScalarField &F1) const
 Blending for CDES parameter. More...
 
virtual void correctNut (const volScalarField &S2)
 
virtual void correctNut ()
 
tmp< volScalarFieldr (const volScalarField &nur, const volScalarField &magGradU) const
 
virtual tmp< volScalarFieldS2 (const volTensorField &gradU) const
 Return square of strain rate. More...
 
virtual tmp< volScalarFielddTilda (const volScalarField &magGradU, const volScalarField &CDES) const
 Return length scale. More...
 
virtual tmp< volScalarField::InternalepsilonByk (const volScalarField &F1, const volTensorField &gradU) const
 Return epsilon/k. More...
 
virtual tmp< volScalarField::InternalGbyNu0 (const volTensorField &gradU, const volScalarField &S2) const
 Return (G/nu)_0. More...
 
virtual tmp< volScalarField::InternalGbyNu (const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
 Return G/nu. More...
 
- Protected Member Functions inherited from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >
void setDecayControl (const dictionary &dict)
 Set decay control with kInf and omegaInf. More...
 
virtual tmp< volScalarFieldF1 (const volScalarField &CDkOmega) const
 
virtual tmp< volScalarFieldF2 () const
 
virtual tmp< volScalarFieldF3 () const
 
virtual tmp< volScalarFieldF23 () const
 
tmp< volScalarFieldblend (const volScalarField &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
 Return the blended field. More...
 
tmp< volScalarField::Internalblend (const volScalarField::Internal &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
 Return the internal blended field. More...
 
tmp< volScalarFieldalphaK (const volScalarField &F1) const
 
tmp< volScalarFieldalphaOmega (const volScalarField &F1) const
 
tmp< volScalarField::Internalbeta (const volScalarField::Internal &F1) const
 
tmp< volScalarField::Internalgamma (const volScalarField::Internal &F1) const
 
virtual tmp< volScalarField::InternalPk (const volScalarField::Internal &G) const
 Return k production rate. More...
 
virtual tmp< volScalarField::InternalGbyNu (const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
 Return G/nu. More...
 
virtual tmp< fvScalarMatrixkSource () const
 
virtual tmp< fvScalarMatrixomegaSource () const
 
virtual tmp< fvScalarMatrixQsas (const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
 

Protected Attributes

Switch useSigma_
 Switch to activate grey-area enhanced sigma-(D)DES. More...
 
dimensionedScalar kappa_
 
dimensionedScalar CDESkom_
 DES coefficients. More...
 
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_
 Flag to include the F3 term. More...
 
const volScalarFieldy_
 Wall distance. More...
 
volScalarField k_
 Turbulent kinetic energy field [m^2/s^2]. More...
 
volScalarField omega_
 Specific dissipation rate field [1/s]. More...
 
Switch decayControl_
 Flag to include the decay control. More...
 
dimensionedScalar kInf_
 
dimensionedScalar omegaInf_
 

Detailed Description

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

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

Reference:

    Strelets, M. (2001).
    Detached Eddy Simulation of Massively Separated Flows.
    39th AIAA Aerospace Sciences Meeting and Exhibit, Reno, NV.
Note
The default values of the DES constants implemented are code-specific values calibrated for OpenFOAM using decaying isotropic turbulence, and hence deviate slightly from the values suggested in the reference publication.
Source files

Definition at line 67 of file kOmegaSSTDES.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 171 of file kOmegaSSTDES.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 172 of file kOmegaSSTDES.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 173 of file kOmegaSSTDES.H.

Constructor & Destructor Documentation

◆ kOmegaSSTDES()

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.

Definition at line 163 of file kOmegaSSTDES.C.

References Foam::endl(), dimensioned< Type >::getOrAddToDict(), Foam::type(), and WarningInFunction.

Here is the call graph for this function:

◆ ~kOmegaSSTDES()

virtual ~kOmegaSSTDES ( )
virtualdefault

Destructor.

Member Function Documentation

◆ CDES()

virtual tmp<volScalarField> CDES ( const volScalarField F1) const
inlineprotectedvirtual

Blending for CDES parameter.

Definition at line 109 of file kOmegaSSTDES.H.

References kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::blend(), kOmegaSSTDES< BasicTurbulenceModel >::CDESkeps_, and kOmegaSSTDES< BasicTurbulenceModel >::CDESkom_.

Here is the call graph for this function:

◆ correctNut() [1/2]

void correctNut ( const volScalarField S2)
protectedvirtual

Reimplemented from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >.

Definition at line 35 of file kOmegaSSTDES.C.

◆ correctNut() [2/2]

void correctNut ( )
protectedvirtual

Reimplemented from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >.

Definition at line 46 of file kOmegaSSTDES.C.

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

Here is the call graph for this function:

◆ r()

tmp< volScalarField > r ( const volScalarField nur,
const volScalarField magGradU 
) const
protected

Definition at line 54 of file kOmegaSSTDES.C.

References DimensionedField< Type, GeoMesh >::dimensions(), Foam::max(), Foam::min(), Foam::sqr(), and Foam::tr().

Here is the call graph for this function:

◆ S2()

tmp< volScalarField > S2 ( const volTensorField gradU) const
protectedvirtual

Return square of strain rate.

Reimplemented from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >.

Reimplemented in kOmegaSSTDDES< BasicTurbulenceModel >.

Definition at line 72 of file kOmegaSSTDES.C.

References F1, Foam::fvc::grad(), k, Foam::mag(), Foam::pos(), tmp< T >::ref(), and Foam::sqr().

Here is the call graph for this function:

◆ dTilda()

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

Return length scale.

Reimplemented in kOmegaSSTIDDES< BasicTurbulenceModel >, and kOmegaSSTDDES< BasicTurbulenceModel >.

Definition at line 109 of file kOmegaSSTDES.C.

References Foam::min().

Here is the call graph for this function:

◆ epsilonByk()

tmp< volScalarField::Internal > epsilonByk ( const volScalarField F1,
const volTensorField gradU 
) const
protectedvirtual

Return epsilon/k.

Reimplemented from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >.

Definition at line 120 of file kOmegaSSTDES.C.

References F1, Foam::mag(), and Foam::sqrt().

Here is the call graph for this function:

◆ GbyNu0()

tmp< volScalarField::Internal > GbyNu0 ( const volTensorField gradU,
const volScalarField S2 
) const
protectedvirtual

Return (G/nu)_0.

Reimplemented from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >.

Reimplemented in kOmegaSSTDDES< BasicTurbulenceModel >.

Definition at line 132 of file kOmegaSSTDES.C.

◆ GbyNu()

tmp< volScalarField::Internal > GbyNu ( const volScalarField::Internal GbyNu0,
const volScalarField::Internal F2,
const volScalarField::Internal S2 
) const
protectedvirtual

Return G/nu.

Reimplemented in kOmegaSSTDDES< BasicTurbulenceModel >.

Definition at line 149 of file kOmegaSSTDES.C.

◆ TypeName()

TypeName ( "kOmegaSSTDES< BasicTurbulenceModel >"  )

Runtime type information.

◆ read()

bool read ( )
virtual

Re-read model coefficients if they have changed.

Reimplemented from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >.

Reimplemented in kOmegaSSTDDES< BasicTurbulenceModel >, and kOmegaSSTIDDES< BasicTurbulenceModel >.

Definition at line 242 of file kOmegaSSTDES.C.

References Foam::read().

Here is the call graph for this function:

◆ lengthScaleRAS()

Foam::tmp< Foam::volScalarField > lengthScaleRAS ( ) const
virtual

RAS length scale.

Definition at line 260 of file kOmegaSSTDES.C.

References k, and Foam::sqrt().

Here is the call graph for this function:

◆ lengthScaleLES()

Foam::tmp< Foam::volScalarField > lengthScaleLES ( const volScalarField CDES) const
virtual

LES length scale.

Definition at line 272 of file kOmegaSSTDES.C.

References delta.

◆ LESRegion()

tmp< volScalarField > LESRegion ( ) const
virtual

Return the LES field indicator.

Definition at line 281 of file kOmegaSSTDES.C.

References F1, Foam::fvc::grad(), k, Foam::mag(), Foam::neg(), tmp< T >::New(), IOobject::scopedName(), and U.

Here is the call graph for this function:

Member Data Documentation

◆ useSigma_

Switch useSigma_
protected

Switch to activate grey-area enhanced sigma-(D)DES.

Definition at line 91 of file kOmegaSSTDES.H.

◆ kappa_

dimensionedScalar kappa_
protected

Definition at line 95 of file kOmegaSSTDES.H.

◆ CDESkom_

dimensionedScalar CDESkom_
protected

DES coefficients.

Definition at line 100 of file kOmegaSSTDES.H.

Referenced by kOmegaSSTDES< BasicTurbulenceModel >::CDES().

◆ CDESkeps_

dimensionedScalar CDESkeps_
protected

Definition at line 101 of file kOmegaSSTDES.H.

Referenced by kOmegaSSTDES< BasicTurbulenceModel >::CDES().


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