Templated abstract base class for LES SGS models. More...
Public Types | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
Public Member Functions | |
TypeName ("LES") | |
Runtime type information. More... | |
declareRunTimeSelectionTable (autoPtr, LESModel, 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)) | |
LESModel (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 | ~LESModel () |
Destructor. More... | |
virtual bool | read () |
Read model coefficients if they have changed. More... | |
virtual const dictionary & | coeffDict () const |
Const access to the coefficients dictionary. More... | |
const dimensionedScalar & | kMin () const |
Return the lower allowable limit for k (default: SMALL) More... | |
dimensionedScalar & | kMin () |
Allow kMin to be changed. More... | |
const volScalarField & | delta () const |
Access function to filter width. 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... | |
virtual void | correct () |
Solve the turbulence equations and correct the turbulence viscosity. More... | |
Static Public Member Functions | |
static autoPtr< LESModel > | 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 LES model. More... | |
Protected Member Functions | |
virtual void | printCoeffs (const word &type) |
Print model coefficients. More... | |
Protected Attributes | |
dictionary | LESDict_ |
LES 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 | omegaMin_ |
Lower limit of omega. More... | |
autoPtr< Foam::LESdelta > | delta_ |
Run-time selectable delta model. More... | |
Private Member Functions | |
LESModel (const LESModel &) | |
Disallow default bitwise copy construct. More... | |
void | operator= (const LESModel &) |
Disallow default bitwise assignment. More... | |
Templated abstract base class for LES SGS models.
Definition at line 56 of file LESModel.H.
typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 106 of file LESModel.H.
typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 107 of file LESModel.H.
typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 108 of file LESModel.H.
Disallow default bitwise copy construct.
LESModel | ( | 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.
Definition at line 44 of file LESModel.C.
|
inlinevirtual |
Destructor.
Definition at line 167 of file LESModel.H.
|
protectedvirtual |
Print model coefficients.
Definition at line 31 of file LESModel.C.
References Foam::endl(), Foam::Info, and Foam::type().
|
private |
Disallow default bitwise assignment.
TypeName | ( | "LES" | ) |
Runtime type information.
declareRunTimeSelectionTable | ( | autoPtr | , |
LESModel< BasicTurbulenceModel > | , | ||
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) | |||
) |
|
static |
Return a reference to the selected LES model.
Definition at line 115 of file LESModel.C.
References Foam::constant::atomic::alpha, Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, Foam::Info, dictionary::lookup(), Foam::nl, phi, rho, dictionary::subDict(), and U.
|
virtual |
Read model coefficients if they have changed.
Reimplemented in ReynoldsStress< LESModel< BasicTurbulenceModel > >, eddyViscosity< LESModel< BasicTurbulenceModel > >, linearViscousStress< LESModel< BasicTurbulenceModel > >, SpalartAllmarasDES< BasicTurbulenceModel >, dynamicKEqn< BasicTurbulenceModel >, NicenoKEqn< BasicTurbulenceModel >, SpalartAllmarasIDDES< BasicTurbulenceModel >, Smagorinsky< BasicTurbulenceModel >, SmagorinskyZhang< BasicTurbulenceModel >, DeardorffDiffStress< BasicTurbulenceModel >, WALE< BasicTurbulenceModel >, continuousGasKEqn< BasicTurbulenceModel >, kEqn< BasicTurbulenceModel >, SpalartAllmarasDDES< BasicTurbulenceModel >, dynamicLagrangian< BasicTurbulenceModel >, and LESeddyViscosity< BasicTurbulenceModel >.
Definition at line 168 of file LESModel.C.
References Foam::read(), dictionary::subDictPtr(), and Foam::type().
|
inlinevirtual |
Const access to the coefficients dictionary.
Definition at line 180 of file LESModel.H.
References LESModel< BasicTurbulenceModel >::coeffDict_.
|
inline |
Return the lower allowable limit for k (default: SMALL)
Definition at line 186 of file LESModel.H.
References LESModel< BasicTurbulenceModel >::kMin_.
|
inline |
Allow kMin to be changed.
Definition at line 192 of file LESModel.H.
References LESModel< BasicTurbulenceModel >::kMin_.
|
inline |
Access function to filter width.
Definition at line 198 of file LESModel.H.
References LESModel< BasicTurbulenceModel >::delta_.
Referenced by dynamicLagrangian< BasicTurbulenceModel >::k().
|
inlinevirtual |
Return the effective viscosity.
Definition at line 205 of file LESModel.H.
References IOobject::groupName().
|
inlinevirtual |
Return the effective viscosity on patch.
Definition at line 218 of file LESModel.H.
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Reimplemented in ReynoldsStress< LESModel< BasicTurbulenceModel > >, eddyViscosity< LESModel< BasicTurbulenceModel > >, linearViscousStress< LESModel< BasicTurbulenceModel > >, SpalartAllmarasDES< BasicTurbulenceModel >, dynamicKEqn< BasicTurbulenceModel >, Smagorinsky< BasicTurbulenceModel >, kEqn< BasicTurbulenceModel >, dynamicLagrangian< BasicTurbulenceModel >, WALE< BasicTurbulenceModel >, and DeardorffDiffStress< BasicTurbulenceModel >.
Definition at line 194 of file LESModel.C.
References correct().
|
protected |
LES coefficients dictionary.
Definition at line 66 of file LESModel.H.
|
protected |
Turbulence on/off flag.
Definition at line 69 of file LESModel.H.
|
protected |
Flag to print the model coeffs at run-time.
Definition at line 72 of file LESModel.H.
|
protected |
Model coefficients dictionary.
Definition at line 75 of file LESModel.H.
Referenced by LESModel< BasicTurbulenceModel >::coeffDict().
|
protected |
Lower limit of k.
Definition at line 78 of file LESModel.H.
Referenced by LESModel< BasicTurbulenceModel >::kMin().
|
protected |
Lower limit of omega.
Definition at line 81 of file LESModel.H.
|
protected |
Run-time selectable delta model.
Definition at line 84 of file LESModel.H.
Referenced by LESModel< BasicTurbulenceModel >::delta().
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.