Templated abstract base class for DES models. More...
#include <DESModel.H>


Public Types | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
| Public Types inherited from LESeddyViscosity< BasicTurbulenceModel > | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
| Public Types inherited from eddyViscosity< LESModel< BasicTurbulenceModel > > | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
| Public Types inherited from linearViscousStress< BasicTurbulenceModel > | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
Public Member Functions | |
| DESModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName) | |
| Construct from components. | |
| virtual | ~DESModel ()=default |
| Destructor. | |
| virtual bool | read () |
| Re-read model coefficients if they have changed. | |
| virtual tmp< volScalarField > | LESRegion () const =0 |
| Return the LES field indicator. | |
| virtual tmp< volScalarField > | Ssigma (const volTensorField &gradU) const |
| Return modified strain rate. | |
| virtual tmp< volScalarField > | fd () const |
| Return the shielding function. | |
| Public Member Functions inherited from DESModelBase | |
| ClassName ("DESModelBase") | |
| Runtime type information. | |
| DESModelBase () noexcept=default | |
| Constructor. | |
| virtual | ~DESModelBase ()=default |
| Destructor. | |
| Public Member Functions inherited from LESeddyViscosity< BasicTurbulenceModel > | |
| LESeddyViscosity (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName) | |
| Construct from components. | |
| virtual | ~LESeddyViscosity ()=default |
| Destructor. | |
| Public Member Functions inherited from eddyViscosity< LESModel< BasicTurbulenceModel > > | |
| eddyViscosity (const word &modelName, 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 | ~eddyViscosity ()=default |
| Destructor. | |
| virtual tmp< volScalarField > | nut () const |
| Return the turbulence viscosity. | |
| virtual tmp< volScalarField > | k () const=0 |
| Return the turbulence kinetic energy. | |
| virtual tmp< volSymmTensorField > | R () const |
| Return the Reynolds stress tensor. | |
| virtual void | validate () |
| Validate the turbulence fields after construction. | |
| virtual void | correct ()=0 |
| Solve the turbulence equations and correct the turbulence viscosity. | |
| Public Member Functions inherited from linearViscousStress< BasicTurbulenceModel > | |
| linearViscousStress (const word &modelName, 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 | ~linearViscousStress ()=default |
| Destructor. | |
| virtual tmp< volSymmTensorField > | devRhoReff () const |
| Return the effective stress tensor. | |
| virtual tmp< volSymmTensorField > | devRhoReff (const volVectorField &U) const |
| Return the effective stress tensor based on a given velocity field. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
| Return the source term for the momentum equation. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (const volScalarField &rho, volVectorField &U) const |
| Return the source term for the momentum equation. | |
Static Public Member Functions | |
| static tmp< volScalarField > | Ssigma (const volTensorField &gradU, const dimensionedScalar &coeff) |
| Return modified strain rate. | |
Protected Attributes | |
| dimensionedScalar | Ctrans_ |
| Model-specific transition constant. | |
| Protected Attributes inherited from eddyViscosity< LESModel< BasicTurbulenceModel > > | |
| volScalarField | nut_ |
Additional Inherited Members | |
| Protected Member Functions inherited from eddyViscosity< LESModel< BasicTurbulenceModel > > | |
| virtual void | correctNut ()=0 |
Templated abstract base class for DES models.
Definition at line 51 of file DESModel.H.
| typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 81 of file DESModel.H.
| typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 82 of file DESModel.H.
| typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 83 of file DESModel.H.
| DESModel | ( | const word & | type, |
| const alphaField & | alpha, | ||
| const rhoField & | rho, | ||
| const volVectorField & | U, | ||
| const surfaceScalarField & | alphaRhoPhi, | ||
| const surfaceScalarField & | phi, | ||
| const transportModel & | transport, | ||
| const word & | propertiesName = turbulenceModel::propertiesName ) |
Construct from components.
Definition at line 35 of file DESModel.C.
References alpha, Ctrans_, Foam::dimless, phi, rho, U, and Foam::Zero.
|
virtualdefault |
Destructor.
|
virtual |
Re-read model coefficients if they have changed.
Reimplemented from LESeddyViscosity< BasicTurbulenceModel >.
Definition at line 67 of file DESModel.C.
References Ctrans_, and LESeddyViscosity< BasicTurbulenceModel >::read().

|
pure virtual |
Return the LES field indicator.
Implements DESModelBase.
|
static |
Return modified strain rate.
Definition at line 81 of file DESModel.C.
References Foam::acos(), alpha1, Foam::cos(), Foam::det(), Foam::dimless, Foam::dimTime, Foam::max(), VectorSpace< Form, Cmpt, Mrows *Ncols >::max, Foam::min(), VectorSpace< Form, Cmpt, Mrows *Ncols >::min, Foam::constant::mathematical::pi(), Foam::pow(), Foam::pow3(), Foam::pow4(), Foam::sqr(), Foam::sqrt(), GeometricField< Type, PatchField, GeoMesh >::T(), and Foam::tr().
Referenced by sigma< BasicTurbulenceModel >::correctNut(), and Ssigma().


|
virtual |
Return modified strain rate.
Note: uses Ctrans_ coefficient
Definition at line 152 of file DESModel.C.
References Ctrans_, and Ssigma().

|
virtual |
Return the shielding function.
Implements DESModelBase.
Definition at line 162 of file DESModel.C.
References Foam::dimless, GeometricField< scalar, fvPatchField, volMesh >::New(), IOobjectOption::NO_REGISTER, and Foam::Zero.

|
protected |
Model-specific transition constant.
Definition at line 76 of file DESModel.H.
Referenced by DESModel(), read(), and Ssigma().