Abstract base class for incompressible turbulence models (RAS, LES and laminar). More...
#include <incompressibleTurbulenceModel/turbulenceModel.H>
Abstract base class for incompressible turbulence models (RAS, LES and laminar).
Definition at line 68 of file turbulenceModel.H.
Public Member Functions | |
TypeName ("turbulenceModel") | |
Runtime type information.
| |
declareRunTimeNewSelectionTable (autoPtr, turbulenceModel, turbulenceModel,(const volVectorField &U, const surfaceScalarField &phi, transportModel &lamTransportModel),(U, phi, lamTransportModel)) | |
turbulenceModel (const volVectorField &U, const surfaceScalarField &phi, transportModel &lamTransportModel) | |
Construct from components.
| |
virtual | ~turbulenceModel () |
Destructor.
| |
const volVectorField & | U () const |
Access function to velocity field.
| |
const surfaceScalarField & | phi () const |
Access function to flux field.
| |
transportModel & | transport () const |
Access function to incompressible transport model.
| |
const volScalarField & | nu () const |
Return the laminar viscosity.
| |
virtual tmp< volScalarField > | nut () const =0 |
Return the turbulence viscosity.
| |
virtual tmp< volScalarField > | nuEff () const =0 |
Return the effective viscosity.
| |
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 > | devReff () const =0 |
Return the effective stress tensor including the laminar stress.
| |
virtual tmp< fvVectorMatrix > | divDevReff (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 volVectorField &U, const surfaceScalarField &phi, transportModel &lamTransportModel) |
Return a reference to the selected turbulence model.
| |
Protected Attributes | |
const Time & | runTime_ |
const fvMesh & | mesh_ |
const volVectorField & | U_ |
const surfaceScalarField & | phi_ |
transportModel & | transportModel_ |
turbulenceModel | ( | const volVectorField & | U, |
const surfaceScalarField & | phi, | ||
transportModel & | lamTransportModel | ||
) |
Construct from components.
Definition at line 45 of file turbulenceModel.C.
virtual ~turbulenceModel | ( | ) | [inline, virtual]
|
Destructor.
Definition at line 140 of file turbulenceModel.H.
TypeName | ( | "turbulenceModel" | ) |
Runtime type information.
declareRunTimeNewSelectionTable | ( | autoPtr | , |
turbulenceModel | , | ||
turbulenceModel | , | ||
(const volVectorField &U, const surfaceScalarField &phi, transportModel &lamTransportModel) | , | ||
(U, phi, lamTransportModel) | |||
) |
autoPtr< turbulenceModel > New | ( | const volVectorField & | U, |
const surfaceScalarField & | phi, | ||
transportModel & | lamTransportModel | ||
) | [static]
|
Return a reference to the selected turbulence model.
Definition at line 63 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, IOobject::time(), and U.
const volVectorField& U | ( | ) | const [inline]
|
Access function to velocity field.
Definition at line 147 of file turbulenceModel.H.
References incompressible::turbulenceModel::U_.
Referenced by Smagorinsky2::B(), scaleSimilarity::B(), kOmegaSSTSAS::B(), laminar::B(), SpalartAllmaras::B(), GenEddyVisc::B(), locDynOneEqEddy::correct(), kOmegaSSTSAS::devBeff(), laminar::devBeff(), SpalartAllmaras::devBeff(), GenSGSStress::devBeff(), GenEddyVisc::devBeff(), scaleSimilarity::epsilon(), kOmegaSSTSAS::epsilon(), laminar::epsilon(), SpalartAllmaras::epsilon(), spectEddyVisc::k(), scaleSimilarity::k(), laminar::k(), SpalartAllmaras::k(), Smagorinsky::k(), and incompressible::fixedShearStressFvPatchVectorField::updateCoeffs().
const surfaceScalarField& phi | ( | ) | const [inline]
|
Access function to flux field.
Definition at line 153 of file turbulenceModel.H.
References incompressible::turbulenceModel::phi_.
Referenced by LRRDiffStress::correct(), locDynOneEqEddy::correct(), kOmegaSSTSAS::correct(), SpalartAllmaras::correct(), oneEqEddy::correct(), dynOneEqEddy::correct(), and DeardorffDiffStress::correct().
transportModel& transport | ( | ) | const [inline]
|
Access function to incompressible transport model.
Definition at line 159 of file turbulenceModel.H.
References incompressible::turbulenceModel::transportModel_.
const volScalarField& nu | ( | ) | const [inline]
|
Return the laminar viscosity.
Definition at line 165 of file turbulenceModel.H.
References transportModel::nu(), and incompressible::turbulenceModel::transportModel_.
Referenced by qZeta::correct(), LienLeschzinerLowRe::correct(), LienCubicKELowRe::correct(), LamBremhorstKE::correct(), realizableKE::correct(), LaunderSharmaKE::correct(), SpalartAllmaras::correct(), LRRDiffStress::DBEff(), DeardorffDiffStress::DBEff(), NonlinearKEShih::DepsilonEff(), LienLeschzinerLowRe::DepsilonEff(), LienCubicKELowRe::DepsilonEff(), LienCubicKE::DepsilonEff(), LamBremhorstKE::DepsilonEff(), RNGkEpsilon::DepsilonEff(), realizableKE::DepsilonEff(), LRR::DepsilonEff(), LaunderSharmaKE::DepsilonEff(), LaunderGibsonRSTM::DepsilonEff(), kEpsilon::DepsilonEff(), laminar::devBeff(), GenSGSStress::devBeff(), LRR::devReff(), LaunderGibsonRSTM::devReff(), laminar::devReff(), laminar::divDevBeff(), NonlinearKEShih::DkEff(), LienLeschzinerLowRe::DkEff(), LienCubicKELowRe::DkEff(), LienCubicKE::DkEff(), LamBremhorstKE::DkEff(), kOmega::DkEff(), locDynOneEqEddy::DkEff(), kOmegaSSTSAS::DkEff(), RNGkEpsilon::DkEff(), realizableKE::DkEff(), LaunderSharmaKE::DkEff(), kOmegaSST::DkEff(), kEpsilon::DkEff(), oneEqEddy::DkEff(), dynOneEqEddy::DkEff(), SpalartAllmaras::DnuTildaEff(), kOmega::DomegaEff(), kOmegaSSTSAS::DomegaEff(), kOmegaSST::DomegaEff(), qZeta::DqEff(), LRR::DREff(), LaunderGibsonRSTM::DREff(), qZeta::DzetaEff(), laminar::epsilon(), kOmegaSSTSAS::F1(), kOmegaSSTSAS::F2(), SpalartAllmaras::fv1(), SpalartAllmaras::fv2(), SpalartAllmaras::fv3(), spectEddyVisc::k(), laminar::nuEff(), laminar::nuSgs(), and laminar::nut().
virtual tmp<volScalarField> nut | ( | ) | const [pure virtual]
|
Return the turbulence viscosity.
Implemented in incompressible::LESModel, kEpsilon, kOmega, kOmegaSST, LamBremhorstKE, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LienCubicKE, LienCubicKELowRe, LienLeschzinerLowRe, LRR, NonlinearKEShih, qZeta, incompressible::RASModel, realizableKE, RNGkEpsilon, SpalartAllmaras, and incompressible::laminar.
virtual tmp<volScalarField> nuEff | ( | ) | const [pure virtual]
|
Return the effective viscosity.
Implemented in laminar, incompressible::LESModel, laminar, incompressible::RASModel, and incompressible::laminar.
virtual tmp<volScalarField> k | ( | ) | const [pure virtual]
|
Return the turbulence kinetic energy.
Implemented in dynOneEqEddy, GenEddyVisc, GenSGSStress, homogeneousDynSmagorinsky, kOmegaSSTSAS, laminar, incompressible::LESModel, locDynOneEqEddy, mixedSmagorinsky, oneEqEddy, scaleSimilarity, Smagorinsky, SpalartAllmaras, spectEddyVisc, kEpsilon, kOmega, kOmegaSST, LamBremhorstKE, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LienCubicKE, LienCubicKELowRe, LienLeschzinerLowRe, LRR, NonlinearKEShih, qZeta, incompressible::RASModel, realizableKE, RNGkEpsilon, SpalartAllmaras, and incompressible::laminar.
virtual tmp<volScalarField> epsilon | ( | ) | const [pure virtual]
|
Return the turbulence kinetic energy dissipation rate.
Implemented in GenEddyVisc, GenSGSStress, kOmegaSSTSAS, laminar, incompressible::LESModel, mixedSmagorinsky, scaleSimilarity, SpalartAllmaras, kEpsilon, kOmega, kOmegaSST, LamBremhorstKE, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LienCubicKE, LienCubicKELowRe, LienLeschzinerLowRe, LRR, NonlinearKEShih, qZeta, incompressible::RASModel, realizableKE, RNGkEpsilon, SpalartAllmaras, and incompressible::laminar.
virtual tmp<volSymmTensorField> R | ( | ) | const [pure virtual]
|
Return the Reynolds stress tensor.
Implemented in incompressible::LESModel, kEpsilon, kOmega, kOmegaSST, LamBremhorstKE, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LienCubicKE, LienCubicKELowRe, LienLeschzinerLowRe, LRR, NonlinearKEShih, qZeta, incompressible::RASModel, realizableKE, RNGkEpsilon, SpalartAllmaras, and incompressible::laminar.
virtual tmp<volSymmTensorField> devReff | ( | ) | const [pure virtual]
|
Return the effective stress tensor including the laminar stress.
Implemented in incompressible::LESModel, kEpsilon, kOmega, kOmegaSST, LamBremhorstKE, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LienCubicKE, LienCubicKELowRe, LienLeschzinerLowRe, LRR, NonlinearKEShih, qZeta, incompressible::RASModel, realizableKE, RNGkEpsilon, SpalartAllmaras, and incompressible::laminar.
virtual tmp<fvVectorMatrix> divDevReff | ( | 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 incompressible::LESModel, kEpsilon, kOmega, kOmegaSST, LamBremhorstKE, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LienCubicKE, LienCubicKELowRe, LienLeschzinerLowRe, LRR, NonlinearKEShih, qZeta, incompressible::RASModel, realizableKE, RNGkEpsilon, SpalartAllmaras, and incompressible::laminar.
Definition at line 114 of file turbulenceModel.C.
virtual bool read | ( | ) | [pure virtual]
|
Read turbulenceProperties dictionary.
Implemented in DeardorffDiffStress, dynOneEqEddy, GenEddyVisc, GenSGSStress, homogeneousDynSmagorinsky, kOmegaSSTSAS, laminar, incompressible::LESModel, locDynOneEqEddy, LRRDiffStress, mixedSmagorinsky, oneEqEddy, scaleSimilarity, Smagorinsky, Smagorinsky2, SpalartAllmaras, SpalartAllmarasIDDES, spectEddyVisc, kEpsilon, kOmega, kOmegaSST, LamBremhorstKE, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LienCubicKE, LienCubicKELowRe, LienLeschzinerLowRe, LRR, NonlinearKEShih, qZeta, incompressible::RASModel, realizableKE, RNGkEpsilon, SpalartAllmaras, and incompressible::laminar.
const Time& runTime_ [protected]
|
Definition at line 75 of file turbulenceModel.H.
Referenced by laminar::B(), GenSGSStress::devBeff(), qZeta::devReff(), NonlinearKEShih::devReff(), LienLeschzinerLowRe::devReff(), LienCubicKELowRe::devReff(), LienCubicKE::devReff(), LamBremhorstKE::devReff(), kOmega::devReff(), RNGkEpsilon::devReff(), realizableKE::devReff(), LRR::devReff(), LaunderSharmaKE::devReff(), LaunderGibsonRSTM::devReff(), laminar::devReff(), kOmegaSST::devReff(), kEpsilon::devReff(), SpalartAllmaras::devReff(), laminar::epsilon(), SpalartAllmaras::epsilon(), laminar::k(), SpalartAllmaras::k(), laminar::nuSgs(), laminar::nut(), qZeta::R(), NonlinearKEShih::R(), LienLeschzinerLowRe::R(), LienCubicKELowRe::R(), LienCubicKE::R(), LamBremhorstKE::R(), kOmega::R(), RNGkEpsilon::R(), realizableKE::R(), LaunderSharmaKE::R(), laminar::R(), kOmegaSST::R(), kEpsilon::R(), and SpalartAllmaras::R().
const fvMesh& mesh_ [protected]
|
Definition at line 76 of file turbulenceModel.H.
Referenced by laminar::B(), LienLeschzinerLowRe::correct(), LienCubicKELowRe::correct(), LamBremhorstKE::correct(), kOmegaSSTSAS::correct(), LRR::correct(), LaunderGibsonRSTM::correct(), kOmegaSST::correct(), SpalartAllmaras::correct(), GenSGSStress::devBeff(), qZeta::devReff(), NonlinearKEShih::devReff(), LienLeschzinerLowRe::devReff(), LienCubicKELowRe::devReff(), LienCubicKE::devReff(), LamBremhorstKE::devReff(), kOmega::devReff(), RNGkEpsilon::devReff(), realizableKE::devReff(), LRR::devReff(), LaunderSharmaKE::devReff(), LaunderGibsonRSTM::devReff(), laminar::devReff(), kOmegaSST::devReff(), kEpsilon::devReff(), SpalartAllmaras::devReff(), kOmega::epsilon(), laminar::epsilon(), kOmegaSST::epsilon(), SpalartAllmaras::epsilon(), laminar::k(), SpalartAllmaras::k(), laminar::nuSgs(), laminar::nut(), qZeta::R(), NonlinearKEShih::R(), LienLeschzinerLowRe::R(), LienCubicKELowRe::R(), LienCubicKE::R(), LamBremhorstKE::R(), kOmega::R(), RNGkEpsilon::R(), realizableKE::R(), LaunderSharmaKE::R(), laminar::R(), kOmegaSST::R(), kEpsilon::R(), and SpalartAllmaras::R().
const volVectorField& U_ [protected]
|
Definition at line 78 of file turbulenceModel.H.
Referenced by qZeta::correct(), NonlinearKEShih::correct(), LienLeschzinerLowRe::correct(), LienCubicKELowRe::correct(), LienCubicKE::correct(), LamBremhorstKE::correct(), kOmega::correct(), RNGkEpsilon::correct(), realizableKE::correct(), LRR::correct(), LaunderSharmaKE::correct(), LaunderGibsonRSTM::correct(), kOmegaSST::correct(), kEpsilon::correct(), SpalartAllmaras::correct(), qZeta::devReff(), NonlinearKEShih::devReff(), LienLeschzinerLowRe::devReff(), LienCubicKELowRe::devReff(), LienCubicKE::devReff(), LamBremhorstKE::devReff(), kOmega::devReff(), RNGkEpsilon::devReff(), realizableKE::devReff(), LRR::devReff(), LaunderSharmaKE::devReff(), LaunderGibsonRSTM::devReff(), laminar::devReff(), kOmegaSST::devReff(), kEpsilon::devReff(), SpalartAllmaras::devReff(), laminar::epsilon(), laminar::k(), qZeta::R(), NonlinearKEShih::R(), LienLeschzinerLowRe::R(), LamBremhorstKE::R(), kOmega::R(), RNGkEpsilon::R(), realizableKE::R(), LaunderSharmaKE::R(), laminar::R(), kOmegaSST::R(), kEpsilon::R(), SpalartAllmaras::R(), and incompressible::turbulenceModel::U().
const surfaceScalarField& phi_ [protected]
|
Definition at line 79 of file turbulenceModel.H.
Referenced by qZeta::correct(), NonlinearKEShih::correct(), LienLeschzinerLowRe::correct(), LienCubicKELowRe::correct(), LienCubicKE::correct(), LamBremhorstKE::correct(), kOmega::correct(), RNGkEpsilon::correct(), realizableKE::correct(), LRR::correct(), LaunderSharmaKE::correct(), LaunderGibsonRSTM::correct(), kOmegaSST::correct(), kEpsilon::correct(), SpalartAllmaras::correct(), and incompressible::turbulenceModel::phi().
transportModel& transportModel_ [protected]
|
Definition at line 81 of file turbulenceModel.H.
Referenced by incompressible::turbulenceModel::nu(), and incompressible::turbulenceModel::transport().