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

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

Inheritance diagram for kOmega< BasicTurbulenceModel >:
Inheritance graph
[legend]
Collaboration diagram for kOmega< 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 ("kOmega")
 
 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)
 
virtual ~kOmega ()=default
 
virtual bool read ()
 
tmp< volScalarFieldDkEff () const
 
tmp< volScalarFieldDomegaEff () const
 
virtual tmp< volScalarFieldk () const
 
virtual tmp< volScalarFieldomega () const
 
virtual void correct ()
 
- 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)
 
virtual ~eddyViscosity ()=default
 
virtual tmp< volScalarFieldnut () const
 
virtual tmp< scalarFieldnut (const label patchi) const
 
virtual tmp< volScalarFieldk () const=0
 
virtual tmp< volSymmTensorFieldR () const
 
virtual void validate ()
 
- 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)
 
virtual ~linearViscousStress ()=default
 
virtual tmp< volSymmTensorFielddevRhoReff () const
 
virtual tmp< volSymmTensorFielddevRhoReff (const volVectorField &U) const
 
virtual tmp< fvVectorMatrixdivDevRhoReff (volVectorField &U) const
 
virtual tmp< fvVectorMatrixdivDevRhoReff (const volScalarField &rho, volVectorField &U) const
 
- Public Member Functions inherited from RASModel< BasicTurbulenceModel >
 TypeName ("RAS")
 
 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)
 
virtual ~RASModel ()=default
 
const dimensionedScalarkMin () const
 
const dimensionedScalarepsilonMin () const
 
const dimensionedScalaromegaMin () const
 
dimensionedScalarkMin ()
 
dimensionedScalarepsilonMin ()
 
dimensionedScalaromegaMin ()
 
virtual const dictionarycoeffDict () const
 
virtual tmp< volScalarFieldnuEff () const
 
virtual tmp< scalarFieldnuEff (const label patchi) const
 
virtual tmp< volScalarFieldepsilon () const
 

Protected Member Functions

virtual void correctNut ()
 
- Protected Member Functions inherited from RASModel< BasicTurbulenceModel >
virtual void printCoeffs (const word &type)
 
 RASModel (const RASModel &)=delete
 
void operator= (const RASModel &)=delete
 

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_
 
- Protected Attributes inherited from RASModel< BasicTurbulenceModel >
dictionary RASDict_
 
Switch turbulence_
 
Switch printCoeffs_
 
dictionary coeffDict_
 
dimensionedScalar kMin_
 
dimensionedScalar epsilonMin_
 
dimensionedScalar omegaMin_
 

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)
 

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. (1998).
    Turbulence modeling for CFD
    (Vol. 2, pp. 103-217). La Canada, CA: DCW industries.

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 74 of file kOmega.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 105 of file kOmega.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 106 of file kOmega.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 107 of file kOmega.H.

Constructor & Destructor Documentation

◆ kOmega()

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 
)

Definition at line 50 of file kOmega.C.

References Foam::bound(), and Foam::type().

Here is the call graph for this function:

◆ ~kOmega()

virtual ~kOmega ( )
virtualdefault

Member Function Documentation

◆ correctNut()

void correctNut
protectedvirtual

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 36 of file kOmega.C.

References optionList::correct(), and options::New().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "kOmega< BasicTurbulenceModel >"  )

◆ read()

bool read
virtual

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 157 of file kOmega.C.

References Foam::read().

Here is the call graph for this function:

◆ DkEff()

tmp<volScalarField> DkEff ( ) const
inline

◆ DomegaEff()

tmp<volScalarField> DomegaEff ( ) const
inline

◆ k()

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Definition at line 166 of file kOmega.H.

References kOmega< BasicTurbulenceModel >::k_.

◆ omega()

virtual tmp<volScalarField> omega ( ) const
inlinevirtual

Reimplemented from RASModel< BasicTurbulenceModel >.

Definition at line 172 of file kOmega.H.

References kOmega< BasicTurbulenceModel >::omega_.

◆ correct()

void correct
virtual

Member Data Documentation

◆ Cmu_

dimensionedScalar Cmu_
protected

Definition at line 85 of file kOmega.H.

◆ beta_

dimensionedScalar beta_
protected

Definition at line 86 of file kOmega.H.

◆ gamma_

dimensionedScalar gamma_
protected

Definition at line 87 of file kOmega.H.

◆ alphaK_

dimensionedScalar alphaK_
protected

Definition at line 88 of file kOmega.H.

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

◆ alphaOmega_

dimensionedScalar alphaOmega_
protected

Definition at line 89 of file kOmega.H.

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

◆ k_

volScalarField k_
protected

Definition at line 94 of file kOmega.H.

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

◆ omega_

volScalarField omega_
protected

Definition at line 95 of file kOmega.H.

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


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