Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
kOmegaSST< TurbulenceModel, BasicTurbulenceModel > Class Template Reference

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

Inheritance diagram for kOmegaSST< TurbulenceModel, BasicTurbulenceModel >:
Inheritance graph
[legend]
Collaboration diagram for kOmegaSST< TurbulenceModel, BasicTurbulenceModel >:
Collaboration graph
[legend]

Public Types

typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 
- Public Types inherited from TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel >
typedef Alpha alphaField
 
typedef Rho rhoField
 
typedef TransportModel transportModel
 

Public Member Functions

 kOmegaSST (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 ~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...
 
- Public Member Functions inherited from TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel >
 declareRunTimeNewSelectionTable (autoPtr, TurbulenceModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName),(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName))
 
 TurbulenceModel (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
 Construct. More...
 
virtual ~TurbulenceModel ()
 Destructor. More...
 
const alphaFieldalpha () const
 Access function to phase fraction. More...
 
const transportModeltransport () const
 Access function to incompressible transport model. More...
 
tmp< volScalarFieldnu () const
 Return the laminar viscosity. More...
 
tmp< scalarFieldnu (const label patchi) const
 Return the laminar viscosity on patchi. 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
 
void correctNut (const volScalarField &S2, const volScalarField &F2)
 
virtual void correctNut ()
 
virtual tmp< volScalarFieldepsilonByk (const volScalarField &F1, const volScalarField &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 &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_
 
- Protected Attributes inherited from TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel >
const alphaFieldalpha_
 
const transportModeltransport_
 

Additional Inherited Members

- Static Public Member Functions inherited from TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel >
static autoPtr< TurbulenceModelNew (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName)
 Return a reference to the selected turbulence model. More...
 

Detailed Description

template<class TurbulenceModel, class BasicTurbulenceModel>
class Foam::kOmegaSST< TurbulenceModel, BasicTurbulenceModel >

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

Member Typedef Documentation

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 229 of file kOmegaSSTBase.H.

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 230 of file kOmegaSSTBase.H.

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 231 of file kOmegaSSTBase.H.

Constructor & Destructor Documentation

kOmegaSST ( 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 201 of file kOmegaSSTBase.C.

References Foam::bound().

Here is the call graph for this function:

virtual ~kOmegaSST ( )
inlinevirtual

Destructor.

Reimplemented in kOmegaSST< BasicTurbulenceModel >.

Definition at line 251 of file kOmegaSSTBase.H.

Member Function Documentation

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

Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().

Here is the caller graph for this function:

tmp<volScalarField> F3 ( ) const
protected
tmp<volScalarField> F23 ( ) const
protected
tmp<volScalarField> blend ( const volScalarField F1,
const dimensionedScalar psi1,
const dimensionedScalar psi2 
) const
inlineprotected
tmp<volScalarField> alphaK ( const volScalarField F1) const
inlineprotected

Definition at line 181 of file kOmegaSSTBase.H.

Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::DkEff().

Here is the caller graph for this function:

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

Definition at line 186 of file kOmegaSSTBase.H.

Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::DomegaEff().

Here is the caller graph for this function:

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

Definition at line 191 of file kOmegaSSTBase.H.

Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().

Here is the caller graph for this function:

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

Definition at line 196 of file kOmegaSSTBase.H.

Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().

Here is the caller graph for this function:

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

Definition at line 115 of file kOmegaSSTBase.C.

References optionList::correct(), Foam::max(), options::New(), and Foam::sqrt().

Here is the call graph for this function:

void correctNut ( )
protectedvirtual

Reimplemented in kOmegaSSTSato< BasicTurbulenceModel >.

Definition at line 131 of file kOmegaSSTBase.C.

References kOmegaSST< TurbulenceModel, BasicTurbulenceModel >::epsilonByk(), Foam::fvc::grad(), Foam::magSqr(), and Foam::symm().

Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().

Here is the call graph for this function:

Here is the caller graph for this function:

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

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

Reimplemented in kOmegaSSTDES< BasicTurbulenceModel >.

Definition at line 139 of file kOmegaSSTBase.C.

Referenced by kOmegaSST< TurbulenceModel, BasicTurbulenceModel >::correctNut(), and kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().

Here is the caller graph for this function:

tmp< fvScalarMatrix > kSource ( ) const
protectedvirtual

Definition at line 150 of file kOmegaSSTBase.C.

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

Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().

Here is the caller graph for this function:

tmp< fvScalarMatrix > omegaSource ( ) const
protectedvirtual

Definition at line 165 of file kOmegaSSTBase.C.

References Foam::dimTime, Foam::dimVolume, and kOmegaSST< TurbulenceModel, BasicTurbulenceModel >::Qsas().

Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().

Here is the call graph for this function:

Here is the caller graph for this function:

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

Re-read model coefficients if they have changed.

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

Definition at line 377 of file kOmegaSSTBase.C.

References Foam::read().

Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::~kOmegaSST().

Here is the call graph for this function:

Here is the caller graph for this function:

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

Return the effective diffusivity for k.

Definition at line 261 of file kOmegaSSTBase.H.

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

Return the effective diffusivity for omega.

Definition at line 270 of file kOmegaSSTBase.H.

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Return the turbulence kinetic energy.

Definition at line 283 of file kOmegaSSTBase.H.

virtual tmp<volScalarField> epsilon ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 289 of file kOmegaSSTBase.H.

virtual tmp<volScalarField> omega ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 308 of file kOmegaSSTBase.H.

void correct ( )
virtual

Member Data Documentation

dimensionedScalar alphaK1_
protected

Definition at line 132 of file kOmegaSSTBase.H.

dimensionedScalar alphaK2_
protected

Definition at line 133 of file kOmegaSSTBase.H.

dimensionedScalar alphaOmega1_
protected

Definition at line 135 of file kOmegaSSTBase.H.

dimensionedScalar alphaOmega2_
protected

Definition at line 136 of file kOmegaSSTBase.H.

dimensionedScalar gamma1_
protected

Definition at line 138 of file kOmegaSSTBase.H.

dimensionedScalar gamma2_
protected

Definition at line 139 of file kOmegaSSTBase.H.

dimensionedScalar beta1_
protected

Definition at line 141 of file kOmegaSSTBase.H.

dimensionedScalar beta2_
protected

Definition at line 142 of file kOmegaSSTBase.H.

dimensionedScalar betaStar_
protected

Definition at line 144 of file kOmegaSSTBase.H.

dimensionedScalar a1_
protected

Definition at line 146 of file kOmegaSSTBase.H.

dimensionedScalar b1_
protected

Definition at line 147 of file kOmegaSSTBase.H.

dimensionedScalar c1_
protected

Definition at line 148 of file kOmegaSSTBase.H.

Switch F3_
protected

Definition at line 150 of file kOmegaSSTBase.H.

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.

volScalarField k_
protected
volScalarField omega_
protected

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