Launder and Sharma low-Reynolds k-epsilon turbulence model for incompressible and compressible and combusting flows including rapid distortion theory (RDT) based compression term. 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 ("LaunderSharmaKE") | |
Runtime type information. More... | |
LaunderSharmaKE (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 | ~LaunderSharmaKE () |
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 > | fMu () const |
tmp< volScalarField > | f2 () const |
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 | C1_ |
dimensionedScalar | C2_ |
dimensionedScalar | C3_ |
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... | |
Launder and Sharma low-Reynolds k-epsilon turbulence model for incompressible and compressible and combusting flows including rapid distortion theory (RDT) based compression term.
Launder, B. E., & Sharma, B. I. (1974). Application of the energy-dissipation model of turbulence to the calculation of flow near a spinning disc. Letters in heat and mass transfer, 1(2), 131-137. 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
LaunderSharmaKECoeffs { Cmu 0.09; C1 1.44; C2 1.92; C3 -0.33; alphah 1.0; // only for compressible alphahk 1.0; // only for compressible alphaEps 0.76923; }
Definition at line 84 of file LaunderSharmaKE.H.
typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 127 of file LaunderSharmaKE.H.
typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 128 of file LaunderSharmaKE.H.
typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 129 of file LaunderSharmaKE.H.
LaunderSharmaKE | ( | 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 100 of file LaunderSharmaKE.C.
References Foam::bound().
|
inlinevirtual |
Destructor.
Definition at line 153 of file LaunderSharmaKE.H.
References LaunderSharmaKE< BasicTurbulenceModel >::read().
|
protected |
Definition at line 40 of file LaunderSharmaKE.C.
References Foam::exp(), nu, and Foam::sqr().
|
protected |
Definition at line 47 of file LaunderSharmaKE.C.
References Foam::exp(), Foam::min(), nu, and Foam::sqr().
|
protectedvirtual |
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 56 of file LaunderSharmaKE.C.
References optionList::correct(), options::New(), and Foam::sqr().
|
protectedvirtual |
Definition at line 67 of file LaunderSharmaKE.C.
References Foam::dimTime, and Foam::dimVolume.
|
protectedvirtual |
Definition at line 82 of file LaunderSharmaKE.C.
References Foam::dimTime, and Foam::dimVolume.
TypeName | ( | "LaunderSharmaKE< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Re-read model coefficients if they have changed.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 217 of file LaunderSharmaKE.C.
References Foam::read().
Referenced by LaunderSharmaKE< BasicTurbulenceModel >::~LaunderSharmaKE().
|
inline |
Return the effective diffusivity for k.
Definition at line 163 of file LaunderSharmaKE.H.
References nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.
|
inline |
Return the effective diffusivity for epsilon.
Definition at line 176 of file LaunderSharmaKE.H.
References nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.
|
inlinevirtual |
Return the turbulence kinetic energy.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 189 of file LaunderSharmaKE.H.
References LaunderSharmaKE< BasicTurbulenceModel >::k_.
|
inlinevirtual |
Return the turbulence kinetic energy dissipation rate.
Definition at line 195 of file LaunderSharmaKE.H.
References LaunderSharmaKE< BasicTurbulenceModel >::correct(), and LaunderSharmaKE< BasicTurbulenceModel >::epsilon_.
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 238 of file LaunderSharmaKE.C.
References Foam::fvc::absolute(), Foam::bound(), tmp< T >::clear(), 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::fvc::magSqrGradGrad(), options::New(), nu, nut, phi, tmp< T >::ref(), Foam::solve(), Foam::fvm::Sp(), Foam::sqrt(), Foam::fvm::SuSp(), and Foam::twoSymm().
Referenced by LaunderSharmaKE< BasicTurbulenceModel >::epsilon().
|
protected |
Definition at line 101 of file LaunderSharmaKE.H.
|
protected |
Definition at line 102 of file LaunderSharmaKE.H.
|
protected |
Definition at line 103 of file LaunderSharmaKE.H.
|
protected |
Definition at line 104 of file LaunderSharmaKE.H.
|
protected |
Definition at line 105 of file LaunderSharmaKE.H.
|
protected |
Definition at line 106 of file LaunderSharmaKE.H.
|
protected |
Definition at line 111 of file LaunderSharmaKE.H.
Referenced by LaunderSharmaKE< BasicTurbulenceModel >::k().
|
protected |
Definition at line 112 of file LaunderSharmaKE.H.
Referenced by LaunderSharmaKE< BasicTurbulenceModel >::epsilon().