Renormalization group k-epsilon turbulence model for incompressible and compressible flows. More...
#include <RNGkEpsilon.H>


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< volScalarField > | DkEff () const |
| Return the effective diffusivity for k. | |
| tmp< volScalarField > | DepsilonEff () const |
| Return the effective diffusivity for epsilon. | |
| virtual tmp< volScalarField > | k () const |
| Return the turbulence kinetic energy. | |
| virtual tmp< volScalarField > | epsilon () 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< volScalarField > | nut () const |
| Return the turbulence viscosity. | |
| virtual tmp< volSymmTensorField > | R () 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< volSymmTensorField > | devRhoReff () const |
| Return the effective stress tensor. | |
| virtual tmp< volSymmTensorField > | devRhoReff (const volVectorField &U) const |
| Return the effective stress tensor based on a given velocity field. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
| Return the source term for the momentum equation. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (const volScalarField &rho, volVectorField &U) const |
| Return the source term for the momentum equation. | |
Protected Member Functions | |
| virtual void | correctNut () |
| virtual tmp< fvScalarMatrix > | kSource () const |
| virtual tmp< fvScalarMatrix > | epsilonSource () 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_ |
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;
}
Definition at line 84 of file RNGkEpsilon.H.
| typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 132 of file RNGkEpsilon.H.
| typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 133 of file RNGkEpsilon.H.
| typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 134 of file RNGkEpsilon.H.
| 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.

|
virtualdefault |
Destructor.
|
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().


|
protectedvirtual |
Definition at line 47 of file RNGkEpsilon.C.
References Foam::dimTime, Foam::dimVolume, k_, and tmp< T >::New().
Referenced by correct().


|
protectedvirtual |
Definition at line 58 of file RNGkEpsilon.C.
References Foam::dimTime, Foam::dimVolume, epsilon_, and tmp< T >::New().
Referenced by correct().


| TypeName | ( | "RNGkEpsilon< BasicTurbulenceModel >" | ) |
Runtime type information.
|
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().


|
inline |
Return the effective diffusivity for k.
Definition at line 177 of file RNGkEpsilon.H.
Referenced by correct().

|
inline |
Return the effective diffusivity for epsilon.
Definition at line 192 of file RNGkEpsilon.H.
Referenced by correct().

|
inlinevirtual |
Return the turbulence kinetic energy.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 207 of file RNGkEpsilon.H.
|
inlinevirtual |
Return the turbulence kinetic energy dissipation rate.
Definition at line 215 of file RNGkEpsilon.H.
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 227 of file RNGkEpsilon.C.
References Foam::fvc::absolute(), alpha, beta_, Foam::bound(), C1_, C2_, C3_, tmp< T >::clear(), Cmu_, correct(), correctNut(), Foam::fvm::ddt(), DepsilonEff(), Foam::devTwoSymm(), Foam::fvc::div(), Foam::fvm::div(), divU, DkEff(), eddyViscosity< RASModel< BasicTurbulenceModel > >::eddyViscosity(), epsilon_, epsilonSource(), eta0_, fvOptions, Foam::fvc::grad(), k_, kSource(), Foam::fvm::laplacian(), Foam::mag(), options::New(), eddyViscosity< RASModel< BasicTurbulenceModel > >::nut(), nut, eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_, phi, eddyViscosity< RASModel< BasicTurbulenceModel > >::R(), R, tmp< T >::ref(), rho, IOobject::scopedName(), solve(), Foam::fvm::Sp(), Foam::sqr(), Foam::sqrt(), Foam::fvm::SuSp(), Foam::type(), and U.

|
protected |
Definition at line 107 of file RNGkEpsilon.H.
Referenced by correct(), correctNut(), read(), and RNGkEpsilon().
|
protected |
Definition at line 108 of file RNGkEpsilon.H.
Referenced by correct(), read(), and RNGkEpsilon().
|
protected |
Definition at line 109 of file RNGkEpsilon.H.
Referenced by correct(), read(), and RNGkEpsilon().
|
protected |
Definition at line 110 of file RNGkEpsilon.H.
Referenced by correct(), read(), and RNGkEpsilon().
|
protected |
Definition at line 111 of file RNGkEpsilon.H.
Referenced by read(), and RNGkEpsilon().
|
protected |
Definition at line 112 of file RNGkEpsilon.H.
Referenced by read(), and RNGkEpsilon().
|
protected |
Definition at line 113 of file RNGkEpsilon.H.
Referenced by correct(), read(), and RNGkEpsilon().
|
protected |
Definition at line 114 of file RNGkEpsilon.H.
Referenced by correct(), read(), and RNGkEpsilon().
|
protected |
Definition at line 119 of file RNGkEpsilon.H.
Referenced by correct(), correctNut(), kSource(), and RNGkEpsilon().
|
protected |
Definition at line 120 of file RNGkEpsilon.H.
Referenced by correct(), correctNut(), epsilonSource(), and RNGkEpsilon().