Loading...
Searching...
No Matches
kOmega< BasicTurbulenceModel > Class Template Reference

Standard high Reynolds-number k-omega turbulence model for incompressible and compressible flows. More...

#include <kOmega.H>

Inheritance diagram for kOmega< BasicTurbulenceModel >:
Collaboration diagram for kOmega< BasicTurbulenceModel >:

Public Types

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

Public Member Functions

 TypeName ("kOmega")
 Runtime type information.
 kOmega (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.
virtual ~kOmega ()=default
 Destructor.
virtual bool read ()
 Read RASProperties dictionary.
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k.
tmp< volScalarFieldDomegaEff () const
 Return the effective diffusivity for omega.
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy.
virtual tmp< volScalarFieldomega () const
 Return the turbulence specific dissipation rate.
virtual void correct ()
 Solve the turbulence equations and correct the turbulence viscosity.
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.
virtual ~eddyViscosity ()=default
 Destructor.
virtual tmp< volScalarFieldnut () const
 Return the turbulence viscosity.
virtual tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor.
virtual void validate ()
 Validate the turbulence fields after construction.
Public Member Functions inherited from linearViscousStress< 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.
virtual ~linearViscousStress ()=default
 Destructor.
virtual tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor.
virtual tmp< volSymmTensorFielddevRhoReff (const volVectorField &U) const
 Return the effective stress tensor based on a given velocity field.
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 Return the source term for the momentum equation.
virtual tmp< fvVectorMatrixdivDevRhoReff (const volScalarField &rho, volVectorField &U) const
 Return the source term for the momentum equation.

Protected Member Functions

virtual void correctNut ()

Protected Attributes

dimensionedScalar Cmu_
dimensionedScalar beta_
dimensionedScalar gamma_
dimensionedScalar alphaK_
dimensionedScalar alphaOmega_
volScalarField k_
volScalarField omega_
Protected Attributes inherited from eddyViscosity< RASModel< BasicTurbulenceModel > >
volScalarField nut_

Detailed Description

template<class BasicTurbulenceModel>
class Foam::RASModels::kOmega< BasicTurbulenceModel >

Standard high Reynolds-number k-omega turbulence model for incompressible and compressible flows.

References:

    Wilcox, D. C. (1988).
    Reassessment of the scale-determining
    equation for advanced turbulence models.
    AIAA Journal, 26(11), 1299-1310.
    DOI:10.2514/3.10041

The default model coefficients are

    kOmegaCoeffs
    {
        Cmu         0.09;  // Equivalent to betaStar
        alpha       0.52;
        beta        0.072;
        alphak      0.5;
        alphaOmega  0.5;
    }
Source files

Definition at line 76 of file kOmega.H.

Member Typedef Documentation

◆ alphaField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 107 of file kOmega.H.

◆ rhoField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 108 of file kOmega.H.

◆ transportModel

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 109 of file kOmega.H.

Constructor & Destructor Documentation

◆ kOmega()

template<class BasicTurbulenceModel>
kOmega ( 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 49 of file kOmega.C.

References alpha, alphaK_, alphaOmega_, beta_, Foam::bound(), Cmu_, eddyViscosity< RASModel< BasicTurbulenceModel > >::eddyViscosity(), gamma_, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, k_, omega_, phi, rho, timeName, and U.

Here is the call graph for this function:

◆ ~kOmega()

template<class BasicTurbulenceModel>
virtual ~kOmega ( )
virtualdefault

Destructor.

Member Function Documentation

◆ correctNut()

template<class BasicTurbulenceModel>
void correctNut ( )
protectedvirtual

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 36 of file kOmega.C.

References optionList::correct(), k_, options::New(), eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_, and omega_.

Referenced by correct().

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

◆ TypeName()

template<class BasicTurbulenceModel>
TypeName ( "kOmega< BasicTurbulenceModel >" )

Runtime type information.

◆ read()

template<class BasicTurbulenceModel>
bool read ( )
virtual

Read RASProperties dictionary.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 157 of file kOmega.C.

References alphaK_, alphaOmega_, beta_, Cmu_, eddyViscosity< RASModel< BasicTurbulenceModel > >::eddyViscosity(), gamma_, and read().

Referenced by read().

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

◆ DkEff()

template<class BasicTurbulenceModel>
tmp< volScalarField > DkEff ( ) const
inline

Return the effective diffusivity for k.

Definition at line 152 of file kOmega.H.

Referenced by correct().

Here is the caller graph for this function:

◆ DomegaEff()

template<class BasicTurbulenceModel>
tmp< volScalarField > DomegaEff ( ) const
inline

Return the effective diffusivity for omega.

Definition at line 167 of file kOmega.H.

Referenced by correct().

Here is the caller graph for this function:

◆ k()

template<class BasicTurbulenceModel>
virtual tmp< volScalarField > k ( ) const
inlinevirtual

Return the turbulence kinetic energy.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 182 of file kOmega.H.

◆ omega()

template<class BasicTurbulenceModel>
virtual tmp< volScalarField > omega ( ) const
inlinevirtual

Return the turbulence specific dissipation rate.

Definition at line 190 of file kOmega.H.

◆ correct()

Member Data Documentation

◆ Cmu_

template<class BasicTurbulenceModel>
dimensionedScalar Cmu_
protected

Definition at line 87 of file kOmega.H.

Referenced by correct(), kOmega(), and read().

◆ beta_

template<class BasicTurbulenceModel>
dimensionedScalar beta_
protected

Definition at line 88 of file kOmega.H.

Referenced by correct(), kOmega(), and read().

◆ gamma_

template<class BasicTurbulenceModel>
dimensionedScalar gamma_
protected

Definition at line 89 of file kOmega.H.

Referenced by correct(), kOmega(), and read().

◆ alphaK_

template<class BasicTurbulenceModel>
dimensionedScalar alphaK_
protected

Definition at line 90 of file kOmega.H.

Referenced by kOmega(), and read().

◆ alphaOmega_

template<class BasicTurbulenceModel>
dimensionedScalar alphaOmega_
protected

Definition at line 91 of file kOmega.H.

Referenced by kOmega(), and read().

◆ k_

template<class BasicTurbulenceModel>
volScalarField k_
protected

Definition at line 96 of file kOmega.H.

Referenced by correct(), correctNut(), and kOmega().

◆ omega_

template<class BasicTurbulenceModel>
volScalarField omega_
protected

Definition at line 97 of file kOmega.H.

Referenced by correct(), correctNut(), and kOmega().


The documentation for this class was generated from the following files:
  • src/TurbulenceModels/turbulenceModels/RAS/kOmega/kOmega.H
  • src/TurbulenceModels/turbulenceModels/RAS/kOmega/kOmega.C