kOmegaSST< MomentumTransportModel, BasicMomentumTransportModel > Class Template Reference

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

Inheritance diagram for kOmegaSST< MomentumTransportModel, BasicMomentumTransportModel >:
Collaboration diagram for kOmegaSST< MomentumTransportModel, BasicMomentumTransportModel >:

Public Types

typedef BasicMomentumTransportModel::alphaField alphaField
 
typedef BasicMomentumTransportModel::rhoField rhoField
 

Public Member Functions

 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...
 
virtual bool read ()
 Re-read model coefficients if they have changed. 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

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 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 MomentumTransportModel, class BasicMomentumTransportModel>
class Foam::kOmegaSST< MomentumTransportModel, BasicMomentumTransportModel >

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

Turbulence model described in:

    Menter, F. R. & Esch, T. (2001).
    Elements of Industrial Heat Transfer Prediction.
    16th Brazilian Congress of Mechanical Engineering (COBEM).

with updated coefficients from

    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.

but with the consistent production terms from the 2001 paper as form in the 2003 paper is a typo, see

    http://turbmodels.larc.nasa.gov/sst.html

and the addition of the optional F3 term for rough walls from

    Hellsten, A. (1998).
    "Some Improvements in Menter’s k-omega-SST turbulence model"
    29th AIAA Fluid Dynamics Conference, AIAA-98-2554.

Note that this implementation is written in terms of alpha diffusion coefficients rather than the more traditional sigma (alpha = 1/sigma) so that the blending can be applied to all coefficuients in a consistent manner. The paper suggests that sigma is blended but this would not be consistent with the blending of the k-epsilon and k-omega models.

Also note that the error in the last term of equation (2) relating to sigma has been corrected.

Wall-functions are applied in this implementation by using equations (14) to specify the near-wall omega as appropriate.

The blending functions (15) and (16) are not currently used because of the uncertainty in their origin, range of applicability and that if y+ becomes sufficiently small blending u_tau in this manner clearly becomes nonsense.

The default model coefficients are

    kOmegaSSTCoeffs
    {
        alphaK1     0.85;
        alphaK2     1.0;
        alphaOmega1 0.5;
        alphaOmega2 0.856;
        beta1       0.075;
        beta2       0.0828;
        betaStar    0.09;
        gamma1      5/9;
        gamma2      0.44;
        a1          0.31;
        b1          1.0;
        c1          10.0;
        F3          no;
    }
Source files

Definition at line 112 of file kOmegaSSTBase.H.

Member Typedef Documentation

◆ alphaField

typedef BasicMomentumTransportModel::alphaField alphaField

Definition at line 242 of file kOmegaSSTBase.H.

◆ rhoField

typedef BasicMomentumTransportModel::rhoField rhoField

Definition at line 243 of file kOmegaSSTBase.H.

Constructor & Destructor Documentation

◆ kOmegaSST() [1/2]

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.

Definition at line 218 of file kOmegaSSTBase.C.

Referenced by kOmegaSST< eddyViscosity< RASModel< BasicMomentumTransportModel > >, BasicMomentumTransportModel >::Qsas().

Here is the caller graph for this function:

◆ kOmegaSST() [2/2]

kOmegaSST ( const kOmegaSST< MomentumTransportModel, BasicMomentumTransportModel > &  )
delete

Disallow default bitwise copy construction.

◆ ~kOmegaSST()

virtual ~kOmegaSST ( )
inlinevirtual

Destructor.

Reimplemented in kOmegaSST< BasicMomentumTransportModel >.

Definition at line 265 of file kOmegaSSTBase.H.

Member Function Documentation

◆ F1()

virtual tmp<volScalarField> F1 ( const volScalarField CDkOmega) const
protectedvirtual

◆ F2()

virtual tmp<volScalarField> F2 ( ) const
protectedvirtual

◆ F3()

virtual tmp<volScalarField> F3 ( ) const
protectedvirtual

◆ F23()

virtual tmp<volScalarField> F23 ( ) const
protectedvirtual

◆ blend() [1/2]

◆ blend() [2/2]

tmp<volScalarField::Internal> blend ( const volScalarField::Internal F1,
const dimensionedScalar psi1,
const dimensionedScalar psi2 
) const
inlineprotected

Definition at line 172 of file kOmegaSSTBase.H.

◆ alphaK()

tmp<volScalarField> alphaK ( const volScalarField F1) const
inlineprotected

Definition at line 181 of file kOmegaSSTBase.H.

Referenced by kOmegaSST< eddyViscosity< RASModel< BasicMomentumTransportModel > >, BasicMomentumTransportModel >::DkEff().

Here is the caller graph for this function:

◆ alphaOmega()

tmp<volScalarField> alphaOmega ( const volScalarField F1) const
inlineprotected

Definition at line 186 of file kOmegaSSTBase.H.

