Templated abstract base class for LES SGS models. More...
#include <LESModel.H>


Public Types | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
Public Member Functions | |
| TypeName ("LES") | |
| Runtime type information. | |
| 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. | |
| virtual | ~LESModel ()=default |
| Destructor. | |
| virtual bool | read () |
| Read model coefficients if they have changed. | |
| virtual const dictionary & | coeffDict () const |
| Const access to the coefficients dictionary. | |
| const dimensionedScalar & | Ce () const noexcept |
| const dimensionedScalar & | kMin () const |
| Return the lower allowable limit for k (default: SMALL). | |
| dimensionedScalar & | kMin () |
| Allow kMin to be changed. | |
| const volScalarField & | delta () const |
| Access function to filter width. | |
| virtual tmp< volScalarField > | nuEff () const |
| Return the effective viscosity. | |
| virtual tmp< scalarField > | nuEff (const label patchi) const |
| Return the effective viscosity on patch. | |
| virtual tmp< volScalarField > | epsilon () const |
| Return the turbulence kinetic energy dissipation rate. | |
| virtual tmp< volScalarField > | omega () const |
| Return the specific dissipation rate. | |
| virtual void | correct () |
| Solve the turbulence equations and correct the turbulence viscosity. | |
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. | |
Protected Member Functions | |
| virtual void | printCoeffs (const word &type) |
| Print model coefficients. | |
| LESModel (const LESModel &)=delete | |
| No copy construct. | |
| void | operator= (const LESModel &)=delete |
| No copy assignment. | |
Protected Attributes | |
| dictionary | LESDict_ |
| LES coefficients dictionary. | |
| Switch | turbulence_ |
| Turbulence on/off flag. | |
| Switch | printCoeffs_ |
| Flag to print the model coeffs at run-time. | |
| dictionary | coeffDict_ |
| Model coefficients dictionary. | |
| dimensionedScalar | Ce_ |
| Empirical model constant. | |
| dimensionedScalar | kMin_ |
| Lower limit of k. | |
| dimensionedScalar | epsilonMin_ |
| Lower limit of epsilon. | |
| dimensionedScalar | omegaMin_ |
| Lower limit for omega. | |
| autoPtr< Foam::LESdelta > | delta_ |
| Run-time selectable delta model. | |
Templated abstract base class for LES SGS models.
Definition at line 56 of file LESModel.H.
| typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 130 of file LESModel.H.
| typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 131 of file LESModel.H.
| typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 132 of file LESModel.H.
|
protecteddelete |
No 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 39 of file LESModel.C.
References alpha, Ce_, coeffDict_, delta_, Foam::dimless, Foam::dimTime, Foam::dimVelocity, epsilonMin_, kMin_, LESDict_, New(), omegaMin_, phi, printCoeffs_, rho, Foam::sqr(), turbulence_, and U.

|
virtualdefault |
Destructor.
|
protectedvirtual |
Print model coefficients.
Definition at line 27 of file LESModel.C.
References coeffDict_, Foam::endl(), Foam::Info, and printCoeffs_.

|
protecteddelete |
No copy 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 125 of file LESModel.C.
References alpha, dict, IOobject::group(), IOobject::groupName(), IOobjectOption::MUST_READ, IOobjectOption::NO_REGISTER, IOobjectOption::NO_WRITE, phi, rho, dictionary::subDict(), and U.
Referenced by LESModel().


|
virtual |
Read model coefficients if they have changed.
Definition at line 182 of file LESModel.C.
References Ce_, coeffDict_, delta_, kMin_, LESDict_, turbulence_, and Foam::type().

|
inlinevirtual |
Const access to the coefficients dictionary.
Definition at line 215 of file LESModel.H.
|
inlinenoexcept |
Definition at line 221 of file LESModel.H.
|
inline |
Return the lower allowable limit for k (default: SMALL).
Definition at line 229 of file LESModel.H.
|
inline |
Allow kMin to be changed.
Definition at line 237 of file LESModel.H.
|
inline |
Access function to filter width.
Definition at line 245 of file LESModel.H.
|
inlinevirtual |
Return the effective viscosity.
Definition at line 254 of file LESModel.H.
|
inlinevirtual |
Return the effective viscosity on patch.
Definition at line 269 of file LESModel.H.
|
virtual |
Return the turbulence kinetic energy dissipation rate.
Definition at line 206 of file LESModel.C.
References delta, IOobject::groupName(), k, tmp< T >::New(), and Foam::pow().

|
virtual |
Return the specific dissipation rate.
Definition at line 223 of file LESModel.C.
References Foam::dimLength, Foam::dimTime, IOobject::groupName(), tmp< T >::New(), and Foam::sqr().

|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Definition at line 242 of file LESModel.C.
References delta_.
|
protected |
LES coefficients dictionary.
Definition at line 67 of file LESModel.H.
Referenced by LESModel(), and read().
|
protected |
Turbulence on/off flag.
Definition at line 72 of file LESModel.H.
Referenced by LESModel(), and read().
|
protected |
Flag to print the model coeffs at run-time.
Definition at line 77 of file LESModel.H.
Referenced by LESModel(), and printCoeffs().
|
protected |
Model coefficients dictionary.
Definition at line 82 of file LESModel.H.
Referenced by LESModel(), printCoeffs(), and read().
|
protected |
Empirical model constant.
Definition at line 87 of file LESModel.H.
Referenced by LESModel(), and read().
|
protected |
|
protected |
|
protected |
|
protected |
Run-time selectable delta model.
Definition at line 107 of file LESModel.H.
Referenced by correct(), LESModel(), and read().