Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions
kOmegaSSTBase< BasicEddyViscosityModel > Class Template Reference

Base class implementation of the k-omega-SST turbulence model for incompressible and compressible flows. More...

Inheritance diagram for kOmegaSSTBase< BasicEddyViscosityModel >:
Inheritance graph
[legend]
Collaboration diagram for kOmegaSSTBase< BasicEddyViscosityModel >:
Collaboration graph
[legend]

Public Types

typedef BasicEddyViscosityModel::alphaField alphaField
 
typedef BasicEddyViscosityModel::rhoField rhoField
 
typedef BasicEddyViscosityModel::transportModel transportModel
 

Public Member Functions

 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...
 
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...
 

Protected Member Functions

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 void correctNut (const volScalarField &S2)
 
virtual void correctNut ()
 
virtual tmp< fvScalarMatrixkSource () const
 
virtual tmp< fvScalarMatrixomegaSource () const
 
virtual tmp< fvScalarMatrixQsas (const volScalarField &S2, const volScalarField &gamma, const volScalarField &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_
 

Private Member Functions

 kOmegaSSTBase (const kOmegaSSTBase &)
 
kOmegaSSTBaseoperator= (const kOmegaSSTBase &)
 

Detailed Description

template<class BasicEddyViscosityModel>
class Foam::kOmegaSSTBase< BasicEddyViscosityModel >

Base class 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

    kOmegaSSTBaseCoeffs
    {
        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 115 of file kOmegaSSTBase.H.

Member Typedef Documentation

◆ alphaField

typedef BasicEddyViscosityModel::alphaField alphaField

Definition at line 216 of file kOmegaSSTBase.H.

◆ rhoField

typedef BasicEddyViscosityModel::rhoField rhoField

Definition at line 217 of file kOmegaSSTBase.H.

◆ transportModel

typedef BasicEddyViscosityModel::transportModel transportModel

Definition at line 218 of file kOmegaSSTBase.H.

Constructor & Destructor Documentation

◆ kOmegaSSTBase() [1/2]

kOmegaSSTBase ( const kOmegaSSTBase< BasicEddyViscosityModel > &  )
private

◆ kOmegaSSTBase() [2/2]

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.

Definition at line 181 of file kOmegaSSTBase.C.

◆ ~kOmegaSSTBase()

virtual ~kOmegaSSTBase ( )
inlinevirtual

Destructor.

Definition at line 238 of file kOmegaSSTBase.H.

Member Function Documentation

◆ operator=()

kOmegaSSTBase& operator= ( const kOmegaSSTBase< BasicEddyViscosityModel > &  )
private

◆ F1()

tmp< volScalarField > F1 ( const volScalarField CDkOmega) const
protected

◆ F2()

tmp< volScalarField > F2
protected

Definition at line 68 of file kOmegaSSTBase.C.

◆ F3()

tmp< volScalarField > F3
protected

Definition at line 85 of file kOmegaSSTBase.C.

◆ F23()

tmp< volScalarField > F23
protected

Definition at line 98 of file kOmegaSSTBase.C.

◆ blend()

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

◆ alphaK()

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

Definition at line 181 of file kOmegaSSTBase.H.

Referenced by kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::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 kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::DomegaEff().

Here is the caller graph for this function:

◆ beta()

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

Definition at line 191 of file kOmegaSSTBase.H.

◆ gamma()

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

Definition at line 196 of file kOmegaSSTBase.H.

◆ correctNut() [1/2]

void correctNut ( const volScalarField S2)
protectedvirtual

◆ correctNut() [2/2]

void correctNut
protectedvirtual

◆ kSource()

tmp< fvScalarMatrix > kSource
protectedvirtual

Definition at line 131 of file kOmegaSSTBase.C.

◆ omegaSource()

tmp< fvScalarMatrix > omegaSource
protectedvirtual

Definition at line 145 of file kOmegaSSTBase.C.

◆ Qsas()

tmp< fvScalarMatrix > Qsas ( const volScalarField S2,
const volScalarField gamma,
const volScalarField beta 
) const
protectedvirtual

Reimplemented in kOmegaSSTSAS< BasicTurbulenceModel >.

Definition at line 160 of file kOmegaSSTBase.C.

◆ read()

bool read
virtual

◆ DkEff()

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

Return the effective diffusivity for k.

Definition at line 248 of file kOmegaSSTBase.H.

◆ DomegaEff()

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

Return the effective diffusivity for omega.

Definition at line 257 of file kOmegaSSTBase.H.

◆ k()

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Return the turbulence kinetic energy.

Definition at line 270 of file kOmegaSSTBase.H.

◆ epsilon()

virtual tmp<volScalarField> epsilon ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 276 of file kOmegaSSTBase.H.

◆ omega()

virtual tmp<volScalarField> omega ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 295 of file kOmegaSSTBase.H.

◆ correct()

void correct
virtual

Solve the turbulence equations and correct the turbulence viscosity.

Reimplemented in kOmegaSSTSato< BasicTurbulenceModel >, and kOmegaSSTDES< BasicTurbulenceModel >.

Definition at line 385 of file kOmegaSSTBase.C.

Field Documentation

◆ alphaK1_

dimensionedScalar alphaK1_
protected

◆ alphaK2_

dimensionedScalar alphaK2_
protected

◆ alphaOmega1_

dimensionedScalar alphaOmega1_
protected

◆ alphaOmega2_

dimensionedScalar alphaOmega2_
protected

◆ gamma1_

dimensionedScalar gamma1_
protected

◆ gamma2_

dimensionedScalar gamma2_
protected

◆ beta1_

dimensionedScalar beta1_
protected

◆ beta2_

dimensionedScalar beta2_
protected

◆ betaStar_

dimensionedScalar betaStar_
protected

◆ a1_

dimensionedScalar a1_
protected

Definition at line 146 of file kOmegaSSTBase.H.

◆ b1_

dimensionedScalar b1_
protected

Definition at line 147 of file kOmegaSSTBase.H.

◆ c1_

dimensionedScalar c1_
protected

Definition at line 148 of file kOmegaSSTBase.H.

◆ F3_

Switch F3_
protected

Definition at line 150 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 158 of file kOmegaSSTBase.H.

◆ k_

volScalarField k_
protected

◆ omega_

volScalarField omega_
protected

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