Referenced by kOmegaSST< eddyViscosity< RASModel< BasicMomentumTransportModel > >, BasicMomentumTransportModel >::DomegaEff().

Here is the caller graph for this function:

◆ beta()

◆ gamma()

◆ correctNut() [1/2]

void correctNut ( const volScalarField S2,
const volScalarField F2 
)
protectedvirtual

Reimplemented in kOmegaSSTSato< BasicMomentumTransportModel >.

Definition at line 119 of file kOmegaSSTBase.C.

◆ correctNut() [2/2]

◆ Pk()

◆ 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 154 of file kOmegaSSTBase.C.

Referenced by kOmegaSST< eddyViscosity< RASModel< BasicMomentumTransportModel > >, BasicMomentumTransportModel >::gamma(), and kOmegaSST< eddyViscosity< RASModel< BasicMomentumTransportModel > >, BasicMomentumTransportModel >::Pk().

Here is the caller graph for this function:

◆ kSource()

tmp< fvScalarMatrix > kSource ( ) const
protectedvirtual

Definition at line 165 of file kOmegaSSTBase.C.

Referenced by kOmegaSST< eddyViscosity< RASModel< BasicMomentumTransportModel > >, BasicMomentumTransportModel >::gamma().

Here is the caller graph for this function:

◆ omegaSource()

◆ Qsas()

◆ read()

◆ DkEff()

tmp<volScalarField> DkEff ( const volScalarField F1) const
inline

Return the effective diffusivity for k.

Definition at line 275 of file kOmegaSSTBase.H.

◆ DomegaEff()

tmp<volScalarField> DomegaEff ( const volScalarField F1) const
inline

Return the effective diffusivity for omega.

Definition at line 285 of file kOmegaSSTBase.H.

◆ k()

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Return the turbulence kinetic energy.

Definition at line 295 of file kOmegaSSTBase.H.

◆ epsilon()

virtual tmp<volScalarField> epsilon ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 301 of file kOmegaSSTBase.H.

◆ omega()

virtual tmp<volScalarField> omega ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 312 of file kOmegaSSTBase.H.

◆ correct()

void correct ( )
virtual

Solve the turbulence equations and correct the turbulence viscosity.

Reimplemented in kOmegaSSTLM< BasicMomentumTransportModel >, and kOmegaSSTSato< BasicMomentumTransportModel >.

Definition at line 420 of file kOmegaSSTBase.C.

Referenced by kOmegaSST< eddyViscosity< RASModel< BasicMomentumTransportModel > >, BasicMomentumTransportModel >::omega().

Here is the caller graph for this function:

◆ operator=()

void operator= ( const kOmegaSST< MomentumTransportModel, BasicMomentumTransportModel > &  )
delete

Disallow default bitwise assignment.

Referenced by kOmegaSST< eddyViscosity< RASModel< BasicMomentumTransportModel > >, BasicMomentumTransportModel >::omega().

Here is the caller graph for this function:

Member Data Documentation

◆ alphaK1_

dimensionedScalar alphaK1_
protected

Definition at line 122 of file kOmegaSSTBase.H.

◆ alphaK2_

dimensionedScalar alphaK2_
protected

Definition at line 123 of file kOmegaSSTBase.H.

◆ alphaOmega1_

dimensionedScalar alphaOmega1_
protected

Definition at line 125 of file kOmegaSSTBase.H.

◆ alphaOmega2_

dimensionedScalar alphaOmega2_
protected

Definition at line 126 of file kOmegaSSTBase.H.

◆ gamma1_

dimensionedScalar gamma1_
protected

Definition at line 128 of file kOmegaSSTBase.H.

◆ gamma2_

dimensionedScalar gamma2_
protected

Definition at line 129 of file kOmegaSSTBase.H.

◆ beta1_

dimensionedScalar beta1_
protected

Definition at line 131 of file kOmegaSSTBase.H.

◆ beta2_

dimensionedScalar beta2_
protected

Definition at line 132 of file kOmegaSSTBase.H.

◆ betaStar_

dimensionedScalar betaStar_
protected

Definition at line 134 of file kOmegaSSTBase.H.

◆ a1_

dimensionedScalar a1_
protected

Definition at line 136 of file kOmegaSSTBase.H.

◆ b1_

dimensionedScalar b1_
protected

Definition at line 137 of file kOmegaSSTBase.H.

◆ c1_

dimensionedScalar c1_
protected

Definition at line 138 of file kOmegaSSTBase.H.

◆ F3_

Switch F3_
protected

Definition at line 140 of file kOmegaSSTBase.H.

◆ y_

const volScalarField& y_
protected

Wall distance.

Note: different to wall distance in parent RASModel which is for near-wall cells only

Definition at line 148 of file kOmegaSSTBase.H.

◆ k_

◆ omega_


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