Base class for all incompressible flow LES SGS models. More...
#include <incompressibleLESModels/LESModel.H>
Base class for all incompressible flow LES SGS models.
This class defines the basic interface for an incompressible flow SGS model, and encapsulates data of value to all possible models. In particular this includes references to all the dependent fields (U, phi), the physical viscosity nu, and the LESProperties dictionary, which contains the model selection and model coefficients.
Definition at line 71 of file LESModel.H.
Public Member Functions | |
TypeName ("LESModel") | |
Runtime type information.
| |
declareRunTimeSelectionTable (autoPtr, LESModel, dictionary,(const volVectorField &U, const surfaceScalarField &phi, transportModel &lamTransportModel),(U, phi, lamTransportModel)) | |
LESModel (const word &type, const volVectorField &U, const surfaceScalarField &phi, transportModel &lamTransportModel) | |
Construct from components.
| |
virtual | ~LESModel () |
Destructor.
| |
const dictionary & | coeffDict () const |
Const access to the coefficients dictionary,.
| |
virtual const volScalarField & | delta () const |
Access function to filter width.
| |
const dimensionedScalar & | k0 () const |
Return the value of k0 which k is not allowed to be less than.
| |
dimensionedScalar & | k0 () |
Allow k0 to be changed.
| |
virtual tmp< volScalarField > | k () const =0 |
Return the SGS turbulent kinetic energy.
| |
virtual tmp< volScalarField > | epsilon () const =0 |
Return the SGS turbulent dissipation.
| |
virtual tmp< volScalarField > | nuSgs () const =0 |
Return the SGS viscosity.
| |
virtual tmp< volScalarField > | nuEff () const |
Return the effective viscosity.
| |
virtual tmp< volSymmTensorField > | B () const =0 |
Return the sub-grid stress tensor.
| |
virtual tmp< volSymmTensorField > | devBeff () const =0 |
Return the deviatoric part of the effective sub-grid.
| |
virtual tmp< fvVectorMatrix > | divDevBeff (volVectorField &U) const =0 |
Returns div(dev(Beff)).
| |
virtual tmp< volScalarField > | nut () const |
Return the turbulence viscosity.
| |
virtual tmp< volSymmTensorField > | R () const |
Return the Reynolds stress tensor.
| |
virtual tmp< volSymmTensorField > | devReff () const |
Return the effective stress tensor including the laminar stress.
| |
virtual tmp< fvVectorMatrix > | divDevReff (volVectorField &U) const |
Return the source term for the momentum equation.
| |
void | correct () |
Correct Eddy-Viscosity and related properties.
| |
virtual void | correct (const tmp< volTensorField > &gradU) |
Correct Eddy-Viscosity and related properties.
| |
virtual bool | read ()=0 |
Read LESProperties dictionary.
| |
Static Public Member Functions | |
static autoPtr< LESModel > | New (const volVectorField &U, const surfaceScalarField &phi, transportModel &lamTransportModel) |
Return a reference to the selected LES model.
| |
Protected Member Functions | |
virtual void | printCoeffs () |
Print model coefficients.
| |
Protected Attributes | |
Switch | printCoeffs_ |
dictionary | coeffDict_ |
dimensionedScalar | k0_ |
autoPtr< LESdelta > | delta_ |
LESModel | ( | const word & | type, |
const volVectorField & | U, | ||
const surfaceScalarField & | phi, | ||
transportModel & | lamTransportModel | ||
) |
virtual ~LESModel | ( | ) | [inline, virtual]
|
Destructor.
Definition at line 152 of file LESModel.H.
void printCoeffs | ( | ) | [protected, virtual]
|
Print model coefficients.
Definition at line 44 of file LESModel.C.
References Foam::endl(), Foam::Info, and Foam::type().
TypeName | ( | "LESModel" | ) |
Runtime type information.
declareRunTimeSelectionTable | ( | autoPtr | , |
LESModel | , | ||
dictionary | , | ||
(const volVectorField &U, const surfaceScalarField &phi, transportModel &lamTransportModel) | , | ||
(U, phi, lamTransportModel) | |||
) |
autoPtr< LESModel > New | ( | const volVectorField & | U, |
const surfaceScalarField & | phi, | ||
transportModel & | lamTransportModel | ||
) | [static]
|
Return a reference to the selected LES model.
Definition at line 94 of file LESModel.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 dictionary& coeffDict | ( | ) | const [inline]
|
Const access to the coefficients dictionary,.
which provides info. about choice of models, and all related data (particularly model coefficients).
Definition at line 161 of file LESModel.H.
References incompressible::LESModel::coeffDict_.
Referenced by spectEddyVisc::read(), SpalartAllmarasIDDES::read(), Smagorinsky2::read(), scaleSimilarity::read(), LRRDiffStress::read(), locDynOneEqEddy::read(), kOmegaSSTSAS::read(), homogeneousDynSmagorinsky::read(), SpalartAllmaras::read(), Smagorinsky::read(), oneEqEddy::read(), GenSGSStress::read(), GenEddyVisc::read(), dynOneEqEddy::read(), and DeardorffDiffStress::read().
virtual const volScalarField& delta | ( | ) | const [inline, virtual]
|
Access function to filter width.
Reimplemented in SpalartAllmarasIDDES.
Definition at line 167 of file LESModel.H.
References incompressible::LESModel::delta_.
Referenced by Smagorinsky2::B(), LRRDiffStress::correct(), locDynOneEqEddy::correct(), homogeneousDynSmagorinsky::correct(), oneEqEddy::correct(), dynOneEqEddy::correct(), DeardorffDiffStress::correct(), SpalartAllmarasDDES::dTilda(), SpalartAllmaras::dTilda(), GenSGSStress::epsilon(), GenEddyVisc::epsilon(), spectEddyVisc::k(), and Smagorinsky::k().
const dimensionedScalar& k0 | ( | ) | const [inline]
|
Return the value of k0 which k is not allowed to be less than.
Definition at line 173 of file LESModel.H.
References incompressible::LESModel::k0_.
Referenced by LRRDiffStress::correct(), locDynOneEqEddy::correct(), kOmegaSSTSAS::correct(), oneEqEddy::correct(), dynOneEqEddy::correct(), and DeardorffDiffStress::correct().
dimensionedScalar& k0 | ( | ) | [inline]
|
Allow k0 to be changed.
Definition at line 179 of file LESModel.H.
References incompressible::LESModel::k0_.
virtual tmp<volScalarField> k | ( | ) | const [pure virtual]
|
Return the SGS turbulent kinetic energy.
Implements incompressible::turbulenceModel.
Implemented in dynOneEqEddy, GenEddyVisc, GenSGSStress, homogeneousDynSmagorinsky, kOmegaSSTSAS, laminar, locDynOneEqEddy, mixedSmagorinsky, oneEqEddy, scaleSimilarity, Smagorinsky, SpalartAllmaras, and spectEddyVisc.
virtual tmp<volScalarField> epsilon | ( | ) | const [pure virtual]
|
Return the SGS turbulent dissipation.
Implements incompressible::turbulenceModel.
Implemented in GenEddyVisc, GenSGSStress, kOmegaSSTSAS, laminar, mixedSmagorinsky, scaleSimilarity, and SpalartAllmaras.
virtual tmp<volScalarField> nuSgs | ( | ) | const [pure virtual]
|
Return the SGS viscosity.
Implemented in GenEddyVisc, GenSGSStress, kOmegaSSTSAS, laminar, mixedSmagorinsky, and SpalartAllmaras.
Referenced by incompressible::LESModel::nuEff(), and incompressible::LESModel::nut().
virtual tmp<volScalarField> nuEff | ( | ) | const [inline, virtual]
|
Return the effective viscosity.
Implements incompressible::turbulenceModel.
Reimplemented in laminar.
Definition at line 195 of file LESModel.H.
References nu, and incompressible::LESModel::nuSgs().
Referenced by LRRDiffStress::correct(), DeardorffDiffStress::correct(), kOmegaSSTSAS::devBeff(), SpalartAllmaras::devBeff(), GenEddyVisc::devBeff(), kOmegaSSTSAS::divDevBeff(), SpalartAllmaras::divDevBeff(), GenEddyVisc::divDevBeff(), kOmegaSSTSAS::epsilon(), SpalartAllmaras::epsilon(), and spectEddyVisc::k().
virtual tmp<volSymmTensorField> B | ( | ) | const [pure virtual]
|
Return the sub-grid stress tensor.
Implemented in GenEddyVisc, GenSGSStress, kOmegaSSTSAS, laminar, mixedSmagorinsky, scaleSimilarity, Smagorinsky2, and SpalartAllmaras.
Referenced by incompressible::LESModel::R().
virtual tmp<volSymmTensorField> devBeff | ( | ) | const [pure virtual]
|
Return the deviatoric part of the effective sub-grid.
turbulence stress tensor including the laminar stress
Implemented in GenEddyVisc, GenSGSStress, kOmegaSSTSAS, laminar, mixedSmagorinsky, scaleSimilarity, and SpalartAllmaras.
Referenced by incompressible::LESModel::devReff().
virtual tmp<fvVectorMatrix> divDevBeff | ( | volVectorField & | U ) | const [pure virtual]
|
Returns div(dev(Beff)).
This is the additional term due to the filtering of the NSE.
Implemented in GenEddyVisc, GenSGSStress, kOmegaSSTSAS, laminar, mixedSmagorinsky, scaleSimilarity, Smagorinsky2, and SpalartAllmaras.
Referenced by incompressible::LESModel::divDevReff().
virtual tmp<volScalarField> nut | ( | ) | const [inline, virtual]
|
Return the turbulence viscosity.
Implements incompressible::turbulenceModel.
Definition at line 218 of file LESModel.H.
References incompressible::LESModel::nuSgs().
virtual tmp<volSymmTensorField> R | ( | ) | const [inline, virtual]
|
Return the Reynolds stress tensor.
Implements incompressible::turbulenceModel.
Definition at line 224 of file LESModel.H.
References incompressible::LESModel::B().
virtual tmp<volSymmTensorField> devReff | ( | ) | const [inline, virtual]
|
Return the effective stress tensor including the laminar stress.
Implements incompressible::turbulenceModel.
Definition at line 230 of file LESModel.H.
References incompressible::LESModel::devBeff().
virtual tmp<fvVectorMatrix> divDevReff | ( | volVectorField & | U ) | const [inline, virtual]
|
Return the source term for the momentum equation.
Definition at line 236 of file LESModel.H.
References incompressible::LESModel::divDevBeff().
void correct | ( | ) | [virtual]
|
Correct Eddy-Viscosity and related properties.
This calls correct(const tmp<volTensorField>& gradU) by supplying gradU calculated locally.
Implements incompressible::turbulenceModel.
Definition at line 152 of file LESModel.C.
References correct(), and Foam::fvc::grad().
Referenced by spectEddyVisc::correct(), mixedSmagorinsky::correct(), LRRDiffStress::correct(), locDynOneEqEddy::correct(), kOmegaSSTSAS::correct(), homogeneousDynSmagorinsky::correct(), SpalartAllmaras::correct(), Smagorinsky::correct(), oneEqEddy::correct(), GenEddyVisc::correct(), dynOneEqEddy::correct(), and DeardorffDiffStress::correct().
void correct | ( | const tmp< volTensorField > & | gradU ) | [virtual]
|
Correct Eddy-Viscosity and related properties.
Reimplemented in DeardorffDiffStress, dynOneEqEddy, GenEddyVisc, homogeneousDynSmagorinsky, kOmegaSSTSAS, locDynOneEqEddy, LRRDiffStress, mixedSmagorinsky, oneEqEddy, scaleSimilarity, Smagorinsky, SpalartAllmaras, and spectEddyVisc.
Definition at line 145 of file LESModel.C.
References correct().
bool read | ( | ) | [pure virtual]
|
Read LESProperties dictionary.
Implements incompressible::turbulenceModel.
Implemented in DeardorffDiffStress, dynOneEqEddy, GenEddyVisc, GenSGSStress, homogeneousDynSmagorinsky, kOmegaSSTSAS, laminar, locDynOneEqEddy, LRRDiffStress, mixedSmagorinsky, oneEqEddy, scaleSimilarity, Smagorinsky, Smagorinsky2, SpalartAllmaras, SpalartAllmarasIDDES, and spectEddyVisc.
Definition at line 158 of file LESModel.C.
References regIOobject::read(), readIfPresent(), and Foam::type().
Switch printCoeffs_ [protected]
|
Definition at line 81 of file LESModel.H.
dictionary coeffDict_ [protected]
|
Definition at line 82 of file LESModel.H.
Referenced by incompressible::LESModel::coeffDict().
dimensionedScalar k0_ [protected]
|
Definition at line 84 of file LESModel.H.
Referenced by incompressible::LESModel::k0().
Definition at line 86 of file LESModel.H.
Referenced by incompressible::LESModel::delta().