Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
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)
 
virtual ~kOmegaSSTBase ()=default
 
virtual bool read ()
 
tmp< volScalarFieldDkEff (const volScalarField &F1) const
 
tmp< volScalarFieldDomegaEff (const volScalarField &F1) const
 
virtual tmp< volScalarFieldk () const
 
virtual tmp< volScalarFieldomega () const
 
virtual void correct ()
 

Protected Member Functions

void setDecayControl (const dictionary &dict)
 
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)
 
virtual void correctNut ()
 
virtual tmp< volScalarField::InternalPk (const volScalarField::Internal &G) const
 
virtual tmp< volScalarField::InternalepsilonByk (const volScalarField &F1, const volTensorField &gradU) const
 
virtual tmp< volScalarField::InternalGbyNu (const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
 
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_
 
volScalarField k_
 
volScalarField omega_
 
Switch decayControl_
 
dimensionedScalar kInf_
 
dimensionedScalar omegaInf_
 

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.

and the optional decay control from:

    Spalart, P. R. and Rumsey, C. L. (2007).
    Effective Inflow Conditions for Turbulence Models in Aerodynamic
    Calculations
    AIAA Journal, 45(10), 2544 - 2553.

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 coefficients 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;

        // Optional decay control
        decayControl    yes;
        kInf            \<far-field k value\>;
        omegaInf        \<far-field omega value\>;
    }
Source files

Definition at line 124 of file kOmegaSSTBase.H.

Member Typedef Documentation

◆ alphaField

typedef BasicEddyViscosityModel::alphaField alphaField

Definition at line 286 of file kOmegaSSTBase.H.

◆ rhoField

typedef BasicEddyViscosityModel::rhoField rhoField

Definition at line 287 of file kOmegaSSTBase.H.

◆ transportModel

typedef BasicEddyViscosityModel::transportModel transportModel

Definition at line 288 of file kOmegaSSTBase.H.

Constructor & Destructor Documentation

◆ kOmegaSSTBase()

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 
)

Definition at line 218 of file kOmegaSSTBase.C.

◆ ~kOmegaSSTBase()

virtual ~kOmegaSSTBase ( )
virtualdefault

Member Function Documentation

◆ setDecayControl()

void setDecayControl ( const dictionary dict)
protected

Definition at line 426 of file kOmegaSSTBase.C.

◆ F1()

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

◆ F2()

tmp< volScalarField > F2
protectedvirtual

Definition at line 65 of file kOmegaSSTBase.C.

◆ F3()

tmp< volScalarField > F3
protectedvirtual

Definition at line 82 of file kOmegaSSTBase.C.

◆ F23()

tmp< volScalarField > F23
protectedvirtual

Definition at line 95 of file kOmegaSSTBase.C.

◆ blend() [1/2]

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

◆ blend() [2/2]

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

Definition at line 204 of file kOmegaSSTBase.H.

◆ alphaK()

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

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

Referenced by kOmegaSSTBase< eddyViscosity< RASModel< BasicTurbulenceModel > > >::DomegaEff().

Here is the caller graph for this function:

◆ beta()

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

Definition at line 224 of file kOmegaSSTBase.H.

◆ gamma()

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

Definition at line 236 of file kOmegaSSTBase.H.

◆ correctNut() [1/2]

void correctNut ( const volScalarField S2)
protectedvirtual

◆ correctNut() [2/2]

void correctNut
protectedvirtual

◆ Pk()

tmp< volScalarField::Internal > Pk ( const volScalarField::Internal G) const
protectedvirtual

Reimplemented in kOmegaSSTLM< BasicTurbulenceModel >.

Definition at line 130 of file kOmegaSSTBase.C.

◆ epsilonByk()

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

◆ GbyNu()

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

Reimplemented in kOmegaSSTDES< BasicTurbulenceModel >.

Definition at line 152 of file kOmegaSSTBase.C.

◆ kSource()

tmp< fvScalarMatrix > kSource
protectedvirtual

Definition at line 168 of file kOmegaSSTBase.C.

◆ omegaSource()

tmp< fvScalarMatrix > omegaSource
protectedvirtual

Definition at line 182 of file kOmegaSSTBase.C.

◆ Qsas()

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

Reimplemented in kOmegaSSTSAS< BasicTurbulenceModel >.

Definition at line 197 of file kOmegaSSTBase.C.

◆ read()

bool read
virtual

◆ DkEff()

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

Definition at line 317 of file kOmegaSSTBase.H.

◆ DomegaEff()

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

Definition at line 326 of file kOmegaSSTBase.H.

◆ k()

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Definition at line 339 of file kOmegaSSTBase.H.

◆ omega()

virtual tmp<volScalarField> omega ( ) const
inlinevirtual

Definition at line 345 of file kOmegaSSTBase.H.

◆ correct()

void correct
virtual

Member Data 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

Definition at line 155 of file kOmegaSSTBase.H.

◆ a1_

dimensionedScalar a1_
protected

Definition at line 157 of file kOmegaSSTBase.H.

◆ b1_

dimensionedScalar b1_
protected

Definition at line 158 of file kOmegaSSTBase.H.

◆ c1_

dimensionedScalar c1_
protected

Definition at line 159 of file kOmegaSSTBase.H.

◆ F3_

Switch F3_
protected

Definition at line 162 of file kOmegaSSTBase.H.

◆ y_

const volScalarField& y_
protected

Definition at line 170 of file kOmegaSSTBase.H.

◆ k_

volScalarField k_
protected

◆ omega_

volScalarField omega_
protected

◆ decayControl_

Switch decayControl_
protected

Definition at line 179 of file kOmegaSSTBase.H.

◆ kInf_

dimensionedScalar kInf_
protected

Definition at line 180 of file kOmegaSSTBase.H.

◆ omegaInf_

dimensionedScalar omegaInf_
protected

Definition at line 181 of file kOmegaSSTBase.H.


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