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 | 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 BasicTurbulenceModel::alphaField alphaField |
Definition at line 109 of file LESModel.H.
| typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 110 of file LESModel.H.
| typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 111 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, | ||
| const word & | propertiesName | ||
| ) |
Construct from components.
Definition at line 44 of file LESModel.C.
References LESModel< BasicTurbulenceModel >::New().

|
inlinevirtual |
Destructor.
Definition at line 170 of file LESModel.H.
References LESModel< BasicTurbulenceModel >::read().

|
protectedvirtual |
Print model coefficients.
Definition at line 31 of file LESModel.C.
References Foam::endl(), and Foam::Info.

| 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 126 of file LESModel.C.
References Foam::constant::atomic::alpha, TimePaths::constant(), IOobject::db(), Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, IOobject::group(), Foam::Info, dictionary::lookup(), Foam::nl, phi, rho, dictionary::subDict(), IOobject::time(), and U.
Referenced by LESModel< BasicTurbulenceModel >::LESModel().


|
virtual |
Read model coefficients if they have changed.
Reimplemented in SpalartAllmarasDES< BasicTurbulenceModel >, dynamicKEqn< BasicTurbulenceModel >, NicenoKEqn< BasicTurbulenceModel >, Smagorinsky< BasicTurbulenceModel >, SpalartAllmarasIDDES< BasicTurbulenceModel >, SmagorinskyZhang< BasicTurbulenceModel >, DeardorffDiffStress< BasicTurbulenceModel >, WALE< BasicTurbulenceModel >, continuousGasKEqn< BasicTurbulenceModel >, kEqn< BasicTurbulenceModel >, dynamicLagrangian< BasicTurbulenceModel >, ReynoldsStress< LESModel< BasicTurbulenceModel > >, LESeddyViscosity< BasicTurbulenceModel >, eddyViscosity< LESModel< BasicTurbulenceModel > >, and linearViscousStress< LESModel< BasicTurbulenceModel > >.
Definition at line 179 of file LESModel.C.
References Foam::blockMeshTools::read().
Referenced by LESModel< BasicTurbulenceModel >::~LESModel().


|
inlinevirtual |
Const access to the coefficients dictionary.
Definition at line 183 of file LESModel.H.
References LESModel< BasicTurbulenceModel >::coeffDict_.
|
inline |
Return the lower allowable limit for k (default: SMALL)
Definition at line 189 of file LESModel.H.
References LESModel< BasicTurbulenceModel >::kMin_.
|
inline |
Allow kMin to be changed.
Definition at line 195 of file LESModel.H.
References LESModel< BasicTurbulenceModel >::kMin_.
|
inline |
Access function to filter width.
Definition at line 201 of file LESModel.H.
References LESModel< BasicTurbulenceModel >::delta_.
Referenced by dynamicLagrangian< BasicTurbulenceModel >::k().

|
inlinevirtual |
Return the effective viscosity.
Definition at line 208 of file LESModel.H.
References IOobject::groupName(), nu, and nut.

|
inlinevirtual |
Return the effective viscosity on patch.
Definition at line 221 of file LESModel.H.
References LESModel< BasicTurbulenceModel >::correct(), nu, and nut.

|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Reimplemented in SpalartAllmarasDES< BasicTurbulenceModel >, dynamicKEqn< BasicTurbulenceModel >, Smagorinsky< BasicTurbulenceModel >, kEqn< BasicTurbulenceModel >, dynamicLagrangian< BasicTurbulenceModel >, WALE< BasicTurbulenceModel >, ReynoldsStress< LESModel< BasicTurbulenceModel > >, DeardorffDiffStress< BasicTurbulenceModel >, eddyViscosity< LESModel< BasicTurbulenceModel > >, and linearViscousStress< LESModel< BasicTurbulenceModel > >.
Definition at line 202 of file LESModel.C.
References correct.
Referenced by LESModel< BasicTurbulenceModel >::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< BasicTurbulenceModel >::coeffDict().
|
protected |
Lower limit of k.
Definition at line 78 of file LESModel.H.
Referenced by LESModel< BasicTurbulenceModel >::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< BasicTurbulenceModel >::delta().
1.8.13