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

Renormalization group k-epsilon turbulence model for incompressible and compressible flows. More...

#include <RNGkEpsilon.H>

Inheritance diagram for RNGkEpsilon< BasicTurbulenceModel >:
Collaboration diagram for RNGkEpsilon< 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 ("RNGkEpsilon")
 Runtime type information.
 RNGkEpsilon (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 ~RNGkEpsilon ()=default
 Destructor.
virtual bool read ()
 Re-read model coefficients if they have changed.
tmp< volScalarFieldDkEff () const
 Return the effective diffusivity for k.
tmp< volScalarFieldDepsilonEff () const
 Return the effective diffusivity for epsilon.
virtual tmp< volScalarFieldk () const
 Return the turbulence kinetic energy.
virtual tmp< volScalarFieldepsilon () const
 Return the turbulence kinetic energy 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 ()
virtual tmp< fvScalarMatrixkSource () const
virtual tmp< fvScalarMatrixepsilonSource () const

Protected Attributes

dimensionedScalar Cmu_
dimensionedScalar C1_
dimensionedScalar C2_
dimensionedScalar C3_
dimensionedScalar sigmak_
dimensionedScalar sigmaEps_
dimensionedScalar eta0_
dimensionedScalar beta_
volScalarField k_
volScalarField epsilon_
Protected Attributes inherited from eddyViscosity< RASModel< BasicTurbulenceModel > >
volScalarField nut_

Detailed Description

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

Renormalization group k-epsilon turbulence model for incompressible and compressible flows.

Reference:

    Yakhot, V., Orszag, S. A., Thangam, S.,
    Gatski, T. B., & Speziale, C. G. (1992).
    Development of turbulence models for shear flows
    by a double expansion technique.
    Physics of Fluids A: Fluid Dynamics (1989-1993), 4(7), 1510-1520.

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

    RNGkEpsilonCoeffs
    {
        Cmu         0.0845;
        C1          1.42;
        C2          1.68;
        C3          -0.33;
        sigmak      0.71942;
        sigmaEps    0.71942;
        eta0        4.38;
        beta        0.012;
    }
Source files

Definition at line 84 of file RNGkEpsilon.H.

Member Typedef Documentation

◆ alphaField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 132 of file RNGkEpsilon.H.

◆ rhoField

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 133 of file RNGkEpsilon.H.

◆ transportModel

template<class BasicTurbulenceModel>
typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 134 of file RNGkEpsilon.H.

Constructor & Destructor Documentation

◆ RNGkEpsilon()

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

References alpha, beta_, Foam::bound(), C1_, C2_, C3_, Cmu_, eddyViscosity< RASModel< BasicTurbulenceModel > >::eddyViscosity(), epsilon_, eta0_, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, k_, phi, rho, sigmaEps_, sigmak_, timeName, and U.

Here is the call graph for this function:

◆ ~RNGkEpsilon()

template<class BasicTurbulenceModel>
virtual ~RNGkEpsilon ( )
virtualdefault

Destructor.

Member Function Documentation

◆ correctNut()

template<class BasicTurbulenceModel>
void correctNut ( )
protectedvirtual

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 36 of file RNGkEpsilon.C.

References Cmu_, optionList::correct(), epsilon_, k_, options::New(), eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_, and Foam::sqr().

Referenced by correct().

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

◆ kSource()

template<class BasicTurbulenceModel>
tmp< fvScalarMatrix > kSource ( ) const
protectedvirtual

Definition at line 47 of file RNGkEpsilon.C.

References Foam::dimTime, Foam::dimVolume, k_, and tmp< T >::New().

Referenced by correct().

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

◆ epsilonSource()

template<class BasicTurbulenceModel>
tmp< fvScalarMatrix > epsilonSource ( ) const
protectedvirtual

Definition at line 58 of file RNGkEpsilon.C.

References Foam::dimTime, Foam::dimVolume, epsilon_, and tmp< T >::New().

Referenced by correct().

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

◆ TypeName()

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

Runtime type information.

◆ read()

template<class BasicTurbulenceModel>
bool read ( )
virtual

Re-read model coefficients if they have changed.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 206 of file RNGkEpsilon.C.

References beta_, C1_, C2_, C3_, Cmu_, eddyViscosity< RASModel< BasicTurbulenceModel > >::eddyViscosity(), eta0_, read(), sigmaEps_, and sigmak_.

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 177 of file RNGkEpsilon.H.

Referenced by correct().

Here is the caller graph for this function:

◆ DepsilonEff()

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

Return the effective diffusivity for epsilon.

Definition at line 192 of file RNGkEpsilon.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 207 of file RNGkEpsilon.H.

◆ epsilon()

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

Return the turbulence kinetic energy dissipation rate.

Definition at line 215 of file RNGkEpsilon.H.

◆ correct()

Member Data Documentation

◆ Cmu_

template<class BasicTurbulenceModel>
dimensionedScalar Cmu_
protected

Definition at line 107 of file RNGkEpsilon.H.

Referenced by correct(), correctNut(), read(), and RNGkEpsilon().

◆ C1_

template<class BasicTurbulenceModel>
dimensionedScalar C1_
protected

Definition at line 108 of file RNGkEpsilon.H.

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

◆ C2_

template<class BasicTurbulenceModel>
dimensionedScalar C2_
protected

Definition at line 109 of file RNGkEpsilon.H.

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

◆ C3_

template<class BasicTurbulenceModel>
dimensionedScalar C3_
protected

Definition at line 110 of file RNGkEpsilon.H.

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

◆ sigmak_

template<class BasicTurbulenceModel>
dimensionedScalar sigmak_
protected

Definition at line 111 of file RNGkEpsilon.H.

Referenced by read(), and RNGkEpsilon().

◆ sigmaEps_

template<class BasicTurbulenceModel>
dimensionedScalar sigmaEps_
protected

Definition at line 112 of file RNGkEpsilon.H.

Referenced by read(), and RNGkEpsilon().

◆ eta0_

template<class BasicTurbulenceModel>
dimensionedScalar eta0_
protected

Definition at line 113 of file RNGkEpsilon.H.

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

◆ beta_

template<class BasicTurbulenceModel>
dimensionedScalar beta_
protected

Definition at line 114 of file RNGkEpsilon.H.

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

◆ k_

template<class BasicTurbulenceModel>
volScalarField k_
protected

Definition at line 119 of file RNGkEpsilon.H.

Referenced by correct(), correctNut(), kSource(), and RNGkEpsilon().

◆ epsilon_

template<class BasicTurbulenceModel>
volScalarField epsilon_
protected

Definition at line 120 of file RNGkEpsilon.H.

Referenced by correct(), correctNut(), epsilonSource(), and RNGkEpsilon().


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