RNGkEpsilon< BasicTurbulenceModel > Class Template Reference

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

Inheritance diagram for RNGkEpsilon< BasicTurbulenceModel >:
[legend]
Collaboration diagram for RNGkEpsilon< BasicTurbulenceModel >:
[legend]

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. More...
 
 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. More...
 
virtual ~RNGkEpsilon ()=default
 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...
 
- 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 bool read ()=0
 Re-read model coefficients if they have changed. 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< volScalarFieldk () const=0
 Return the turbulence kinetic energy. More...
 
virtual tmp< volSymmTensorFieldR () const
 Return the Reynolds stress tensor. More...
 
virtual void validate ()
 Validate the turbulence fields after construction. More...
 
virtual void correct ()=0
 Solve the turbulence equations and correct the turbulence viscosity. More...
 
- 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. More...
 
virtual ~linearViscousStress ()=default
 Destructor. More...
 
virtual bool read ()=0
 Re-read model coefficients if they have changed. 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...
 
virtual void correct ()=0
 Solve the turbulence equations and correct the turbulence viscosity. More...
 

Protected Member Functions

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

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

Member Typedef Documentation

◆ alphaField

typedef BasicTurbulenceModel::alphaField alphaField

Definition at line 132 of file RNGkEpsilon.H.

◆ rhoField

typedef BasicTurbulenceModel::rhoField rhoField

Definition at line 133 of file RNGkEpsilon.H.

◆ transportModel

typedef BasicTurbulenceModel::transportModel transportModel

Definition at line 134 of file RNGkEpsilon.H.

Constructor & Destructor Documentation

◆ RNGkEpsilon()

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 86 of file RNGkEpsilon.C.

References Foam::bound(), RNGkEpsilon< BasicTurbulenceModel >::epsilon_, RNGkEpsilon< BasicTurbulenceModel >::k_, and Foam::type().

Here is the call graph for this function:

◆ ~RNGkEpsilon()

virtual ~RNGkEpsilon ( )
virtualdefault

Destructor.

Member Function Documentation

◆ correctNut()

void correctNut
protectedvirtual

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 43 of file RNGkEpsilon.C.

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

Here is the call graph for this function:

◆ kSource()

tmp< fvScalarMatrix > kSource
protectedvirtual

Definition at line 54 of file RNGkEpsilon.C.

References Foam::dimTime, and Foam::dimVolume.

◆ epsilonSource()

tmp< fvScalarMatrix > epsilonSource
protectedvirtual

Definition at line 69 of file RNGkEpsilon.C.

References Foam::dimTime, and Foam::dimVolume.

◆ TypeName()

TypeName ( "RNGkEpsilon< BasicTurbulenceModel >"  )

Runtime type information.

◆ read()

bool read
virtual

Re-read model coefficients if they have changed.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

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

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

◆ DepsilonEff()

tmp< volScalarField > DepsilonEff ( ) const
inline

Return the effective diffusivity for epsilon.

Definition at line 180 of file RNGkEpsilon.H.

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

References RNGkEpsilon< BasicTurbulenceModel >::k_.

◆ epsilon()

virtual tmp< volScalarField > epsilon ( ) const
inlinevirtual

Return the turbulence kinetic energy dissipation rate.

Definition at line 199 of file RNGkEpsilon.H.

References RNGkEpsilon< BasicTurbulenceModel >::epsilon_.

◆ correct()

void correct
virtual

Solve the turbulence equations and correct the turbulence viscosity.

Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.

Definition at line 242 of file RNGkEpsilon.C.

References Foam::fvc::absolute(), alpha, Foam::bound(), tmp< T >::clear(), correct(), Foam::fvm::ddt(), Foam::dev(), Foam::fvc::div(), Foam::fvm::div(), divU, fvOptions, Foam::fvc::grad(), Foam::fvm::laplacian(), Foam::mag(), Time::New(), nut, phi, R, tmp< T >::ref(), rho, solve(), Foam::fvm::Sp(), Foam::sqr(), Foam::sqrt(), Foam::fvm::SuSp(), Foam::twoSymm(), Foam::type(), and U.

Here is the call graph for this function:

Member Data Documentation

◆ Cmu_

dimensionedScalar Cmu_
protected

Definition at line 107 of file RNGkEpsilon.H.

◆ C1_

dimensionedScalar C1_
protected

Definition at line 108 of file RNGkEpsilon.H.

◆ C2_

dimensionedScalar C2_
protected

Definition at line 109 of file RNGkEpsilon.H.

◆ C3_

dimensionedScalar C3_
protected

Definition at line 110 of file RNGkEpsilon.H.

◆ sigmak_

dimensionedScalar sigmak_
protected

Definition at line 111 of file RNGkEpsilon.H.

◆ sigmaEps_

dimensionedScalar sigmaEps_
protected

Definition at line 112 of file RNGkEpsilon.H.

◆ eta0_

dimensionedScalar eta0_
protected

Definition at line 113 of file RNGkEpsilon.H.

◆ beta_

dimensionedScalar beta_
protected

Definition at line 114 of file RNGkEpsilon.H.

◆ k_

◆ epsilon_


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