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

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

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

Public Types

typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 
- Public Types inherited from eddyViscosity< RASModel< BasicTurbulenceModel > >
typedef RASModel< BasicTurbulenceModel >::alphaField alphaField
 
typedef RASModel< BasicTurbulenceModel >::rhoField rhoField
 
typedef RASModel< BasicTurbulenceModel >::transportModel transportModel
 
- Public Types inherited from linearViscousStress< RASModel< BasicTurbulenceModel > >
typedef RASModel< BasicTurbulenceModel >::alphaField alphaField
 
typedef RASModel< BasicTurbulenceModel >::rhoField rhoField
 
typedef RASModel< BasicTurbulenceModel >::transportModel transportModel
 
- Public Types inherited from RASModel< BasicTurbulenceModel >
typedef BasicTurbulenceModel::alphaField alphaField
 
typedef BasicTurbulenceModel::rhoField rhoField
 
typedef BasicTurbulenceModel::transportModel transportModel
 

Public Member Functions

 TypeName ("kOmegaSST")
 Runtime type information. More...
 
 kOmegaSST (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
 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 eddyViscosity< RASModel< BasicTurbulenceModel > >
 eddyViscosity (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
 Construct from components. More...
 
virtual ~eddyViscosity ()
 Destructor. More...
 
volScalarFieldevNut ()
 Return non-const access to the turbulence viscosity. More...
 
virtual tmp< volScalarFieldnut () const
 Return the turbulence viscosity. More...
 
virtual tmp< scalarFieldnut (const label patchi) const
 Return the turbulence viscosity on patch. More...
 
virtual tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor. More...
 
- Public Member Functions inherited from linearViscousStress< RASModel< BasicTurbulenceModel > >
 linearViscousStress (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
 Construct from components. More...
 
virtual ~linearViscousStress ()
 Destructor. More...
 
virtual tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor. More...
 
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 Return the source term for the momentum equation. More...
 
virtual tmp< fvVectorMatrixdivDevRhoReff (const volScalarField &rho, volVectorField &U) const
 Return the source term for the momentum equation. More...
 
- Public Member Functions inherited from RASModel< BasicTurbulenceModel >
 TypeName ("RAS")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, RASModel, 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))
 
 RASModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
 Construct from components. More...
 
virtual ~RASModel ()
 Destructor. More...
 
const dimensionedScalarkMin () const
 Return the lower allowable limit for k (default: SMALL) More...
 
const dimensionedScalarepsilonMin () const
 Return the lower allowable limit for epsilon (default: SMALL) More...
 
const dimensionedScalaromegaMin () const
 Return the lower allowable limit for omega (default: SMALL) More...
 
dimensionedScalarkMin ()
 Allow kMin to be changed. More...
 
dimensionedScalarepsilonMin ()
 Allow epsilonMin to be changed. More...
 
dimensionedScalaromegaMin ()
 Allow omegaMin to be changed. More...
 
virtual const dictionarycoeffDict () const
 Const access to the coefficients dictionary. More...
 
virtual tmp< volScalarFieldnuEff () const
 Return the effective viscosity. More...
 
virtual tmp< scalarFieldnuEff (const label patchi) const
 Return the effective viscosity on patch. 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)
 
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 Member Functions inherited from RASModel< BasicTurbulenceModel >
virtual void printCoeffs (const word &type)
 Print model coefficients. More...
 

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 eddyViscosity< RASModel< BasicTurbulenceModel > >
volScalarField nut_
 
- Protected Attributes inherited from RASModel< BasicTurbulenceModel >
dictionary RASDict_
 RAS coefficients dictionary. More...
 
Switch turbulence_
 Turbulence on/off flag. More...
 
Switch printCoeffs_
 Flag to print the model coeffs at run-time. More...
 
dictionary coeffDict_
 Model coefficients dictionary. More...
 
dimensionedScalar kMin_
 Lower limit of k. More...
 
dimensionedScalar epsilonMin_
 Lower limit of epsilon. More...
 
dimensionedScalar omegaMin_
 Lower limit for omega. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from RASModel< BasicTurbulenceModel >
