Templated abstract base class for LES SGS models. More...
Public Types | |
typedef BasicMomentumTransportModel::alphaField | alphaField |
typedef BasicMomentumTransportModel::rhoField | rhoField |
typedef BasicMomentumTransportModel::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),(alpha, rho, U, alphaRhoPhi, phi, transport)) | |
LESModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport) | |
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 > | 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... | |
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 transportModel &transport) |
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< Foam::LESdelta > | delta_ |
Run-time selectable delta model. More... | |
Templated abstract base class for LES SGS models.
Definition at line 56 of file LESModel.H.
typedef BasicMomentumTransportModel::alphaField alphaField |
Definition at line 98 of file LESModel.H.
typedef BasicMomentumTransportModel::rhoField rhoField |
Definition at line 99 of file LESModel.H.
typedef BasicMomentumTransportModel::transportModel transportModel |
Definition at line 100 of file LESModel.H.
LESModel | ( | const word & | type, |
const alphaField & | alpha, | ||
const rhoField & | rho, | ||
const volVectorField & | U, | ||
const surfaceScalarField & | alphaRhoPhi, | ||
const surfaceScalarField & | phi, | ||
const transportModel & | transport | ||
) |
Construct from components.
Definition at line 44 of file LESModel.C.
References LESModel< BasicMomentumTransportModel >::New().
Referenced by LESModel< BasicMomentumTransportModel >::printCoeffs().
Disallow default bitwise copy construction.
|
inlinevirtual |
Destructor.
Definition at line 159 of file LESModel.H.
References LESModel< BasicMomentumTransportModel >::read().
|
protectedvirtual |
Print model coefficients.
Definition at line 31 of file LESModel.C.
References Foam::endl(), Foam::Info, and LESModel< BasicMomentumTransportModel >::LESModel().
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 transportModel &transport) | , | ||
(alpha, rho, U, alphaRhoPhi, phi, transport) | |||
) |
|
static |
Return a reference to the selected LES model.
Definition at line 124 of file LESModel.C.
References alpha(), IOobject::db(), Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, IOobject::group(), Foam::Info, dictionary::lookupBackwardsCompatible(), Foam::nl, phi, rho, dictionary::subDict(), and U.
Referenced by LESModel< BasicMomentumTransportModel >::LESModel().
|
virtual |
Read model coefficients if they have changed.
Reimplemented in SpalartAllmarasDES< BasicMomentumTransportModel >, dynamicKEqn< BasicMomentumTransportModel >, SpalartAllmarasIDDES< BasicMomentumTransportModel >, NicenoKEqn< BasicMomentumTransportModel >, Smagorinsky< BasicMomentumTransportModel >, SmagorinskyZhang< BasicMomentumTransportModel >, DeardorffDiffStress< BasicMomentumTransportModel >, WALE< BasicMomentumTransportModel >, continuousGasKEqn< BasicMomentumTransportModel >, kEqn< BasicMomentumTransportModel >, dynamicLagrangian< BasicMomentumTransportModel >, ReynoldsStress< LESModel< BasicMomentumTransportModel > >, eddyViscosity< LESModel< BasicMomentumTransportModel > >, LESeddyViscosity< BasicMomentumTransportModel >, and linearViscousStress< LESModel< BasicMomentumTransportModel > >.
Definition at line 173 of file LESModel.C.
References Foam::blockMeshTools::read(), and Foam::type().
Referenced by LESModel< BasicMomentumTransportModel >::~LESModel().
|
inlinevirtual |
Const access to the coefficients dictionary.
Definition at line 172 of file LESModel.H.
References LESModel< BasicMomentumTransportModel >::coeffDict_.
|
inline |
Return the lower allowable limit for k (default: small)
Definition at line 178 of file LESModel.H.
References LESModel< BasicMomentumTransportModel >::kMin_.
|
inline |
Allow kMin to be changed.
Definition at line 184 of file LESModel.H.
References LESModel< BasicMomentumTransportModel >::kMin_.
|
inline |
Access function to filter width.
Definition at line 190 of file LESModel.H.
References LESModel< BasicMomentumTransportModel >::delta_.
Referenced by dynamicLagrangian< BasicMomentumTransportModel >::k().
|
inlinevirtual |
Return the effective viscosity.
Definition at line 197 of file LESModel.H.
References IOobject::groupName(), GeometricField< scalar, fvPatchField, volMesh >::New(), and nut.
|
inlinevirtual |
Return the effective viscosity on patch.
Definition at line 207 of file LESModel.H.
References LESModel< BasicMomentumTransportModel >::correct(), nut, and LESModel< BasicMomentumTransportModel >::operator=().
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Reimplemented in SpalartAllmarasDES< BasicMomentumTransportModel >, dynamicKEqn< BasicMomentumTransportModel >, Smagorinsky< BasicMomentumTransportModel >, kEqn< BasicMomentumTransportModel >, ReynoldsStress< LESModel< BasicMomentumTransportModel > >, dynamicLagrangian< BasicMomentumTransportModel >, WALE< BasicMomentumTransportModel >, DeardorffDiffStress< BasicMomentumTransportModel >, eddyViscosity< LESModel< BasicMomentumTransportModel > >, and linearViscousStress< LESModel< BasicMomentumTransportModel > >.
Definition at line 196 of file LESModel.C.
References correct.
Referenced by LESModel< BasicMomentumTransportModel >::nuEff().
|
delete |
Disallow default bitwise assignment.
Referenced by LESModel< BasicMomentumTransportModel >::nuEff().
|
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< BasicMomentumTransportModel >::coeffDict().
|
protected |
Lower limit of k.
Definition at line 78 of file LESModel.H.
Referenced by LESModel< BasicMomentumTransportModel >::kMin().
|
protected |
Lower limit of epsilon.
Definition at line 81 of file LESModel.H.
|
protected |
Lower limit for omega.
Definition at line 84 of file LESModel.H.
|
protected |
Run-time selectable delta model.
Definition at line 87 of file LESModel.H.
Referenced by LESModel< BasicMomentumTransportModel >::delta().