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 >:
[legend]
Collaboration diagram for kOmegaSSTBase< BasicEddyViscosityModel >:
[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 ()=default
 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< 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

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
 Return k production rate. More...
 
virtual tmp< volScalarField::InternalepsilonByk (const volScalarField &F1, const volTensorField &gradU) const
 Return epsilon/k which for standard RAS is betaStar*omega. More...
 
virtual tmp< volScalarField::InternalGbyNu (const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
 Return G/nu. 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_
 Flag to include the F3 term. More...
 
const volScalarFieldy_
 Wall distance. More...
 
volScalarField k_
 
volScalarField omega_
 
Switch decayControl_
 Flag to include the decay control. More...
 
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 128 of file kOmegaSSTBase.H.

Member Typedef Documentation

◆ alphaField

typedef BasicEddyViscosityModel::alphaField alphaField

Definition at line 290 of file kOmegaSSTBase.H.

◆ rhoField

typedef BasicEddyViscosityModel::rhoField rhoField

Definition at line 291 of file kOmegaSSTBase.H.

◆ transportModel

typedef BasicEddyViscosityModel::transportModel transportModel

Definition at line 292 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 
)

Construct from components.

Definition at line 224 of file kOmegaSSTBase.C.

◆ ~kOmegaSSTBase()

virtual ~kOmegaSSTBase ( )
virtualdefault

Destructor.

Member Function Documentation

◆ setDecayControl()

void setDecayControl ( const dictionary dict)
protected

Definition at line 432 of file kOmegaSSTBase.C.

References dict, Foam::endl(), and Foam::Info.

Here is the call graph for this function:

◆ F1()

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

Reimplemented in kOmegaSSTLM< BasicTurbulenceModel >.

Definition at line 42 of file kOmegaSSTBase.C.

References Foam::dimless, Foam::dimTime, Foam::max(), Foam::min(), mu, Foam::pow4(), Foam::sqr(), Foam::sqrt(), and Foam::tanh().

Here is the call graph for this function:

◆ F2()

tmp< volScalarField > F2
protectedvirtual

Definition at line 72 of file kOmegaSSTBase.C.

References Foam::max(), Foam::min(), mu, Foam::sqr(), Foam::sqrt(), and Foam::tanh().

Here is the call graph for this function:

◆ F3()

tmp< volScalarField > F3
protectedvirtual

Definition at line 89 of file kOmegaSSTBase.C.

References Foam::min(), mu, Foam::pow4(), Foam::sqr(), and Foam::tanh().

Here is the call graph for this function:

◆ F23()

tmp< volScalarField > F23
protectedvirtual

Definition at line 102 of file kOmegaSSTBase.C.

References F2, F3, and tmp< T >::ref().

Here is the call graph for this function:

◆ blend() [1/2]

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

Definition at line 197 of file kOmegaSSTBase.H.

References F1, psi1, and psi2.

Referenced by kOmegaSSTBase< BasicEddyViscosityModel >::alphaK(), kOmegaSSTBase< BasicEddyViscosityModel >::alphaOmega(), kOmegaSSTBase< BasicEddyViscosityModel >::beta(), and kOmegaSSTBase< BasicEddyViscosityModel >::gamma().

Here is the caller graph for this function:

◆ blend() [2/2]

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

Definition at line 207 of file kOmegaSSTBase.H.

References F1, psi1, and psi2.

◆ alphaK()

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

Definition at line 217 of file kOmegaSSTBase.H.

References kOmegaSSTBase< BasicEddyViscosityModel >::alphaK1_, kOmegaSSTBase< BasicEddyViscosityModel >::alphaK2_, kOmegaSSTBase< BasicEddyViscosityModel >::blend(), and F1.

Referenced by kOmegaSSTBase< BasicEddyViscosityModel >::DkEff().

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

◆ alphaOmega()

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

Definition at line 222 of file kOmegaSSTBase.H.

References kOmegaSSTBase< BasicEddyViscosityModel >::alphaOmega1_, kOmegaSSTBase< BasicEddyViscosityModel >::alphaOmega2_, kOmegaSSTBase< BasicEddyViscosityModel >::blend(), and F1.

Referenced by kOmegaSSTBase< BasicEddyViscosityModel >::DomegaEff().

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

◆ beta()

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

Definition at line 227 of file kOmegaSSTBase.H.

References kOmegaSSTBase< BasicEddyViscosityModel >::beta1_, kOmegaSSTBase< BasicEddyViscosityModel >::beta2_, kOmegaSSTBase< BasicEddyViscosityModel >::blend(), Time::New(), and Foam::type().

Here is the call graph for this function:

◆ gamma()

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

Definition at line 239 of file kOmegaSSTBase.H.

References kOmegaSSTBase< BasicEddyViscosityModel >::blend(), kOmegaSSTBase< BasicEddyViscosityModel >::gamma1_, kOmegaSSTBase< BasicEddyViscosityModel >::gamma2_, Time::New(), and Foam::type().

Here is the call graph for this function:

◆ correctNut() [1/2]

void correctNut ( const volScalarField S2)
protectedvirtual

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

Definition at line 116 of file kOmegaSSTBase.C.

References Foam::max(), Time::New(), and Foam::sqrt().

Here is the call graph for this function:

◆ correctNut() [2/2]

void correctNut
protectedvirtual

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

Definition at line 129 of file kOmegaSSTBase.C.

References Foam::fvc::grad(), Foam::magSqr(), and Foam::symm().

Here is the call graph for this function:

◆ Pk()

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

Return k production rate.

Reimplemented in kOmegaSSTLM< BasicTurbulenceModel >.

Definition at line 136 of file kOmegaSSTBase.C.

References Foam::min().

Referenced by kOmegaSSTLM< BasicTurbulenceModel >::Pk().

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

◆ epsilonByk()

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

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

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

Definition at line 147 of file kOmegaSSTBase.C.

Referenced by kOmegaSSTLM< BasicTurbulenceModel >::epsilonByk().

Here is the caller graph for this function:

◆ GbyNu()

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

Return G/nu.

Reimplemented in kOmegaSSTDES< BasicTurbulenceModel >.

Definition at line 158 of file kOmegaSSTBase.C.

References F2, Foam::max(), Foam::min(), and Foam::sqrt().

Here is the call graph for this function:

◆ kSource()

tmp< fvScalarMatrix > kSource
protectedvirtual

Definition at line 175 of file kOmegaSSTBase.C.

References Foam::dimTime, and Foam::dimVolume.

◆ omegaSource()

tmp< fvScalarMatrix > omegaSource
protectedvirtual

Definition at line 189 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 203 of file kOmegaSSTBase.C.

References Foam::dimTime, and Foam::dimVolume.

◆ read()

◆ DkEff()

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

Return the effective diffusivity for k.

Definition at line 321 of file kOmegaSSTBase.H.

References kOmegaSSTBase< BasicEddyViscosityModel >::alphaK(), F1, and nu.

Here is the call graph for this function:

◆ DomegaEff()

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

Return the effective diffusivity for omega.

Definition at line 330 of file kOmegaSSTBase.H.

References kOmegaSSTBase< BasicEddyViscosityModel >::alphaOmega(), F1, and nu.

Here is the call graph for this function:

◆ k()

virtual tmp< volScalarField > k ( ) const
inlinevirtual

Return the turbulence kinetic energy.

Definition at line 343 of file kOmegaSSTBase.H.

References kOmegaSSTBase< BasicEddyViscosityModel >::k_.

◆ omega()

virtual tmp< volScalarField > omega ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 349 of file kOmegaSSTBase.H.

References kOmegaSSTBase< BasicEddyViscosityModel >::omega_.

◆ correct()

void correct
virtual

Solve the turbulence equations and correct the turbulence viscosity.

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

Definition at line 484 of file kOmegaSSTBase.C.

References Foam::fvc::absolute(), alpha, beta(), Foam::bound(), tmp< T >::clear(), optionList::constrain(), optionList::correct(), Foam::fvm::ddt(), Foam::dev(), Foam::fvc::div(), Foam::fvm::div(), divU, F1, fvOptions, gamma, Foam::fvc::grad(), Foam::fvm::laplacian(), Foam::magSqr(), Time::New(), nut, phi, tmp< T >::ref(), rho, solve(), Foam::fvm::Sp(), Foam::sqr(), Foam::fvm::SuSp(), Foam::symm(), Foam::twoSymm(), Foam::type(), and U.

Here is the call graph for this function:

Member Data Documentation

◆ alphaK1_

dimensionedScalar alphaK1_
protected

Definition at line 147 of file kOmegaSSTBase.H.

Referenced by kOmegaSSTBase< BasicEddyViscosityModel >::alphaK().

◆ alphaK2_

dimensionedScalar alphaK2_
protected

Definition at line 148 of file kOmegaSSTBase.H.

Referenced by kOmegaSSTBase< BasicEddyViscosityModel >::alphaK().

◆ alphaOmega1_

dimensionedScalar alphaOmega1_
protected

◆ alphaOmega2_

dimensionedScalar alphaOmega2_
protected

◆ gamma1_

dimensionedScalar gamma1_
protected

Definition at line 153 of file kOmegaSSTBase.H.

Referenced by kOmegaSSTBase< BasicEddyViscosityModel >::gamma().

◆ gamma2_

dimensionedScalar gamma2_
protected

Definition at line 154 of file kOmegaSSTBase.H.

Referenced by kOmegaSSTBase< BasicEddyViscosityModel >::gamma().

◆ beta1_

dimensionedScalar beta1_
protected

Definition at line 156 of file kOmegaSSTBase.H.

Referenced by kOmegaSSTBase< BasicEddyViscosityModel >::beta().

◆ beta2_

dimensionedScalar beta2_
protected

Definition at line 157 of file kOmegaSSTBase.H.

Referenced by kOmegaSSTBase< BasicEddyViscosityModel >::beta().

◆ betaStar_

dimensionedScalar betaStar_
protected

Definition at line 159 of file kOmegaSSTBase.H.

◆ a1_

dimensionedScalar a1_
protected

Definition at line 161 of file kOmegaSSTBase.H.

◆ b1_

dimensionedScalar b1_
protected

Definition at line 162 of file kOmegaSSTBase.H.

◆ c1_

dimensionedScalar c1_
protected

Definition at line 163 of file kOmegaSSTBase.H.

◆ F3_

Switch F3_
protected

Flag to include the F3 term.

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

◆ k_

volScalarField k_
protected

Definition at line 176 of file kOmegaSSTBase.H.

Referenced by kOmegaSSTBase< BasicEddyViscosityModel >::k().

◆ omega_

volScalarField omega_
protected

Definition at line 177 of file kOmegaSSTBase.H.

Referenced by kOmegaSSTBase< BasicEddyViscosityModel >::omega().

◆ decayControl_

Switch decayControl_
protected

Flag to include the decay control.

Definition at line 183 of file kOmegaSSTBase.H.

◆ kInf_

dimensionedScalar kInf_
protected

Definition at line 184 of file kOmegaSSTBase.H.

◆ omegaInf_

dimensionedScalar omegaInf_
protected

Definition at line 185 of file kOmegaSSTBase.H.


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