Renormalisation group k-epsilon turbulence model for incompressible and compressible flows. More...
Public Types | |
typedef BasicMomentumTransportModel::alphaField | alphaField |
typedef BasicMomentumTransportModel::rhoField | rhoField |
Public Types inherited from eddyViscosity< RASModel< BasicMomentumTransportModel > > | |
typedef RASModel< BasicMomentumTransportModel > ::alphaField | alphaField |
typedef RASModel< BasicMomentumTransportModel > ::rhoField | rhoField |
Public Types inherited from linearViscousStress< RASModel< BasicMomentumTransportModel > > | |
typedef RASModel< BasicMomentumTransportModel > ::alphaField | alphaField |
typedef RASModel< BasicMomentumTransportModel > ::rhoField | rhoField |
Public Types inherited from RASModel< BasicMomentumTransportModel > | |
typedef BasicMomentumTransportModel::alphaField | alphaField |
typedef BasicMomentumTransportModel::rhoField | rhoField |
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 viscosity &viscosity, const word &type=typeName) | |
Construct from components. More... | |
RNGkEpsilon (const RNGkEpsilon &)=delete | |
Disallow default bitwise copy construction. More... | |
virtual | ~RNGkEpsilon () |
Destructor. More... | |
virtual bool | read () |
Re-read model coefficients if they have changed. More... | |
tmp< volScalarField > | DkEff () const |
Return the effective diffusivity for k. More... | |
tmp< volScalarField > | DepsilonEff () const |
Return the effective diffusivity for epsilon. More... | |
virtual tmp< volScalarField > | k () const |
Return the turbulence kinetic energy. More... | |
virtual tmp< volScalarField > | epsilon () const |
Return the turbulence kinetic energy dissipation rate. More... | |
virtual tmp< volScalarField > | omega () const |
Return the turbulence specific dissipation rate. More... | |
virtual void | correct () |
Solve the turbulence equations and correct the turbulence viscosity. More... | |
void | operator= (const RNGkEpsilon &)=delete |
Disallow default bitwise assignment. More... | |
Public Member Functions inherited from eddyViscosity< RASModel< BasicMomentumTransportModel > > | |
eddyViscosity (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity) | |
Construct from components. More... | |
virtual | ~eddyViscosity () |
Destructor. More... | |
virtual tmp< volScalarField > | nut () const |
Return the turbulence viscosity. More... | |
virtual tmp< scalarField > | nut (const label patchi) const |
Return the turbulence viscosity on patch. More... | |
virtual tmp< volSymmTensorField > | sigma () const |
Return the Reynolds stress tensor [m^2/s^2]. More... | |
virtual void | validate () |
Validate the turbulence fields after construction. More... | |
Public Member Functions inherited from linearViscousStress< RASModel< BasicMomentumTransportModel > > | |
linearViscousStress (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity) | |
Construct from components. More... | |
virtual | ~linearViscousStress () |
Destructor. More... | |
virtual tmp< volSymmTensorField > | devTau () const |
Return the effective stress tensor. More... | |
virtual tmp< fvVectorMatrix > | divDevTau (volVectorField &U) const |
Return the source term for the momentum equation. More... | |
virtual tmp< fvVectorMatrix > | divDevTau (const volScalarField &rho, volVectorField &U) const |
Return the source term for the momentum equation. More... | |
Public Member Functions inherited from RASModel< BasicMomentumTransportModel > | |
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 viscosity &viscosity),(alpha, rho, U, alphaRhoPhi, phi, viscosity)) | |
RASModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity) | |
Construct from components. More... | |
RASModel (const RASModel &)=delete | |
Disallow default bitwise copy construction. More... | |
virtual | ~RASModel () |
Destructor. More... | |
const dimensionedScalar & | kMin () const |
Return the lower allowable limit for k (default: small) More... | |
const dimensionedScalar & | epsilonMin () const |
Return the lower allowable limit for epsilon (default: small) More... | |
const dimensionedScalar & | omegaMin () const |
Return the lower allowable limit for omega (default: small) More... | |
dimensionedScalar & | kMin () |
Allow kMin to be changed. More... | |
dimensionedScalar & | epsilonMin () |
Allow epsilonMin to be changed. More... | |
dimensionedScalar & | omegaMin () |
Allow omegaMin to be changed. More... | |
virtual const dictionary & | coeffDict () const |
Const access to the coefficients dictionary. More... | |
virtual tmp< volScalarField > | nu () const |
Return the laminar viscosity. More... | |
virtual tmp< scalarField > | nu (const label patchi) const |
Return the laminar viscosity on patchi. More... | |
virtual tmp< volScalarField > | nuEff () const |
Return the effective viscosity. More... | |
virtual tmp< scalarField > | nuEff (const label patchi) const |
Return the effective viscosity on patch. More... | |
void | operator= (const RASModel &)=delete |
Disallow default bitwise assignment. More... | |
Protected Member Functions | |
virtual void | correctNut () |
virtual tmp< fvScalarMatrix > | kSource () const |
virtual tmp< fvScalarMatrix > | epsilonSource () const |
Protected Member Functions inherited from RASModel< BasicMomentumTransportModel > | |
virtual void | printCoeffs (const word &type) |
Print model coefficients. 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< BasicMomentumTransportModel > > | |
volScalarField | nut_ |
Protected Attributes inherited from RASModel< BasicMomentumTransportModel > | |
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... | |
autoPtr< laminarModels::generalisedNewtonianViscosityModel > | viscosityModel_ |
Run-time selectable generalised Newtonian viscosity model. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from RASModel< BasicMomentumTransportModel > | |
static autoPtr< RASModel > | New (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity) |
Return a reference to the selected RAS model. More... | |
Renormalisation group k-epsilon turbulence model for incompressible and compressible flows.
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; sigmak 0.71942; sigmaEps 0.71942; eta0 4.38; beta 0.012; }
Definition at line 82 of file RNGkEpsilon.H.
typedef BasicMomentumTransportModel::alphaField alphaField |
Definition at line 122 of file RNGkEpsilon.H.
typedef BasicMomentumTransportModel::rhoField rhoField |
Definition at line 123 of file RNGkEpsilon.H.
RNGkEpsilon | ( | const alphaField & | alpha, |
const rhoField & | rho, | ||
const volVectorField & | U, | ||
const surfaceScalarField & | alphaRhoPhi, | ||
const surfaceScalarField & | phi, | ||
const viscosity & | viscosity, | ||
const word & | type = typeName |
||
) |
Construct from components.
Definition at line 84 of file RNGkEpsilon.C.
References Foam::bound().
Referenced by RNGkEpsilon< BasicMomentumTransportModel >::epsilonSource().
|
delete |
Disallow default bitwise copy construction.
|
inlinevirtual |
Destructor.
Definition at line 149 of file RNGkEpsilon.H.
References RNGkEpsilon< BasicMomentumTransportModel >::read().
|
protectedvirtual |
Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.
Definition at line 41 of file RNGkEpsilon.C.
References dictionary::New(), and Foam::sqr().
|
protectedvirtual |
Definition at line 50 of file RNGkEpsilon.C.
References Foam::dimTime, and Foam::dimVolume.
|
protectedvirtual |
Definition at line 66 of file RNGkEpsilon.C.
References Foam::dimTime, Foam::dimVolume, and RNGkEpsilon< BasicMomentumTransportModel >::RNGkEpsilon().
TypeName | ( | "RNGkEpsilon< BasicMomentumTransportModel >" | ) |
Runtime type information.
|
virtual |
Re-read model coefficients if they have changed.
Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.
Definition at line 216 of file RNGkEpsilon.C.
References Foam::read().
Referenced by RNGkEpsilon< BasicMomentumTransportModel >::~RNGkEpsilon().
|
inline |
Return the effective diffusivity for k.
Definition at line 159 of file RNGkEpsilon.H.
References GeometricField< scalar, fvPatchField, volMesh >::New(), RASModel< BasicMomentumTransportModel >::nu(), and eddyViscosity< RASModel< BasicMomentumTransportModel > >::nut_.
|
inline |
Return the effective diffusivity for epsilon.
Definition at line 169 of file RNGkEpsilon.H.
References GeometricField< scalar, fvPatchField, volMesh >::New(), RASModel< BasicMomentumTransportModel >::nu(), and eddyViscosity< RASModel< BasicMomentumTransportModel > >::nut_.
|
inlinevirtual |
Return the turbulence kinetic energy.
Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.
Definition at line 179 of file RNGkEpsilon.H.
References RNGkEpsilon< BasicMomentumTransportModel >::k_.
|
inlinevirtual |
Return the turbulence kinetic energy dissipation rate.
Definition at line 185 of file RNGkEpsilon.H.
References RNGkEpsilon< BasicMomentumTransportModel >::epsilon_.
|
inlinevirtual |
Return the turbulence specific dissipation rate.
Definition at line 191 of file RNGkEpsilon.H.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), RNGkEpsilon< BasicMomentumTransportModel >::correct(), GeometricField< scalar, fvPatchField, volMesh >::New(), RNGkEpsilon< BasicMomentumTransportModel >::operator=(), and GeometricBoundaryField< Type, PatchField, GeoMesh >::types().
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Implements eddyViscosity< RASModel< BasicMomentumTransportModel > >.
Definition at line 239 of file RNGkEpsilon.C.
References Foam::fvc::absolute(), alpha(), Foam::bound(), tmp< T >::clear(), fvConstraints::constrain(), correct, Foam::fvm::ddt(), Foam::dev(), Foam::fvm::div(), Foam::fvc::div(), divU, fvConstraints, fvModels, Foam::constant::universal::G, Foam::fvc::grad(), Foam::fvm::laplacian(), Foam::mag(), dictionary::New(), nut, phi, R, tmp< T >::ref(), rho, Foam::solve(), Foam::fvm::Sp(), Foam::sqr(), Foam::sqrt(), Foam::fvm::SuSp(), and Foam::twoSymm().
Referenced by RNGkEpsilon< BasicMomentumTransportModel >::omega().
|
delete |
Disallow default bitwise assignment.
Referenced by RNGkEpsilon< BasicMomentumTransportModel >::omega().
|
protected |
Definition at line 97 of file RNGkEpsilon.H.
|
protected |
Definition at line 98 of file RNGkEpsilon.H.
|
protected |
Definition at line 99 of file RNGkEpsilon.H.
|
protected |
Definition at line 100 of file RNGkEpsilon.H.
|
protected |
Definition at line 101 of file RNGkEpsilon.H.
|
protected |
Definition at line 102 of file RNGkEpsilon.H.
|
protected |
Definition at line 103 of file RNGkEpsilon.H.
|
protected |
Definition at line 104 of file RNGkEpsilon.H.
|
protected |
Definition at line 109 of file RNGkEpsilon.H.
Referenced by RNGkEpsilon< BasicMomentumTransportModel >::k().
|
protected |
Definition at line 110 of file RNGkEpsilon.H.
Referenced by RNGkEpsilon< BasicMomentumTransportModel >::epsilon().