Abstract base class for compressible turbulence models (RAS, LES and laminar). More...
#include <compressibleTurbulenceModel/turbulenceModel.H>
Abstract base class for compressible turbulence models (RAS, LES and laminar).
Definition at line 69 of file turbulenceModel.H.
Public Member Functions | |
TypeName ("turbulenceModel") | |
Runtime type information.
| |
declareRunTimeNewSelectionTable (autoPtr, turbulenceModel, turbulenceModel,(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const basicThermo &thermoPhysicalModel),(rho, U, phi, thermoPhysicalModel)) | |
turbulenceModel (const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const basicThermo &thermoPhysicalModel) | |
Construct from components.
| |
virtual | ~turbulenceModel () |
Destructor.
| |
const volScalarField & | rho () const |
Access function to density field.
| |
const volVectorField & | U () const |
Access function to velocity field.
| |
const surfaceScalarField & | phi () const |
Access function to flux field.
| |
const basicThermo & | thermo () const |
Access function to thermophysical model.
| |
const volScalarField & | mu () const |
Return the laminar viscosity.
| |
const volScalarField & | alpha () const |
Return the laminar thermal conductivity.
| |
virtual tmp< volScalarField > | mut () const =0 |
Return the turbulence viscosity.
| |
virtual tmp< volScalarField > | muEff () const =0 |
Return the effective viscosity.
| |
virtual tmp< volScalarField > | alphaEff () const =0 |
Return the effective turbulent thermal diffusivity.
| |
virtual tmp< volScalarField > | k () const =0 |
Return the turbulence kinetic energy.
| |
virtual tmp< volScalarField > | epsilon () const =0 |
Return the turbulence kinetic energy dissipation rate.
| |
virtual tmp< volSymmTensorField > | R () const =0 |
Return the Reynolds stress tensor.
| |
virtual tmp< volSymmTensorField > | devRhoReff () const =0 |
Return the effective stress tensor including the laminar stress.
| |
virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const =0 |
Return the source term for the momentum equation.
| |
virtual void | correct ()=0 |
Solve the turbulence equations and correct the turbulence viscosity.
| |
virtual bool | read ()=0 |
Read turbulenceProperties dictionary.
| |
Static Public Member Functions | |
static autoPtr< turbulenceModel > | New (const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const basicThermo &thermoPhysicalModel) |
Return a reference to the selected turbulence model.
| |
Protected Attributes | |
const Time & | runTime_ |
const fvMesh & | mesh_ |
const volScalarField & | rho_ |
const volVectorField & | U_ |
const surfaceScalarField & | phi_ |
const basicThermo & | thermophysicalModel_ |
turbulenceModel | ( | const volScalarField & | rho, |
const volVectorField & | U, | ||
const surfaceScalarField & | phi, | ||
const basicThermo & | thermoPhysicalModel | ||
) |
Construct from components.
Definition at line 45 of file turbulenceModel.C.
virtual ~turbulenceModel | ( | ) | [inline, virtual]
|
Destructor.
Definition at line 145 of file turbulenceModel.H.
TypeName | ( | "turbulenceModel" | ) |
Runtime type information.
declareRunTimeNewSelectionTable | ( | autoPtr | , |
turbulenceModel | , | ||
turbulenceModel | , | ||
(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const basicThermo &thermoPhysicalModel) | , | ||
(rho, U, phi, thermoPhysicalModel) | |||
) |
autoPtr< turbulenceModel > New | ( | const volScalarField & | rho, |
const volVectorField & | U, | ||
const surfaceScalarField & | phi, | ||
const basicThermo & | thermoPhysicalModel | ||
) | [static]
|
Return a reference to the selected turbulence model.
Definition at line 65 of file turbulenceModel.C.
References TimePaths::constant(), IOobject::db(), Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorIn, Foam::Info, dictionary::lookup(), IOobject::MUST_READ, IOobject::NO_WRITE, phi, rho, IOobject::time(), and U.
const volScalarField& rho | ( | ) | const [inline]
|
Access function to density field.
Definition at line 152 of file turbulenceModel.H.
References turbulenceModel::rho_.
Referenced by mutWallFunctionFvPatchScalarField::calcMut(), mutSpalartAllmarasWallFunctionFvPatchScalarField::calcMut(), mutSpalartAllmarasStandardRoughWallFunctionFvPatchScalarField::calcMut(), mutRoughWallFunctionFvPatchScalarField::calcMut(), mutkWallFunctionFvPatchScalarField::calcMut(), mutSpalartAllmarasWallFunctionFvPatchScalarField::calcUTau(), mutSpalartAllmarasStandardWallFunctionFvPatchScalarField::calcYPlus(), mutSpalartAllmarasStandardRoughWallFunctionFvPatchScalarField::calcYPlus(), alphaSgsJayatillekeWallFunctionFvPatchScalarField::evaluate(), mutWallFunctionFvPatchScalarField::yPlus(), mutSpalartAllmarasWallFunctionFvPatchScalarField::yPlus(), and mutkWallFunctionFvPatchScalarField::yPlus().
const volVectorField& U | ( | ) | const [inline]
|
Access function to velocity field.
Definition at line 158 of file turbulenceModel.H.
References turbulenceModel::U_.
Referenced by mutSpalartAllmarasWallFunctionFvPatchScalarField::calcMut(), mutSpalartAllmarasStandardWallFunctionFvPatchScalarField::calcMut(), mutSpalartAllmarasStandardRoughWallFunctionFvPatchScalarField::calcMut(), mutSpalartAllmarasWallFunctionFvPatchScalarField::calcUTau(), alphaSgsJayatillekeWallFunctionFvPatchScalarField::evaluate(), mutSpalartAllmarasWallFunctionFvPatchScalarField::yPlus(), mutSpalartAllmarasStandardWallFunctionFvPatchScalarField::yPlus(), and mutSpalartAllmarasStandardRoughWallFunctionFvPatchScalarField::yPlus().
const surfaceScalarField& phi | ( | ) | const [inline]
|
Access function to flux field.
Definition at line 164 of file turbulenceModel.H.
References turbulenceModel::phi_.
const basicThermo& thermo | ( | ) | const [inline]
|
Access function to thermophysical model.
Definition at line 170 of file turbulenceModel.H.
References turbulenceModel::thermophysicalModel_.
Referenced by turbulentHeatFluxTemperatureFvPatchScalarField::updateCoeffs().
const volScalarField& mu | ( | ) | const [inline]
|
Return the laminar viscosity.
Definition at line 176 of file turbulenceModel.H.
References basicThermo::mu(), and turbulenceModel::thermophysicalModel_.
Referenced by mutWallFunctionFvPatchScalarField::calcMut(), mutSpalartAllmarasWallFunctionFvPatchScalarField::calcMut(), mutSpalartAllmarasStandardWallFunctionFvPatchScalarField::calcMut(), mutSpalartAllmarasStandardRoughWallFunctionFvPatchScalarField::calcMut(), mutRoughWallFunctionFvPatchScalarField::calcMut(), mutkWallFunctionFvPatchScalarField::calcMut(), mutSpalartAllmarasWallFunctionFvPatchScalarField::calcUTau(), mutSpalartAllmarasStandardWallFunctionFvPatchScalarField::calcYPlus(), mutSpalartAllmarasStandardRoughWallFunctionFvPatchScalarField::calcYPlus(), alphaSgsJayatillekeWallFunctionFvPatchScalarField::evaluate(), mutWallFunctionFvPatchScalarField::yPlus(), mutSpalartAllmarasWallFunctionFvPatchScalarField::yPlus(), and mutkWallFunctionFvPatchScalarField::yPlus().
const volScalarField& alpha | ( | ) | const [inline]
|
Return the laminar thermal conductivity.
Definition at line 182 of file turbulenceModel.H.
References basicThermo::alpha(), and turbulenceModel::thermophysicalModel_.
Referenced by RNGkEpsilon::alphaEff(), realizableKE::alphaEff(), LRR::alphaEff(), LaunderSharmaKE::alphaEff(), LaunderGibsonRSTM::alphaEff(), laminar::alphaEff(), kOmegaSST::alphaEff(), kEpsilon::alphaEff(), SpalartAllmaras::alphaEff(), GenSGSStress::alphaEff(), GenEddyVisc::alphaEff(), PDRkEpsilon::alphaEff(), and alphaSgsJayatillekeWallFunctionFvPatchScalarField::evaluate().
virtual tmp<volScalarField> mut | ( | ) | const [pure virtual]
|
Return the turbulence viscosity.
Implemented in LESModel, kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, RASModel, realizableKE, RNGkEpsilon, SpalartAllmaras, laminar, and PDRkEpsilon.
virtual tmp<volScalarField> muEff | ( | ) | const [pure virtual]
|
virtual tmp<volScalarField> alphaEff | ( | ) | const [pure virtual]
|
Return the effective turbulent thermal diffusivity.
Implemented in GenEddyVisc, GenSGSStress, LESModel, SpalartAllmaras, kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, RASModel, realizableKE, RNGkEpsilon, SpalartAllmaras, laminar, and PDRkEpsilon.
virtual tmp<volScalarField> k | ( | ) | const [pure virtual]
|
Return the turbulence kinetic energy.
Implemented in GenEddyVisc, GenSGSStress, LESModel, SpalartAllmaras, kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, RASModel, realizableKE, RNGkEpsilon, SpalartAllmaras, laminar, and PDRkEpsilon.
virtual tmp<volScalarField> epsilon | ( | ) | const [pure virtual]
|
Return the turbulence kinetic energy dissipation rate.
Implemented in GenEddyVisc, GenSGSStress, LESModel, SpalartAllmaras, kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, RASModel, realizableKE, RNGkEpsilon, SpalartAllmaras, laminar, and PDRkEpsilon.
virtual tmp<volSymmTensorField> R | ( | ) | const [pure virtual]
|
Return the Reynolds stress tensor.
Implemented in LESModel, kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, RASModel, realizableKE, RNGkEpsilon, SpalartAllmaras, laminar, and PDRkEpsilon.
virtual tmp<volSymmTensorField> devRhoReff | ( | ) | const [pure virtual]
|
Return the effective stress tensor including the laminar stress.
Implemented in LESModel, kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, RASModel, realizableKE, RNGkEpsilon, SpalartAllmaras, laminar, and PDRkEpsilon.
virtual tmp<fvVectorMatrix> divDevRhoReff | ( | volVectorField & | U ) | const [pure virtual]
|
Return the source term for the momentum equation.
void correct | ( | ) | [pure virtual]
|
Solve the turbulence equations and correct the turbulence viscosity.
Implemented in LESModel, kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, RASModel, realizableKE, RNGkEpsilon, SpalartAllmaras, laminar, and PDRkEpsilon.
Definition at line 121 of file turbulenceModel.C.
virtual bool read | ( | ) | [pure virtual]
|
Read turbulenceProperties dictionary.
Implemented in DeardorffDiffStress, dynOneEqEddy, GenEddyVisc, GenSGSStress, LESModel, lowReOneEqEddy, oneEqEddy, Smagorinsky, SpalartAllmaras, kEpsilon, kOmegaSST, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LRR, RASModel, realizableKE, RNGkEpsilon, SpalartAllmaras, laminar, and PDRkEpsilon.
const Time& runTime_ [protected]
|
Definition at line 76 of file turbulenceModel.H.
Referenced by SpalartAllmaras::epsilon(), and SpalartAllmaras::k().
const fvMesh& mesh_ [protected]
|
Definition at line 77 of file turbulenceModel.H.
Referenced by kOmegaSST::epsilon(), SpalartAllmaras::epsilon(), and SpalartAllmaras::k().
const volScalarField& rho_ [protected]
|
Definition at line 79 of file turbulenceModel.H.
Referenced by SpalartAllmaras::DnuTildaEff(), and turbulenceModel::rho().
const volVectorField& U_ [protected]
|
Definition at line 80 of file turbulenceModel.H.
Referenced by turbulenceModel::U().
const surfaceScalarField& phi_ [protected]
|
Definition at line 81 of file turbulenceModel.H.
Referenced by turbulenceModel::phi().
const basicThermo& thermophysicalModel_ [protected]
|
Definition at line 83 of file turbulenceModel.H.
Referenced by turbulenceModel::alpha(), turbulenceModel::mu(), and turbulenceModel::thermo().