Abstract base class for incompressible turbulence models. More...
#include <incompressibleRASModels/RASModel.H>
Abstract base class for incompressible turbulence models.
Definition at line 68 of file RASModel.H.
Public Member Functions | |
TypeName ("RASModel") | |
Runtime type information.
| |
declareRunTimeSelectionTable (autoPtr, RASModel, dictionary,(const volVectorField &U, const surfaceScalarField &phi, transportModel &lamTransportModel),(U, phi, lamTransportModel)) | |
RASModel (const word &type, const volVectorField &U, const surfaceScalarField &phi, transportModel &lamTransportModel) | |
Construct from components.
| |
virtual | ~RASModel () |
Destructor.
| |
const dimensionedScalar & | k0 () const |
Return the value of k0 which k is not allowed to be less than.
| |
const dimensionedScalar & | epsilon0 () const |
Return the value of epsilon0 which epsilon is not allowed to be.
| |
const dimensionedScalar & | epsilonSmall () const |
Return the value of epsilonSmall which is added to epsilon when.
| |
const dimensionedScalar & | omega0 () const |
Return the value of omega0 which epsilon is not allowed to be.
| |
const dimensionedScalar & | omegaSmall () const |
Return the value of omegaSmall which is added to epsilon when.
| |
dimensionedScalar & | k0 () |
Allow k0 to be changed.
| |
dimensionedScalar & | epsilon0 () |
Allow epsilon0 to be changed.
| |
dimensionedScalar & | epsilonSmall () |
Allow epsilonSmall to be changed.
| |
dimensionedScalar & | omega0 () |
Allow omega0 to be changed.
| |
dimensionedScalar & | omegaSmall () |
Allow omegaSmall to be changed.
| |
const nearWallDist & | y () const |
Return the near wall distances.
| |
scalar | yPlusLam (const scalar kappa, const scalar E) const |
Calculate y+ at the edge of the laminar sublayer.
| |
const dictionary & | coeffDict () const |
Const access to the coefficients dictionary.
| |
virtual tmp< volScalarField > | nut () const =0 |
Return the turbulence viscosity.
| |
virtual tmp< volScalarField > | nuEff () const |
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 tmp< scalarField > | yPlus (const label patchI, const scalar Cmu) const |
Return yPlus for the given patch.
| |
virtual void | correct ()=0 |
Solve the turbulence equations and correct the turbulence viscosity.
| |
virtual bool | read ()=0 |
Read RASProperties dictionary.
| |
Static Public Member Functions | |
static autoPtr< RASModel > | New (const volVectorField &U, const surfaceScalarField &phi, transportModel &lamTransportModel) |
Return a reference to the selected RAS model.
| |
Protected Member Functions | |
virtual void | printCoeffs () |
Print model coefficients.
| |
Protected Attributes | |
Switch | turbulence_ |
Turbulence on/off flag.
| |
Switch | printCoeffs_ |
Flag to print the model coeffs at run-time.
| |
dictionary | coeffDict_ |
Model coefficients dictionary.
| |
dimensionedScalar | k0_ |
Lower limit of k.
| |
dimensionedScalar | epsilon0_ |
Lower limit of epsilon.
| |
dimensionedScalar | epsilonSmall_ |
Small epsilon value used to avoid divide by zero.
| |
dimensionedScalar | omega0_ |
Lower limit for omega.
| |
dimensionedScalar | omegaSmall_ |
Small omega value used to avoid divide by zero.
| |
nearWallDist | y_ |
Near wall distance boundary field.
|
RASModel | ( | const word & | type, |
const volVectorField & | U, | ||
const surfaceScalarField & | phi, | ||
transportModel & | lamTransportModel | ||
) |
Construct from components.
Definition at line 57 of file RASModel.C.
References dictionary::readIfPresent().
virtual ~RASModel | ( | ) | [inline, virtual]
|
Destructor.
Definition at line 169 of file RASModel.H.
void printCoeffs | ( | ) | [protected, virtual]
|
Print model coefficients.
Definition at line 45 of file RASModel.C.
References Foam::endl(), Foam::Info, and Foam::type().
TypeName | ( | "RASModel" | ) |
Runtime type information.
declareRunTimeSelectionTable | ( | autoPtr | , |
RASModel | , | ||
dictionary | , | ||
(const volVectorField &U, const surfaceScalarField &phi, transportModel &lamTransportModel) | , | ||
(U, phi, lamTransportModel) | |||
) |
autoPtr< RASModel > New | ( | const volVectorField & | U, |
const surfaceScalarField & | phi, | ||
transportModel & | lamTransportModel | ||
) | [static]
|
Return a reference to the selected RAS model.
Definition at line 105 of file RASModel.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 dimensionedScalar& k0 | ( | ) | const [inline]
|
Return the value of k0 which k is not allowed to be less than.
Definition at line 178 of file RASModel.H.
References incompressible::RASModel::k0_.
const dimensionedScalar& epsilon0 | ( | ) | const [inline]
|
Return the value of epsilon0 which epsilon is not allowed to be.
less than
Definition at line 185 of file RASModel.H.
References incompressible::RASModel::epsilon0_.
const dimensionedScalar& epsilonSmall | ( | ) | const [inline]
|
Return the value of epsilonSmall which is added to epsilon when.
calculating nut
Definition at line 192 of file RASModel.H.
References incompressible::RASModel::epsilonSmall_.
const dimensionedScalar& omega0 | ( | ) | const [inline]
|
Return the value of omega0 which epsilon is not allowed to be.
less than
Definition at line 199 of file RASModel.H.
References incompressible::RASModel::omega0_.
const dimensionedScalar& omegaSmall | ( | ) | const [inline]
|
Return the value of omegaSmall which is added to epsilon when.
calculating nut
Definition at line 206 of file RASModel.H.
References incompressible::RASModel::omegaSmall_.
dimensionedScalar& k0 | ( | ) | [inline]
|
Allow k0 to be changed.
Definition at line 212 of file RASModel.H.
References incompressible::RASModel::k0_.
dimensionedScalar& epsilon0 | ( | ) | [inline]
|
Allow epsilon0 to be changed.
Definition at line 218 of file RASModel.H.
References incompressible::RASModel::epsilon0_.
dimensionedScalar& epsilonSmall | ( | ) | [inline]
|
Allow epsilonSmall to be changed.
Definition at line 224 of file RASModel.H.
References incompressible::RASModel::epsilonSmall_.
dimensionedScalar& omega0 | ( | ) | [inline]
|
Allow omega0 to be changed.
Definition at line 230 of file RASModel.H.
References incompressible::RASModel::omega0_.
dimensionedScalar& omegaSmall | ( | ) | [inline]
|
Allow omegaSmall to be changed.
Definition at line 236 of file RASModel.H.
References incompressible::RASModel::omegaSmall_.
const nearWallDist& y | ( | ) | const [inline]
|
Return the near wall distances.
Definition at line 242 of file RASModel.H.
References incompressible::RASModel::y_.
scalar yPlusLam | ( | const scalar | kappa, |
const scalar | E | ||
) | const |
Calculate y+ at the edge of the laminar sublayer.
Definition at line 157 of file RASModel.C.
References kappa(), Foam::log(), and Foam::max().
const dictionary& coeffDict | ( | ) | const [inline]
|
Const access to the coefficients dictionary.
Definition at line 251 of file RASModel.H.
References incompressible::RASModel::coeffDict_.
Referenced by qZeta::read(), NonlinearKEShih::read(), LienLeschzinerLowRe::read(), LienCubicKELowRe::read(), LienCubicKE::read(), LamBremhorstKE::read(), kOmega::read(), RNGkEpsilon::read(), realizableKE::read(), LRR::read(), LaunderSharmaKE::read(), LaunderGibsonRSTM::read(), kOmegaSST::read(), kEpsilon::read(), SpalartAllmaras::read(), incompressible::turbulentMixingLengthFrequencyInletFvPatchScalarField::updateCoeffs(), and incompressible::turbulentMixingLengthDissipationRateInletFvPatchScalarField::updateCoeffs().
virtual tmp<volScalarField> nut | ( | ) | const [pure virtual]
|
Return the turbulence viscosity.
Implements incompressible::turbulenceModel.
Implemented in kEpsilon, kOmega, kOmegaSST, LamBremhorstKE, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LienCubicKE, LienCubicKELowRe, LienLeschzinerLowRe, LRR, NonlinearKEShih, qZeta, realizableKE, RNGkEpsilon, and SpalartAllmaras.
Referenced by incompressible::RASModel::nuEff().
virtual tmp<volScalarField> nuEff | ( | ) | const [inline, virtual]
|
Return the effective viscosity.
Implements incompressible::turbulenceModel.
Reimplemented in laminar.
Definition at line 261 of file RASModel.H.
References nu, and incompressible::RASModel::nut().
Referenced by qZeta::devReff(), NonlinearKEShih::devReff(), LienLeschzinerLowRe::devReff(), LienCubicKELowRe::devReff(), LienCubicKE::devReff(), LamBremhorstKE::devReff(), kOmega::devReff(), RNGkEpsilon::devReff(), realizableKE::devReff(), LaunderSharmaKE::devReff(), kOmegaSST::devReff(), kEpsilon::devReff(), SpalartAllmaras::devReff(), qZeta::divDevReff(), NonlinearKEShih::divDevReff(), LienLeschzinerLowRe::divDevReff(), LienCubicKELowRe::divDevReff(), LienCubicKE::divDevReff(), LamBremhorstKE::divDevReff(), kOmega::divDevReff(), RNGkEpsilon::divDevReff(), realizableKE::divDevReff(), LRR::divDevReff(), LaunderSharmaKE::divDevReff(), LaunderGibsonRSTM::divDevReff(), kOmegaSST::divDevReff(), kEpsilon::divDevReff(), SpalartAllmaras::divDevReff(), and incompressible::fixedShearStressFvPatchVectorField::updateCoeffs().
virtual tmp<volScalarField> k | ( | ) | const [pure virtual]
|
Return the turbulence kinetic energy.
Implements incompressible::turbulenceModel.
Implemented in kEpsilon, kOmega, kOmegaSST, LamBremhorstKE, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LienCubicKE, LienCubicKELowRe, LienLeschzinerLowRe, LRR, NonlinearKEShih, qZeta, realizableKE, RNGkEpsilon, and SpalartAllmaras.
virtual tmp<volScalarField> epsilon | ( | ) | const [pure virtual]
|
Return the turbulence kinetic energy dissipation rate.
Implements incompressible::turbulenceModel.
Implemented in kEpsilon, kOmega, kOmegaSST, LamBremhorstKE, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LienCubicKE, LienCubicKELowRe, LienLeschzinerLowRe, LRR, NonlinearKEShih, qZeta, realizableKE, RNGkEpsilon, and SpalartAllmaras.
virtual tmp<volSymmTensorField> R | ( | ) | const [pure virtual]
|
Return the Reynolds stress tensor.
Implements incompressible::turbulenceModel.
Implemented in kEpsilon, kOmega, kOmegaSST, LamBremhorstKE, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LienCubicKE, LienCubicKELowRe, LienLeschzinerLowRe, LRR, NonlinearKEShih, qZeta, realizableKE, RNGkEpsilon, and SpalartAllmaras.
virtual tmp<volSymmTensorField> devReff | ( | ) | const [pure virtual]
|
Return the effective stress tensor including the laminar stress.
Implements incompressible::turbulenceModel.
Implemented in kEpsilon, kOmega, kOmegaSST, LamBremhorstKE, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LienCubicKE, LienCubicKELowRe, LienLeschzinerLowRe, LRR, NonlinearKEShih, qZeta, realizableKE, RNGkEpsilon, and SpalartAllmaras.
Referenced by incompressible::fixedShearStressFvPatchVectorField::updateCoeffs().
virtual tmp<fvVectorMatrix> divDevReff | ( | volVectorField & | U ) | const [pure virtual]
|
Return the source term for the momentum equation.
Implemented in kEpsilon, kOmega, kOmegaSST, LamBremhorstKE, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LienCubicKE, LienCubicKELowRe, LienLeschzinerLowRe, LRR, NonlinearKEShih, qZeta, realizableKE, RNGkEpsilon, and SpalartAllmaras.
tmp< scalarField > yPlus | ( | const label | patchI, |
const scalar | Cmu | ||
) | const [virtual]
|
Return yPlus for the given patch.
Definition at line 170 of file RASModel.C.
References boundaryField(), Foam::endl(), k(), Foam::nl, nu, Foam::pow(), List< T >::setSize(), fvPatch::size(), Foam::sqrt(), and WarningIn.
void correct | ( | ) | [pure virtual]
|
Solve the turbulence equations and correct the turbulence viscosity.
Implements incompressible::turbulenceModel.
Implemented in kEpsilon, kOmega, kOmegaSST, LamBremhorstKE, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LienCubicKE, LienCubicKELowRe, LienLeschzinerLowRe, LRR, NonlinearKEShih, qZeta, realizableKE, RNGkEpsilon, and SpalartAllmaras.
Definition at line 199 of file RASModel.C.
References correct().
bool read | ( | ) | [pure virtual]
|
Read RASProperties dictionary.
Implements incompressible::turbulenceModel.
Implemented in kEpsilon, kOmega, kOmegaSST, LamBremhorstKE, laminar, LaunderGibsonRSTM, LaunderSharmaKE, LienCubicKE, LienCubicKELowRe, LienLeschzinerLowRe, LRR, NonlinearKEShih, qZeta, realizableKE, RNGkEpsilon, and SpalartAllmaras.
Definition at line 210 of file RASModel.C.
References regIOobject::read(), and Foam::type().
Switch turbulence_ [protected]
|
Turbulence on/off flag.
Definition at line 79 of file RASModel.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(), and SpalartAllmaras::correct().
Switch printCoeffs_ [protected]
|
Flag to print the model coeffs at run-time.
Definition at line 82 of file RASModel.H.
dictionary coeffDict_ [protected]
|
Model coefficients dictionary.
Definition at line 85 of file RASModel.H.
Referenced by incompressible::RASModel::coeffDict().
dimensionedScalar k0_ [protected]
|
Lower limit of k.
Definition at line 88 of file RASModel.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(), and incompressible::RASModel::k0().
dimensionedScalar epsilon0_ [protected]
|
Lower limit of epsilon.
Definition at line 91 of file RASModel.H.
Referenced by qZeta::correct(), NonlinearKEShih::correct(), LienLeschzinerLowRe::correct(), LienCubicKELowRe::correct(), LienCubicKE::correct(), LamBremhorstKE::correct(), RNGkEpsilon::correct(), realizableKE::correct(), LRR::correct(), LaunderSharmaKE::correct(), LaunderGibsonRSTM::correct(), kEpsilon::correct(), and incompressible::RASModel::epsilon0().
dimensionedScalar epsilonSmall_ [protected]
|
Small epsilon value used to avoid divide by zero.
Definition at line 94 of file RASModel.H.
Referenced by incompressible::RASModel::epsilonSmall().
dimensionedScalar omega0_ [protected]
|
Lower limit for omega.
Definition at line 97 of file RASModel.H.
Referenced by kOmega::correct(), kOmegaSST::correct(), and incompressible::RASModel::omega0().
dimensionedScalar omegaSmall_ [protected]
|
Small omega value used to avoid divide by zero.
Definition at line 100 of file RASModel.H.
Referenced by kOmega::correct(), and incompressible::RASModel::omegaSmall().
nearWallDist y_ [protected]
|
Near wall distance boundary field.
Definition at line 103 of file RASModel.H.
Referenced by incompressible::RASModel::y().