Realizable k-epsilon turbulence model for incompressible and compressible flows. More...
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 ("realizableKE") | |
Runtime type information. More... | |
realizableKE (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 | ~realizableKE () |
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 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 () |
Destructor. More... | |
volScalarField & | evNut () |
Return non-const access to the turbulence viscosity. 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 > | R () 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 () |
Destructor. More... | |
virtual tmp< volSymmTensorField > | devRhoReff () const |
Return the effective stress tensor. More... | |
virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
Return the source term for the momentum equation. More... | |
virtual tmp< fvVectorMatrix > | divDevRhoReff (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 () |
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 > | nuEff () const |
Return the effective viscosity. More... | |
virtual tmp< scalarField > | nuEff (const label patchi) const |
Return the effective viscosity on patch. More... | |
Protected Member Functions | |
tmp< volScalarField > | rCmu (const volTensorField &gradU, const volScalarField &S2, const volScalarField &magS) |
virtual void | correctNut (const volTensorField &gradU, const volScalarField &S2, const volScalarField &magS) |
virtual void | correctNut () |
virtual tmp< fvScalarMatrix > | kSource () const |
virtual tmp< fvScalarMatrix > | epsilonSource () const |
Protected Member Functions inherited from RASModel< BasicTurbulenceModel > | |
virtual void | printCoeffs (const word &type) |
Print model coefficients. More... | |
Protected Attributes | |
dimensionedScalar | Cmu_ |
dimensionedScalar | A0_ |
dimensionedScalar | C2_ |
dimensionedScalar | sigmak_ |
dimensionedScalar | sigmaEps_ |
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< RASModel > | New (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... | |
Realizable k-epsilon turbulence model for incompressible and compressible flows.
Shih, T. H., Liou, W. W., Shabbir, A., Yang, Z., & Zhu, J. (1994). A new k-epsilon eddy viscosity model for high Reynolds number turbulent flows: Model development and validation. NASA STI/Recon Technical Report N, 95, 11442. Shih, T. H., Liou, W. W., Shabbir, A., Yang, Z., & Zhu, J. (1995). A New k-epsilon Eddy Viscosity Model for High Reynolds Number Turbulent Flows. Computers and Fluids, 24(3), 227-238.
The default model coefficients are
realizableKECoeffs { Cmu 0.09; A0 4.0; C2 1.9; sigmak 1.0; sigmaEps 1.2; }
Definition at line 81 of file realizableKE.H.
typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 128 of file realizableKE.H.
typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 129 of file realizableKE.H.
typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 130 of file realizableKE.H.
realizableKE | ( | 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 132 of file realizableKE.C.
References Foam::bound().
Referenced by realizableKE< BasicTurbulenceModel >::epsilonSource().
|
inlinevirtual |
Destructor.
Definition at line 153 of file realizableKE.H.
References realizableKE< BasicTurbulenceModel >::read().
|
protected |
Definition at line 41 of file realizableKE.C.
References Foam::acos(), tmp< T >::clear(), realizableKE< BasicTurbulenceModel >::correctNut(), Foam::cos(), Foam::dev(), Foam::magSqr(), Foam::max(), Foam::min(), Foam::skew(), Foam::sqrt(), and Foam::symm().
|
protectedvirtual |
Definition at line 74 of file realizableKE.C.
References optionList::correct(), options::New(), and Foam::sqr().
|
protectedvirtual |
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 89 of file realizableKE.C.
References Foam::dev(), Foam::fvc::grad(), Foam::magSqr(), Foam::sqrt(), and Foam::symm().
Referenced by realizableKE< BasicTurbulenceModel >::rCmu().
|
protectedvirtual |
Definition at line 99 of file realizableKE.C.
References Foam::dimTime, and Foam::dimVolume.
|
protectedvirtual |
Definition at line 114 of file realizableKE.C.
References Foam::dimTime, Foam::dimVolume, and realizableKE< BasicTurbulenceModel >::realizableKE().
TypeName | ( | "realizableKE< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Re-read model coefficients if they have changed.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 239 of file realizableKE.C.
References Foam::read().
Referenced by realizableKE< BasicTurbulenceModel >::~realizableKE().
|
inline |
Return the effective diffusivity for k.
Definition at line 163 of file realizableKE.H.
References nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.
|
inline |
Return the effective diffusivity for epsilon.
Definition at line 176 of file realizableKE.H.
References nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.
|
inlinevirtual |
Return the turbulence kinetic energy.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 189 of file realizableKE.H.
References realizableKE< BasicTurbulenceModel >::k_.
|
inlinevirtual |
Return the turbulence kinetic energy dissipation rate.
Definition at line 195 of file realizableKE.H.
References realizableKE< BasicTurbulenceModel >::correct(), and realizableKE< BasicTurbulenceModel >::epsilon_.
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 259 of file realizableKE.C.
References Foam::fvc::absolute(), Foam::bound(), optionList::constrain(), correct, optionList::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::magSqr(), Foam::max(), options::New(), nu, nut, phi, tmp< T >::ref(), Foam::solve(), Foam::fvm::Sp(), Foam::sqrt(), Foam::fvm::SuSp(), Foam::symm(), and Foam::twoSymm().
Referenced by realizableKE< BasicTurbulenceModel >::epsilon().
|
protected |
Definition at line 92 of file realizableKE.H.
|
protected |
Definition at line 93 of file realizableKE.H.
|
protected |
Definition at line 94 of file realizableKE.H.
|
protected |
Definition at line 95 of file realizableKE.H.
|
protected |
Definition at line 96 of file realizableKE.H.
|
protected |
Definition at line 101 of file realizableKE.H.
Referenced by realizableKE< BasicTurbulenceModel >::k().
|
protected |
Definition at line 102 of file realizableKE.H.
Referenced by realizableKE< BasicTurbulenceModel >::epsilon().