Renormalisation group k-epsilon turbulence model for incompressible flows. More...
#include <incompressibleRASModels/RNGkEpsilon.H>
Renormalisation group k-epsilon turbulence model for incompressible flows.
The default model coefficients correspond to the following:
RNGkEpsilonCoeffs
{
Cmu 0.0845;
C1 1.42;
C2 1.68;
sigmak 0.71942;
sigmaEps0.71942;
eta04.38;
beta0.012;
}
Definition at line 67 of file RNGkEpsilon.H.
Inheritance diagram for RNGkEpsilon:
Collaboration diagram for RNGkEpsilon:Public Member Functions | |
| TypeName ("RNGkEpsilon") | |
| Runtime type information.
| |
| RNGkEpsilon (const volVectorField &U, const surfaceScalarField &phi, transportModel &transport) | |
| Construct from components.
| |
| virtual | ~RNGkEpsilon () |
| Destructor.
| |
| virtual tmp< volScalarField > | nut () const |
| Return the turbulence viscosity.
| |
| 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 tmp< volSymmTensorField > | R () const |
| Return the Reynolds stress tensor.
| |
| virtual tmp< volSymmTensorField > | devReff () const |
| Return the effective stress tensor including the laminar stress.
| |
| virtual tmp< fvVectorMatrix > | divDevReff (volVectorField &U) const |
| Return the source term for the momentum equation.
| |
| virtual void | correct () |
| Solve the turbulence equations and correct the turbulence viscosity.
| |
| virtual bool | read () |
| Read RASProperties dictionary.
| |
| RNGkEpsilon | ( | const volVectorField & | U, |
| const surfaceScalarField & | phi, | ||
| transportModel & | transport | ||
| ) |
| virtual ~RNGkEpsilon | ( | ) | [inline, virtual]
|
Destructor.
Definition at line 108 of file RNGkEpsilon.H.
| TypeName | ( | "RNGkEpsilon" | ) |
Runtime type information.
| virtual tmp<volScalarField> nut | ( | ) | const [inline, virtual]
|
Return the turbulence viscosity.
Implements incompressible::RASModel.
Definition at line 115 of file RNGkEpsilon.H.
| tmp<volScalarField> DkEff | ( | ) | const [inline]
|
Return the effective diffusivity for k.
Definition at line 121 of file RNGkEpsilon.H.
References incompressible::turbulenceModel::nu().
Referenced by RNGkEpsilon::correct().
| tmp<volScalarField> DepsilonEff | ( | ) | const [inline]
|
Return the effective diffusivity for epsilon.
Definition at line 130 of file RNGkEpsilon.H.
References incompressible::turbulenceModel::nu().
Referenced by RNGkEpsilon::correct().
| virtual tmp<volScalarField> k | ( | ) | const [inline, virtual]
|
Return the turbulence kinetic energy.
Implements incompressible::RASModel.
Definition at line 139 of file RNGkEpsilon.H.
| virtual tmp<volScalarField> epsilon | ( | ) | const [inline, virtual]
|
Return the turbulence kinetic energy dissipation rate.
Implements incompressible::RASModel.
Definition at line 145 of file RNGkEpsilon.H.
| tmp< volSymmTensorField > R | ( | ) | const [virtual]
|
Return the Reynolds stress tensor.
Implements incompressible::RASModel.
Definition at line 167 of file RNGkEpsilon.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), Foam::fvc::grad(), Foam::I, incompressible::turbulenceModel::mesh_, IOobject::NO_READ, IOobject::NO_WRITE, incompressible::turbulenceModel::runTime_, Time::timeName(), Foam::twoSymm(), and incompressible::turbulenceModel::U_.
Referenced by RNGkEpsilon::correct().
| tmp< volSymmTensorField > devReff | ( | ) | const [virtual]
|
Return the effective stress tensor including the laminar stress.
Implements incompressible::RASModel.
Definition at line 188 of file RNGkEpsilon.C.
References Foam::dev(), Foam::fvc::grad(), incompressible::turbulenceModel::mesh_, IOobject::NO_READ, IOobject::NO_WRITE, incompressible::RASModel::nuEff(), incompressible::turbulenceModel::runTime_, Time::timeName(), Foam::twoSymm(), and incompressible::turbulenceModel::U_.
| tmp< fvVectorMatrix > divDevReff | ( | volVectorField & | U ) | const [virtual]
|
Return the source term for the momentum equation.
Implements incompressible::RASModel.
Definition at line 208 of file RNGkEpsilon.C.
References Foam::dev(), Foam::fvc::div(), Foam::fvc::grad(), Foam::fvm::laplacian(), incompressible::RASModel::nuEff(), and Foam::T().
| void correct | ( | ) | [virtual]
|
Solve the turbulence equations and correct the turbulence viscosity.
fvm::SuSp(R*G/k_, epsilon_)
Implements incompressible::RASModel.
Definition at line 239 of file RNGkEpsilon.C.
References Foam::bound(), GeometricField< Type, PatchField, GeoMesh >::boundaryField(), GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), Foam::fvm::ddt(), RNGkEpsilon::DepsilonEff(), Foam::fvm::div(), RNGkEpsilon::DkEff(), incompressible::RASModel::epsilon0_, Foam::fvc::grad(), incompressible::RASModel::k0_, Foam::fvm::laplacian(), Foam::magSqr(), incompressible::turbulenceModel::phi_, RNGkEpsilon::R(), Foam::solve(), Foam::fvm::Sp(), Foam::sqr(), Foam::sqrt(), Foam::symm(), incompressible::RASModel::turbulence_, and incompressible::turbulenceModel::U_.
| bool read | ( | ) | [virtual]
|
Read RASProperties dictionary.
Implements incompressible::RASModel.
Definition at line 218 of file RNGkEpsilon.C.
References incompressible::RASModel::coeffDict(), and dimensioned< Type >::readIfPresent().