kOmega< BasicTurbulenceModel > Class Template Reference

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

Inheritance diagram for kOmega< BasicTurbulenceModel >:
[legend]
Collaboration diagram for kOmega< BasicTurbulenceModel >:
[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")
 Runtime type information. More...
 
 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. More...
 
virtual ~kOmega ()=default
 Destructor. More...
 
virtual bool read ()
 Read RASProperties dictionary. More...
 
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k. More...
 
tmp< volScalarFieldDomegaEff () 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 specific 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 ()=default
 Destructor. 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...
 
virtual void validate ()
 Validate the turbulence fields after construction. 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 ()=default
 Destructor. More...
 
virtual tmp< volSymmTensorFielddevRhoReff () const
 Return the effective stress tensor. More...
 
virtual tmp< volSymmTensorFielddevRhoReff (const volVectorField &U) const
 Return the effective stress tensor based on a given velocity field. 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 ()=default
 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...
 
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate. More...
 

Protected Member Functions

virtual void correctNut ()
 
- Protected Member Functions inherited from RASModel< BasicTurbulenceModel >
virtual void printCoeffs (const word &type)
 Print model coefficients. More...
 
 RASModel (const RASModel &)=delete
 No copy construct. More...
 
void operator= (const RASModel &)=delete
 No copy assignment. More...
 

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_
 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::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 78 of file kOmega.H.

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 109 of file kOmega.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 110 of file kOmega.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

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

Construct from components.

Definition at line 57 of file kOmega.C.

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

Here is the call graph for this function:

◆ ~kOmega()

virtual ~kOmega ( )
virtualdefault

Destructor.

Member Function Documentation

◆ correctNut()

void correctNut ( )
protectedvirtual

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 43 of file kOmega.C.

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

Here is the call graph for this function:

◆ TypeName()

TypeName ( "kOmega< BasicTurbulenceModel >"  )

Runtime type information.

◆ read()

bool read ( )
virtual

Read RASProperties dictionary.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 164 of file kOmega.C.

References Foam::read().

Here is the call graph for this function:

◆ DkEff()

tmp<volScalarField> DkEff ( ) const
inline

Return the effective diffusivity for k.

Definition at line 144 of file kOmega.H.

References kOmega< BasicTurbulenceModel >::alphaK_, nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.

◆ DomegaEff()

tmp<volScalarField> DomegaEff ( ) const
inline

Return the effective diffusivity for omega.

Definition at line 157 of file kOmega.H.

References kOmega< BasicTurbulenceModel >::alphaOmega_, nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.

◆ k()

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Return the turbulence kinetic energy.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 170 of file kOmega.H.

References kOmega< BasicTurbulenceModel >::k_.

◆ omega()

virtual tmp<volScalarField> omega ( ) const
inlinevirtual

Return the turbulence specific dissipation rate.

Reimplemented from RASModel< BasicTurbulenceModel >.

Definition at line 176 of file kOmega.H.

References kOmega< BasicTurbulenceModel >::omega_.

◆ correct()

void correct ( )
virtual

Solve the turbulence equations and correct the turbulence viscosity.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 182 of file kOmega.C.

References Foam::fvc::absolute(), Foam::constant::atomic::alpha, Foam::bound(), tmp< T >::clear(), correct(), Foam::fvm::ddt(), Foam::dev(), Foam::fvm::div(), Foam::fvc::div(), divU, fvOptions, Foam::constant::universal::G, Foam::fvc::grad(), Foam::fvm::laplacian(), options::New(), nut, phi, tmp< T >::ref(), rho, Foam::solve(), Foam::fvm::Sp(), Foam::fvm::SuSp(), Foam::twoSymm(), and U.

Here is the call graph for this function:

Member Data Documentation

◆ Cmu_

dimensionedScalar Cmu_
protected

Definition at line 89 of file kOmega.H.

◆ beta_

dimensionedScalar beta_
protected

Definition at line 90 of file kOmega.H.

◆ gamma_

dimensionedScalar gamma_
protected

Definition at line 91 of file kOmega.H.

◆ alphaK_

dimensionedScalar alphaK_
protected

Definition at line 92 of file kOmega.H.

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

◆ alphaOmega_

dimensionedScalar alphaOmega_
protected

Definition at line 93 of file kOmega.H.

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

◆ k_

volScalarField k_
protected

Definition at line 98 of file kOmega.H.

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

◆ omega_

volScalarField omega_
protected

Definition at line 99 of file kOmega.H.

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


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