static autoPtr< RASModelNew (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 RAS model. More...
 

Detailed Description

template<class BasicTurbulenceModel>
class Foam::RASModels::kOmegaSST< 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 120 of file kOmegaSST.H.

Member Typedef Documentation

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 224 of file kOmegaSST.H.

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 225 of file kOmegaSST.H.

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 226 of file kOmegaSST.H.

Constructor & Destructor Documentation

kOmegaSST ( const alphaField alpha,
const rhoField rho,
const volVectorField U,
const surfaceScalarField alphaRhoPhi,
const surfaceScalarField phi,
const transportModel transport,
const word propertiesName = turbulenceModel::propertiesName,
const word type = typeName 
)

Construct from components.

Definition at line 180 of file kOmegaSST.C.

References Foam::bound().

Here is the call graph for this function:

virtual ~kOmegaSST ( )
inlinevirtual

Destructor.

Definition at line 250 of file kOmegaSST.H.

References kOmegaSST< BasicTurbulenceModel >::read().

Here is the call graph for this function:

Member Function Documentation

tmp<volScalarField> F1 ( const volScalarField CDkOmega) const
protected
tmp<volScalarField> F2 ( ) const
protected
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 186 of file kOmegaSST.H.

References kOmegaSST< BasicTurbulenceModel >::blend().

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

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 191 of file kOmegaSST.H.

References kOmegaSST< BasicTurbulenceModel >::blend().

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

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 196 of file kOmegaSST.H.

References kOmegaSST< BasicTurbulenceModel >::blend().

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

Here is the call graph for this function:

Here is the caller graph for this function:

tmp<volScalarField> gamma ( const volScalarField F1) const
inlineprotected
void correctNut ( const volScalarField S2)
protected

Definition at line 111 of file kOmegaSST.C.

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

Here is the call graph for this function:

void correctNut ( )
protectedvirtual

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Reimplemented in kOmegaSSTSato< BasicTurbulenceModel >.

Definition at line 123 of file kOmegaSST.C.

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

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

Here is the call graph for this function:

Here is the caller graph for this function:

tmp< fvScalarMatrix > kSource ( ) const
protectedvirtual

Definition at line 130 of file kOmegaSST.C.

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

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

Here is the caller graph for this function:

tmp< fvScalarMatrix > omegaSource ( ) const
protectedvirtual

Definition at line 144 of file kOmegaSST.C.

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

Referenced by kOmegaSST< 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

Reimplemented in kOmegaSSTSAS< BasicTurbulenceModel >.

Definition at line 159 of file kOmegaSST.C.

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

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

Here is the caller graph for this function:

TypeName ( "kOmegaSST< BasicTurbulenceModel >"  )

Runtime type information.

bool read ( )
virtual

Re-read model coefficients if they have changed.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

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

Definition at line 362 of file kOmegaSST.C.

References Foam::read().

Referenced by kOmegaSST< 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 260 of file kOmegaSST.H.

References kOmegaSST< BasicTurbulenceModel >::alphaK(), nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.

Here is the call graph for this function:

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

Return the effective diffusivity for omega.

Definition at line 269 of file kOmegaSST.H.

References kOmegaSST< BasicTurbulenceModel >::alphaOmega(), nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.

Here is the call graph for this function:

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Return the turbulence kinetic energy.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 282 of file kOmegaSST.H.

References kOmegaSST< BasicTurbulenceModel >::k_.

virtual tmp<volScalarField> epsilon ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 288 of file kOmegaSST.H.

References GeometricField< Type, PatchField, GeoMesh >::boundaryField().

Here is the call graph for this function:

virtual tmp<volScalarField> omega ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 307 of file kOmegaSST.H.

References kOmegaSST< BasicTurbulenceModel >::correct(), and kOmegaSST< BasicTurbulenceModel >::omega_.

Here is the call graph for this function:

void correct ( )
virtual

Member Data Documentation

dimensionedScalar alphaK1_
protected

Definition at line 137 of file kOmegaSST.H.

dimensionedScalar alphaK2_
protected

Definition at line 138 of file kOmegaSST.H.

dimensionedScalar alphaOmega1_
protected

Definition at line 140 of file kOmegaSST.H.

dimensionedScalar alphaOmega2_
protected

Definition at line 141 of file kOmegaSST.H.

dimensionedScalar gamma1_
protected

Definition at line 143 of file kOmegaSST.H.

dimensionedScalar gamma2_
protected

Definition at line 144 of file kOmegaSST.H.

dimensionedScalar beta1_
protected

Definition at line 146 of file kOmegaSST.H.

dimensionedScalar beta2_
protected

Definition at line 147 of file kOmegaSST.H.

dimensionedScalar betaStar_
protected

Definition at line 149 of file kOmegaSST.H.

dimensionedScalar a1_
protected

Definition at line 151 of file kOmegaSST.H.

dimensionedScalar b1_
protected

Definition at line 152 of file kOmegaSST.H.

dimensionedScalar c1_
protected

Definition at line 153 of file kOmegaSST.H.

Switch F3_
protected

Definition at line 155 of file kOmegaSST.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 163 of file kOmegaSST.H.

volScalarField k_
protected

Definition at line 165 of file kOmegaSST.H.

Referenced by kOmegaSST< BasicTurbulenceModel >::k().

volScalarField omega_
protected

Definition at line 166 of file kOmegaSST.H.

Referenced by kOmegaSST< BasicTurbulenceModel >::omega().


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