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 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 ("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 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< volScalarFieldomega () const
 Return the specific dissipation rate. More...
 

Protected Member Functions

virtual void correctNut ()
 
virtual tmp< fvScalarMatrixkSource () const
 
virtual tmp< fvScalarMatrixepsilonSource () const
 
- 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 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_
 
- 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::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 87 of file RNGkEpsilon.C.

References Foam::bound(), 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 optionList::correct(), options::New(), and Foam::sqr().

Here is the call graph for this function:

◆ kSource()

tmp< fvScalarMatrix > kSource ( ) const
protectedvirtual

Definition at line 54 of file RNGkEpsilon.C.

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

◆ epsilonSource()

tmp< fvScalarMatrix > epsilonSource ( ) const
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.

Reimplemented from RASModel< BasicTurbulenceModel >.

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(), 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(), Foam::mag(), options::New(), nut, phi, R, tmp< T >::ref(), rho, Foam::solve(), Foam::fvm::Sp(), Foam::sqr(), Foam::sqrt(), 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 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_

volScalarField k_
protected

Definition at line 119 of file RNGkEpsilon.H.

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

◆ epsilon_

volScalarField epsilon_
protected

Definition at line 120 of file RNGkEpsilon.H.

Referenced by RNGkEpsilon< BasicTurbulenceModel >::epsilon().


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