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
 
- Public Types inherited from eddyViscosity< RASModel< BasicMomentumTransportModel > >
typedef BasicMomentumTransportModel::alphaField alphaField
 
typedef BasicMomentumTransportModel::rhoField rhoField
 
- Public Types inherited from linearViscousStress< BasicMomentumTransportModel >
typedef BasicMomentumTransportModel::alphaField alphaField
 
typedef BasicMomentumTransportModel::rhoField rhoField
 

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 viscosity &viscosity, 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 tmp< volScalarFieldomega () const
 Return the turbulence specific 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 viscosity &viscosity)
 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< BasicMomentumTransportModel >
 linearViscousStress (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
 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...
 

Protected Member Functions

virtual void correctNut ()
 
virtual tmp< fvScalarMatrixkSource () const
 
virtual tmp< fvScalarMatrixepsilonSource () const
 

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_
 

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;
        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.

Constructor & Destructor Documentation

◆ kEpsilon() [1/2]

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

Construct from components.

Definition at line 82 of file kEpsilon.C.

◆ kEpsilon() [2/2]

kEpsilon ( const kEpsilon< BasicMomentumTransportModel > &  )
delete

Disallow default bitwise copy construction.

◆ ~kEpsilon()

virtual ~kEpsilon ( )
inlinevirtual

Destructor.

Definition at line 142 of file kEpsilon.H.

Member Function Documentation

◆ correctNut()

void correctNut
protectedvirtual

Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.

Reimplemented in LaheyKEpsilon< BasicMomentumTransportModel >, and continuousGasKEpsilon< BasicMomentumTransportModel >.

Definition at line 41 of file kEpsilon.C.

References dictionary::New(), and Foam::sqr().

Referenced by continuousGasKEpsilon< BasicMomentumTransportModel >::correctNut().

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

◆ kSource()

◆ epsilonSource()

◆ TypeName()

TypeName ( "kEpsilon< BasicMomentumTransportModel >"  )

Runtime type information.

◆ read()

bool read
virtual

Re-read model coefficients if they have changed.

Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.

Reimplemented in LaheyKEpsilon< BasicMomentumTransportModel >, continuousGasKEpsilon< BasicMomentumTransportModel >, buoyantKEpsilon< BasicMomentumTransportModel >, and PDRkEpsilon.

Definition at line 197 of file kEpsilon.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 152 of file kEpsilon.H.

References GeometricField< Type, PatchField, GeoMesh >::New(), and eddyViscosity< RASModel< BasicMomentumTransportModel > >::nut_.

Here is the call graph for this function:

◆ DepsilonEff()

tmp<volScalarField> DepsilonEff ( ) const
inline

Return the effective diffusivity for epsilon.

Definition at line 162 of file kEpsilon.H.

References GeometricField< Type, PatchField, GeoMesh >::New(), and eddyViscosity< RASModel< BasicMomentumTransportModel > >::nut_.

Here is the call graph for this function:

◆ k()

virtual tmp<volScalarField> k ( ) const
inlinevirtual

Return the turbulence kinetic energy.

Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.

Definition at line 172 of file kEpsilon.H.

References kEpsilon< BasicMomentumTransportModel >::k_.

◆ epsilon()

virtual tmp<volScalarField> epsilon ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 178 of file kEpsilon.H.

References kEpsilon< BasicMomentumTransportModel >::epsilon_.

◆ omega()

virtual tmp<volScalarField> omega ( ) const
inlinevirtual

Return the turbulence specific dissipation rate.

Definition at line 184 of file kEpsilon.H.

References kEpsilon< BasicMomentumTransportModel >::Cmu_, kEpsilon< BasicMomentumTransportModel >::epsilon_, kEpsilon< BasicMomentumTransportModel >::k_, and GeometricField< Type, PatchField, GeoMesh >::New().

Here is the call graph for this function:

◆ correct()

◆ operator=()

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

Disallow default bitwise assignment.

Member Data Documentation

◆ Cmu_

dimensionedScalar Cmu_
protected

Definition at line 93 of file kEpsilon.H.

Referenced by kEpsilon< BasicMomentumTransportModel >::omega().

◆ 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_

◆ epsilon_


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