Templated abstract base class for LES SGS models. More...
Public Types | |
typedef BasicMomentumTransportModel::alphaField | alphaField |
typedef BasicMomentumTransportModel::rhoField | rhoField |
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 viscosity &viscosity),(alpha, rho, U, alphaRhoPhi, phi, viscosity)) | |
LESModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity) | |
Construct from components. More... | |
LESModel (const LESModel &)=delete | |
Disallow default bitwise copy construction. 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 > | nu () const |
Return the laminar viscosity. More... | |
virtual tmp< scalarField > | nu (const label patchi) const |
Return the laminar viscosity on patchi. 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 | predict () |
Predict the turbulence transport coefficients if possible. More... | |
virtual void | correct () |
Solve the turbulence equations and correct the turbulence viscosity. More... | |
void | operator= (const LESModel &)=delete |
Disallow default bitwise assignment. 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 viscosity &viscosity) |
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 | epsilonMin_ |
Lower limit of epsilon. More... | |
dimensionedScalar | omegaMin_ |
Lower limit for omega. More... | |
autoPtr< laminarModels::generalisedNewtonianViscosityModel > | viscosityModel_ |
Run-time selectable generalised Newtonian viscosity model. More... | |
autoPtr< Foam::LESdelta > | delta_ |
Run-time selectable delta model. More... | |
Templated abstract base class for LES SGS models.
with support for generalised Newtonian viscosity models including strain-rate dependency.
Definition at line 60 of file LESModel.H.
typedef BasicMomentumTransportModel::alphaField alphaField |
Definition at line 106 of file LESModel.H.
typedef BasicMomentumTransportModel::rhoField rhoField |
Definition at line 107 of file LESModel.H.
LESModel | ( | const word & | type, |
const alphaField & | alpha, | ||
const rhoField & | rho, | ||
const volVectorField & | U, | ||
const surfaceScalarField & | alphaRhoPhi, | ||
const surfaceScalarField & | phi, | ||
const viscosity & | viscosity | ||
) |
Construct from components.
Definition at line 44 of file LESModel.C.
Disallow default bitwise copy construction.
|
inlinevirtual |
Destructor.
Definition at line 166 of file LESModel.H.
|
protectedvirtual |
Print model coefficients.
Definition at line 32 of file LESModel.C.
References Foam::endl(), and Foam::Info.
Referenced by DeardorffDiffStress< BasicMomentumTransportModel >::DeardorffDiffStress().
TypeName | ( | "LES" | ) |
Runtime type information.
declareRunTimeSelectionTable | ( | autoPtr | , |
LESModel< BasicMomentumTransportModel > | , | ||
dictionary | , | ||
(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity) | , | ||
(alpha, rho, U, alphaRhoPhi, phi, viscosity) | |||
) |
|
static |
Return a reference to the selected LES model.
Definition at line 144 of file LESModel.C.
References alpha(), Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, IOobject::group(), Foam::Info, dictionary::lookupBackwardsCompatible(), Foam::nl, momentumTransportModel::readModelDict(), rho, dictionary::subDict(), and U.
|
virtual |
Read model coefficients if they have changed.
Reimplemented in ReynoldsStress< LESModel< BasicMomentumTransportModel > >, and DeardorffDiffStress< BasicMomentumTransportModel >.
Definition at line 194 of file LESModel.C.
References Foam::blockMeshTools::read(), and Foam::type().
|
inlinevirtual |
Const access to the coefficients dictionary.
Definition at line 176 of file LESModel.H.
References LESModel< BasicMomentumTransportModel >::coeffDict_.
|
inline |
Return the lower allowable limit for k (default: small)
Definition at line 182 of file LESModel.H.
References LESModel< BasicMomentumTransportModel >::kMin_.
|
inline |
Allow kMin to be changed.
Definition at line 188 of file LESModel.H.
References LESModel< BasicMomentumTransportModel >::kMin_.
|
inline |
Access function to filter width.
Definition at line 194 of file LESModel.H.
References LESModel< BasicMomentumTransportModel >::delta_.
|
inlinevirtual |
Return the laminar viscosity.
Definition at line 200 of file LESModel.H.
References LESModel< BasicMomentumTransportModel >::viscosityModel_.
Referenced by LESModel< BasicMomentumTransportModel >::nuEff().
|
inlinevirtual |
Return the laminar viscosity on patchi.
Definition at line 206 of file LESModel.H.
References patchi, and LESModel< BasicMomentumTransportModel >::viscosityModel_.
|
inlinevirtual |
Return the effective viscosity.
Definition at line 212 of file LESModel.H.
References GeometricField< Type, PatchField, GeoMesh >::New(), LESModel< BasicMomentumTransportModel >::nu(), and nut.
|
inlinevirtual |
Return the effective viscosity on patch.
Definition at line 222 of file LESModel.H.
References LESModel< BasicMomentumTransportModel >::nu(), and nut.
|
inlinevirtual |
Predict the turbulence transport coefficients if possible.
without solving turbulence transport model equations
Definition at line 229 of file LESModel.H.
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Reimplemented in ReynoldsStress< LESModel< BasicMomentumTransportModel > >, and DeardorffDiffStress< BasicMomentumTransportModel >.
Definition at line 217 of file LESModel.C.
References Foam::MULES::correct().
|
delete |
Disallow default bitwise assignment.
|
protected |
LES coefficients dictionary.
Definition at line 70 of file LESModel.H.
|
protected |
Turbulence on/off flag.
Definition at line 73 of file LESModel.H.
|
protected |
Flag to print the model coeffs at run-time.
Definition at line 76 of file LESModel.H.
|
protected |
Model coefficients dictionary.
Definition at line 79 of file LESModel.H.
Referenced by LESModel< BasicMomentumTransportModel >::coeffDict().
|
protected |
Lower limit of k.
Definition at line 82 of file LESModel.H.
Referenced by LESModel< BasicMomentumTransportModel >::kMin().
|
protected |
Lower limit of epsilon.
Definition at line 85 of file LESModel.H.
|
protected |
Lower limit for omega.
Definition at line 88 of file LESModel.H.
|
protected |
Run-time selectable generalised Newtonian viscosity model.
Definition at line 92 of file LESModel.H.
Referenced by LESModel< BasicMomentumTransportModel >::nu().
|
protected |
Run-time selectable delta model.
Definition at line 95 of file LESModel.H.
Referenced by LESModel< BasicMomentumTransportModel >::delta().