All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
kEpsilon< BasicMomentumTransportModel > Class Template Reference

Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distortion theory (RDT) based compression term. More...

Inheritance diagram for kEpsilon< BasicMomentumTransportModel >:
Collaboration diagram for kEpsilon< BasicMomentumTransportModel >:

Public Types

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

Public Member Functions

 TypeName ("kEpsilon")
 Runtime type information. More...
 
 kEpsilon (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
 Construct from components. More...
 
 kEpsilon (const kEpsilon &)=delete
 Disallow default bitwise copy construction. More...
 
virtual ~kEpsilon ()
 Destructor. More...
 
virtual bool read ()
 Re-read model coefficients if they have changed. More...
 
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k. More...
 
tmp< volScalarFieldDepsilonEff () const
 Return the effective diffusivity for epsilon. More...
 
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy. More...
 
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy dissipation rate. More...
 
virtual void correct ()
 Solve the turbulence equations and correct the turbulence viscosity. More...
 
void operator= (const kEpsilon &)=delete
 Disallow default bitwise assignment. More...
 
- Public Member Functions inherited from eddyViscosity< RASModel< BasicMomentumTransportModel > >
 eddyViscosity (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Construct from components. More...
 
virtual ~eddyViscosity ()
 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< volSymmTensorFieldsigma () const
 Return the Reynolds stress tensor [m^2/s^2]. More...
 
virtual void validate ()
 Validate the turbulence fields after construction. More...
 
- Public Member Functions inherited from linearViscousStress< RASModel< BasicMomentumTransportModel > >
 linearViscousStress (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Construct from components. More...
 
virtual ~linearViscousStress ()
 Destructor. More...
 
virtual tmp< volSymmTensorFielddevTau () const
 Return the effective stress tensor. More...
 
virtual tmp< fvVectorMatrixdivDevTau (volVectorField &U) const
 Return the source term for the momentum equation. More...
 
virtual tmp< fvVectorMatrixdivDevTau (const volScalarField &rho, volVectorField &U) const
 Return the source term for the momentum equation. More...
 
- Public Member Functions inherited from RASModel< BasicMomentumTransportModel >
 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),(alpha, rho, U, alphaRhoPhi, phi, transport))
 
 RASModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Construct from components. More...
 
 RASModel (const RASModel &)=delete
 Disallow default bitwise copy construction. 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...
 
void operator= (const RASModel &)=delete
 Disallow default bitwise assignment. More...
 

Protected Member Functions

virtual void correctNut ()
 
virtual tmp< fvScalarMatrixkSource () const
 
virtual tmp< fvScalarMatrixepsilonSource () const
 
- Protected Member Functions inherited from RASModel< BasicMomentumTransportModel >
virtual void printCoeffs (const word &type)
 Print model coefficients. More...
 

Protected Attributes

dimensionedScalar Cmu_
 
dimensionedScalar C1_
 
dimensionedScalar C2_
 
dimensionedScalar C3_
 
dimensionedScalar sigmak_
 
dimensionedScalar sigmaEps_
 
volScalarField k_
 
volScalarField epsilon_
 
- Protected Attributes inherited from eddyViscosity< RASModel< BasicMomentumTransportModel > >
volScalarField nut_
 
- Protected Attributes inherited from RASModel< BasicMomentumTransportModel >
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< BasicMomentumTransportModel >
static autoPtr< RASModelNew (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
 Return a reference to the selected RAS model. More...
 

Detailed Description

template<class BasicMomentumTransportModel>
class Foam::RASModels::kEpsilon< BasicMomentumTransportModel >

Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distortion theory (RDT) based compression term.

Reference:

    Standard model:
        Launder, B. E., & Spalding, D. B. (1972).
        Lectures in mathematical models of turbulence.

        Launder, B. E., & Spalding, D. B. (1974).
        The numerical computation of turbulent flows.
        Computer methods in applied mechanics and engineering,
        3(2), 269-289.

    For the RDT-based compression term:
        El Tahry, S. H. (1983).
        k-epsilon equation for compressible reciprocating engine flows.
        Journal of Energy, 7(4), 345-353.

The default model coefficients are

    kEpsilonCoeffs
    {
        Cmu         0.09;
        C1          1.44;
        C2          1.92;
        C3          -0.33;
        sigmak      1.0;
        sigmaEps    1.3;
    }
Source files

Definition at line 83 of file kEpsilon.H.

Member Typedef Documentation

◆ alphaField

typedef BasicMomentumTransportModel::alphaField alphaField

Definition at line 115 of file kEpsilon.H.

◆ rhoField

typedef BasicMomentumTransportModel::rhoField rhoField

Definition at line 116 of file kEpsilon.H.

◆ transportModel

typedef BasicMomentumTransportModel::transportModel transportModel

Definition at line 117 of file kEpsilon.H.

Constructor & Destructor Documentation

◆ kEpsilon() [1/2]

kEpsilon ( const alphaField alpha,
const rhoField rho,
const volVectorField U,
const surfaceScalarField alphaRhoPhi,
const surfaceScalarField phi,
const transportModel transport,
const word type = typeName 
)

Construct from components.

Definition at line 82 of file kEpsilon.C.

Referenced by kEpsilon< compressible::momentumTransportModel >::epsilonSource().

Here is the caller graph for this function:

◆ kEpsilon() [2/2]

kEpsilon ( const kEpsilon< BasicMomentumTransportModel > &  )
delete

Disallow default bitwise copy construction.

◆ ~kEpsilon()

virtual ~kEpsilon ( )
inlinevirtual

Destructor.

Definition at line 143 of file kEpsilon.H.

Member Function Documentation

◆ correctNut()

void correctNut ( )
protectedvirtual

◆ kSource()

tmp< fvScalarMatrix > kSource ( ) const
protectedvirtual

◆ epsilonSource()

tmp< fvScalarMatrix > epsilonSource ( ) const
protectedvirtual

◆ TypeName()

TypeName ( "kEpsilon< BasicMomentumTransportModel >"  )

Runtime type information.

◆ read()

bool read ( )
virtual

◆ DkEff()

tmp<volScalarField> DkEff ( ) const
inline

Return the effective diffusivity for k.

Definition at line 153 of file kEpsilon.H.

◆ DepsilonEff()

tmp<volScalarField> DepsilonEff ( ) const
inline

Return the effective diffusivity for epsilon.

Definition at line 163 of file kEpsilon.H.

◆ k()

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Return the turbulence kinetic energy.

Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.

Definition at line 173 of file kEpsilon.H.

◆ epsilon()

virtual tmp<volScalarField> epsilon ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 179 of file kEpsilon.H.

◆ correct()

void correct ( )
virtual

Solve the turbulence equations and correct the turbulence viscosity.

Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.

Reimplemented in LaheyKEpsilon< BasicMomentumTransportModel >, and PDRkEpsilon.

Definition at line 217 of file kEpsilon.C.

Referenced by LaheyKEpsilon< BasicMomentumTransportModel >::correct(), and kEpsilon< compressible::momentumTransportModel >::epsilon().

Here is the caller graph for this function:

◆ operator=()

void operator= ( const kEpsilon< BasicMomentumTransportModel > &  )
delete

Disallow default bitwise assignment.

Referenced by kEpsilon< compressible::momentumTransportModel >::epsilon().

Here is the caller graph for this function:

Member Data Documentation

◆ Cmu_

dimensionedScalar Cmu_
protected

Definition at line 93 of file kEpsilon.H.

◆ C1_

dimensionedScalar C1_
protected

Definition at line 94 of file kEpsilon.H.

◆ C2_

dimensionedScalar C2_
protected

Definition at line 95 of file kEpsilon.H.

◆ C3_

dimensionedScalar C3_
protected

Definition at line 96 of file kEpsilon.H.

◆ sigmak_

dimensionedScalar sigmak_
protected

Definition at line 97 of file kEpsilon.H.

◆ sigmaEps_

dimensionedScalar sigmaEps_
protected

Definition at line 98 of file kEpsilon.H.

◆ k_

volScalarField k_
protected

Definition at line 102 of file kEpsilon.H.

Referenced by kEpsilon< compressible::momentumTransportModel >::k().

◆ epsilon_

volScalarField epsilon_
protected

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