kOmegaSSTDES< BasicMomentumTransportModel > Class Template Reference

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

Inheritance diagram for kOmegaSSTDES< BasicMomentumTransportModel >:
Collaboration diagram for kOmegaSSTDES< BasicMomentumTransportModel >:

Public Types

typedef BasicMomentumTransportModel::alphaField alphaField
 
typedef BasicMomentumTransportModel::rhoField rhoField
 
- Public Types inherited from kOmegaSST< LESeddyViscosity< BasicMomentumTransportModel >, BasicMomentumTransportModel >
typedef BasicMomentumTransportModel::alphaField alphaField
 
typedef BasicMomentumTransportModel::rhoField rhoField
 

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 viscosity &viscosity, const word &type=typeName)
 Construct from components. More...
 
virtual ~kOmegaSSTDES ()
 Destructor. More...
 
virtual bool read ()
 Read model coefficients if they have changed. More...
 
- Public Member Functions inherited from kOmegaSST< LESeddyViscosity< BasicMomentumTransportModel >, BasicMomentumTransportModel >
 kOmegaSST (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
 Construct from components. More...
 
 kOmegaSST (const kOmegaSST &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~kOmegaSST ()
 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...
 
virtual void correct ()
 Solve the turbulence equations and correct the turbulence viscosity. More...
 
void operator= (const kOmegaSST &)=delete
 Disallow default bitwise assignment. More...
 

Protected Member Functions

tmp< volScalarField::InternalLt () const
 Return the turbulent length-scale. More...
 
virtual tmp< volScalarField::InternalFDES (const volScalarField::Internal &F1, const volScalarField::Internal &F2) const
 The DES dissipation-rate multiplier with options zonal filtering. More...
 
virtual tmp< volScalarField::InternalepsilonByk (const volScalarField::Internal &F1, const volScalarField::Internal &F2) const
 Return epsilon/k which for standard RAS is betaStar*omega. More...
 
- Protected Member Functions inherited from kOmegaSST< LESeddyViscosity< BasicMomentumTransportModel >, BasicMomentumTransportModel >
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
 
tmp< volScalarField::Internalblend (const volScalarField::Internal &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
 
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 void correctNut (const volScalarField &S2, const volScalarField &F2)
 
virtual void correctNut ()
 
virtual tmp< volScalarField::InternalPk (const volScalarField::Internal &G) const
 Return k production rate. More...
 
virtual tmp< volScalarField::InternalepsilonByk (const volScalarField::Internal &F1, const volScalarField::Internal &F2) const
 Return epsilon/k which for standard RAS is betaStar*omega. 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

dimensionedScalar CDES_
 DES coefficient. More...
 
direction FSST_
 Zonal filter choice. More...
 
- Protected Attributes inherited from kOmegaSST< LESeddyViscosity< BasicMomentumTransportModel >, BasicMomentumTransportModel >
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_
 

Detailed Description

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

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

DES model described in:

    Menter, F. R., Kuntz, M., and Langtry, R. (2003).
    Ten Years of Industrial Experience with the SST Turbulence Model.
    Turbulence, Heat and Mass Transfer 4, ed: K. Hanjalic, Y. Nagano,
    & M. Tummers, Begell House, Inc., 625 - 632.

Optional support for zonal filtering based on F1 or F2 is provided as described in the paper.

For further details of the implementation of the base k-omega-SST model see Foam::kOmegaSST.

See also
Foam::kOmegaSST
Source files

Definition at line 71 of file kOmegaSSTDES.H.

Member Typedef Documentation

◆ alphaField

typedef BasicMomentumTransportModel::alphaField alphaField

Definition at line 120 of file kOmegaSSTDES.H.

◆ rhoField

typedef BasicMomentumTransportModel::rhoField rhoField

Definition at line 121 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 viscosity viscosity,
const word type = typeName 
)

Construct from components.

Definition at line 85 of file kOmegaSSTDES.C.

Referenced by kOmegaSSTDES< BasicMomentumTransportModel >::epsilonByk().

Here is the caller graph for this function:

◆ ~kOmegaSSTDES()

virtual ~kOmegaSSTDES ( )
inlinevirtual

Destructor.

Definition at line 144 of file kOmegaSSTDES.H.

References kOmegaSSTDES< BasicMomentumTransportModel >::read().

Here is the call graph for this function:

Member Function Documentation

◆ Lt()

tmp< volScalarField::Internal > Lt ( ) const
protected

Return the turbulent length-scale.

Definition at line 39 of file kOmegaSSTDES.C.

References kOmegaSSTDES< BasicMomentumTransportModel >::FDES(), and Foam::sqrt().

Here is the call graph for this function:

◆ FDES()

tmp< volScalarField::Internal > FDES ( const volScalarField::Internal F1,
const volScalarField::Internal F2 
) const
protectedvirtual

The DES dissipation-rate multiplier with options zonal filtering.

based on either F1 or F2

Definition at line 47 of file kOmegaSSTDES.C.

References delta, kOmegaSSTDES< BasicMomentumTransportModel >::epsilonByk(), Foam::exit(), F1, Foam::FatalError, FatalErrorInFunction, and Foam::max().

Referenced by kOmegaSSTDES< BasicMomentumTransportModel >::Lt().

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

◆ epsilonByk()

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

Return epsilon/k which for standard RAS is betaStar*omega.

Definition at line 72 of file kOmegaSSTDES.C.

References kOmegaSSTDES< BasicMomentumTransportModel >::kOmegaSSTDES().

Referenced by kOmegaSSTDES< BasicMomentumTransportModel >::FDES().

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

◆ TypeName()

TypeName ( "kOmegaSSTDES< BasicMomentumTransportModel >"  )

Runtime type information.

◆ read()

bool read ( )
virtual

Read model coefficients if they have changed.

Reimplemented from kOmegaSST< LESeddyViscosity< BasicMomentumTransportModel >, BasicMomentumTransportModel >.

Definition at line 131 of file kOmegaSSTDES.C.

References Foam::read().

Referenced by kOmegaSSTDES< BasicMomentumTransportModel >::~kOmegaSSTDES().

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

Member Data Documentation

◆ CDES_

dimensionedScalar CDES_
protected

DES coefficient.

Definition at line 87 of file kOmegaSSTDES.H.

◆ FSST_

direction FSST_
protected

Zonal filter choice.

  • 0: no filtering
  • 1: (1 - F1)
  • 2: (1 - F2)

Definition at line 94 of file kOmegaSSTDES.H.


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