Loading...
Searching...
No Matches
SpalartAllmarasDDES< BasicTurbulenceModel > Class Template Reference

SpalartAllmaras DDES turbulence model for incompressible and compressible flows. More...

#include <SpalartAllmarasDDES.H>

Inheritance diagram for SpalartAllmarasDDES< BasicTurbulenceModel >:
Collaboration diagram for SpalartAllmarasDDES< BasicTurbulenceModel >:

Public Types

enum class  shieldingMode { standard , ZDES2020 }
 Shielding modes. More...
typedef BasicTurbulenceModel::alphaField alphaField
typedef BasicTurbulenceModel::rhoField rhoField
typedef BasicTurbulenceModel::transportModel transportModel
Public Types inherited from SpalartAllmarasDES< BasicTurbulenceModel >
typedef BasicTurbulenceModel::alphaField alphaField
typedef BasicTurbulenceModel::rhoField rhoField
typedef BasicTurbulenceModel::transportModel transportModel

Public Member Functions

 TypeName ("SpalartAllmarasDDES")
 Runtime type information.
 SpalartAllmarasDDES (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.
virtual ~SpalartAllmarasDDES ()=default
 Destructor.
virtual bool read ()
 Read from dictionary.
virtual tmp< volScalarFieldfd () const
 Return the shielding function.
Public Member Functions inherited from SpalartAllmarasDES< BasicTurbulenceModel >
 TypeName ("SpalartAllmarasDES")
 Runtime type information.
 SpalartAllmarasDES (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.
virtual ~SpalartAllmarasDES ()=default
 Destructor.
virtual tmp< volScalarFieldLESRegion () const
 Return the LES field indicator.

Static Public Attributes

static const Enum< shieldingModeshieldingModeNames
 Shielding mode names.

Protected Member Functions

virtual tmp< volScalarFieldStilda (const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU, const volScalarField &dTilda) const
 Return the production term.
virtual tmp< volScalarFielddTilda (const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
 Return the length scale.
Protected Member Functions inherited from SpalartAllmarasDES< BasicTurbulenceModel >
virtual tmp< volScalarFieldpsi (const volScalarField &chi, const volScalarField &fv1) const
 Return the low Reynolds number correction function.
virtual tmp< volScalarFieldlengthScaleLES (const volScalarField &chi, const volScalarField &fv1) const
 Return the LES length scale.
virtual void correctNut ()

Protected Attributes

shieldingMode shielding_
 Shielding mode.
dimensionedScalar Cd1_
dimensionedScalar Cd2_
dimensionedScalar Cd3_
dimensionedScalar Cd4_
dimensionedScalar betaZDES_
Switch usefP2_
Protected Attributes inherited from SpalartAllmarasDES< BasicTurbulenceModel >
Switch useSigma_
 Switch to activate grey-area enhanced sigma-(D)DES.
dimensionedScalar CDES_
 DES coefficient.
Switch lowReCorrection_
 Flag for low Reynolds number correction.
dimensionedScalar fwStar_

Detailed Description

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

SpalartAllmaras DDES turbulence model for incompressible and compressible flows.

Reference:

    Spalart, P. R., Deck, S., Shur, M. L., Squires,
    K. D., Strelets, M. K., & Travin, A. (2006).
    A new version of detached-eddy simulation,
    resistant to ambiguous grid densities.
    Theoretical and computational fluid dynamics, 20(3), 181-195.
    DOI:10.1007/s00162-006-0015-0

Reference for enhanced shielding function formulation:

    Deck, S., Renard, N. (2020).
    Towards an enhanced protection of attached boundary layers in hybrid
    RANS/LES methods.
    Journal of Computational Physics, 400, 108970.
    DOI:10.1016/j.jcp.2019.108970
Source files

Definition at line 75 of file SpalartAllmarasDDES.H.

Member Typedef Documentation

◆ alphaField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 163 of file SpalartAllmarasDDES.H.

◆ rhoField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 164 of file SpalartAllmarasDDES.H.

◆ transportModel

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 165 of file SpalartAllmarasDDES.H.

Member Enumeration Documentation

◆ shieldingMode

template<class BasicTurbulenceModel>
enum class shieldingMode
strong

Shielding modes.

Enumerator
standard 
ZDES2020 

Definition at line 86 of file SpalartAllmarasDDES.H.

Constructor & Destructor Documentation

◆ SpalartAllmarasDDES()

template<class BasicTurbulenceModel>
SpalartAllmarasDDES ( 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 )

◆ ~SpalartAllmarasDDES()

template<class BasicTurbulenceModel>
virtual ~SpalartAllmarasDDES ( )
virtualdefault

Destructor.

Member Function Documentation

◆ Stilda()

template<class BasicTurbulenceModel>
tmp< volScalarField > Stilda ( const volScalarField & chi,
const volScalarField & fv1,
const volTensorField & gradU,
const volScalarField & dTilda ) const
protectedvirtual

Return the production term.

Reimplemented from SpalartAllmarasDES< BasicTurbulenceModel >.

Definition at line 114 of file SpalartAllmarasDDES.C.

References dTilda(), SpalartAllmarasDES< BasicTurbulenceModel >::lengthScaleLES(), Foam::mag(), Foam::max(), Foam::pos(), Foam::sqr(), Stilda(), and SpalartAllmarasDES< BasicTurbulenceModel >::useSigma_.

Referenced by Stilda().

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

◆ dTilda()

template<class BasicTurbulenceModel>
tmp< volScalarField > dTilda ( const volScalarField & chi,
const volScalarField & fv1,
const volTensorField & gradU ) const
protectedvirtual

Return the length scale.

Reimplemented from SpalartAllmarasDES< BasicTurbulenceModel >.

Definition at line 154 of file SpalartAllmarasDDES.C.

References Foam::dimLength, SpalartAllmarasDES< BasicTurbulenceModel >::lengthScaleLES(), Foam::mag(), Foam::max(), and Foam::Zero.

Referenced by Stilda().

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

◆ TypeName()

template<class BasicTurbulenceModel>
TypeName ( "SpalartAllmarasDDES< BasicTurbulenceModel >" )

Runtime type information.

◆ read()

template<class BasicTurbulenceModel>
bool read ( )
virtual

Read from dictionary.

Reimplemented from SpalartAllmarasDES< BasicTurbulenceModel >.

Definition at line 313 of file SpalartAllmarasDDES.C.

References betaZDES_, Cd1_, Cd2_, Cd3_, Cd4_, SpalartAllmarasDES< BasicTurbulenceModel >::read(), shielding_, shieldingModeNames, and usefP2_.

Here is the call graph for this function:

◆ fd()

template<class BasicTurbulenceModel>
tmp< volScalarField > fd ( ) const
virtual

Return the shielding function.

Definition at line 339 of file SpalartAllmarasDDES.C.

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

Here is the call graph for this function:

Member Data Documentation

◆ shieldingModeNames

template<class BasicTurbulenceModel>
const Foam::Enum< typename Foam::LESModels::SpalartAllmarasDDES< BasicTurbulenceModel >::shieldingMode > shieldingModeNames
static

Shielding mode names.

Definition at line 95 of file SpalartAllmarasDDES.H.

Referenced by read(), and SpalartAllmarasDDES().

◆ shielding_

template<class BasicTurbulenceModel>
shieldingMode shielding_
protected

Shielding mode.

Definition at line 125 of file SpalartAllmarasDDES.H.

Referenced by read(), and SpalartAllmarasDDES().

◆ Cd1_

template<class BasicTurbulenceModel>
dimensionedScalar Cd1_
protected

Definition at line 129 of file SpalartAllmarasDDES.H.

Referenced by read(), and SpalartAllmarasDDES().

◆ Cd2_

template<class BasicTurbulenceModel>
dimensionedScalar Cd2_
protected

Definition at line 130 of file SpalartAllmarasDDES.H.

Referenced by read(), and SpalartAllmarasDDES().

◆ Cd3_

template<class BasicTurbulenceModel>
dimensionedScalar Cd3_
protected

Definition at line 131 of file SpalartAllmarasDDES.H.

Referenced by read(), and SpalartAllmarasDDES().

◆ Cd4_

template<class BasicTurbulenceModel>
dimensionedScalar Cd4_
protected

Definition at line 132 of file SpalartAllmarasDDES.H.

Referenced by read(), and SpalartAllmarasDDES().

◆ betaZDES_

template<class BasicTurbulenceModel>
dimensionedScalar betaZDES_
protected

Definition at line 133 of file SpalartAllmarasDDES.H.

Referenced by read(), and SpalartAllmarasDDES().

◆ usefP2_

template<class BasicTurbulenceModel>
Switch usefP2_
protected

Definition at line 134 of file SpalartAllmarasDDES.H.

Referenced by read(), and SpalartAllmarasDDES().


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