SpalartAllmaras IDDES turbulence model for incompressible and compressible flows. More...
#include <SpalartAllmarasIDDES.H>


Public Types | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
| Public Types inherited from SpalartAllmarasDES< BasicTurbulenceModel > | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
Public Member Functions | |
| TypeName ("SpalartAllmarasIDDES") | |
| Runtime type information. | |
| SpalartAllmarasIDDES (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName) | |
| Construct from components. | |
| virtual | ~SpalartAllmarasIDDES ()=default |
| Destructor. | |
| virtual bool | read () |
| Re-read model coefficients if they have changed. | |
| virtual tmp< volScalarField > | fd () const |
| Return the shielding function. | |
| Public Member Functions inherited from SpalartAllmarasDES< BasicTurbulenceModel > | |
| TypeName ("SpalartAllmarasDES") | |
| Runtime type information. | |
| SpalartAllmarasDES (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName) | |
| Construct from components. | |
| virtual | ~SpalartAllmarasDES ()=default |
| Destructor. | |
| virtual tmp< volScalarField > | LESRegion () const |
| Return the LES field indicator. | |
Protected Member Functions | |
| virtual tmp< volScalarField > | dTilda (const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const |
| Return the length scale. | |
| Protected Member Functions inherited from SpalartAllmarasDES< BasicTurbulenceModel > | |
| virtual tmp< volScalarField > | psi (const volScalarField &chi, const volScalarField &fv1) const |
| Return the low Reynolds number correction function. | |
| virtual tmp< volScalarField > | lengthScaleLES (const volScalarField &chi, const volScalarField &fv1) const |
| Return the LES length scale. | |
| virtual tmp< volScalarField > | Stilda (const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU, const volScalarField &dTilda) const |
| Return the production term. | |
| virtual void | correctNut () |
Protected Attributes | |
| dimensionedScalar | Cdt1_ |
| dimensionedScalar | Cdt2_ |
| dimensionedScalar | Cl_ |
| dimensionedScalar | Ct_ |
| Switch | fe_ |
| const IDDESDelta & | IDDESDelta_ |
| IDDES delta. | |
| Protected Attributes inherited from SpalartAllmarasDES< BasicTurbulenceModel > | |
| Switch | useSigma_ |
| Switch to activate grey-area enhanced sigma-(D)DES. | |
| dimensionedScalar | CDES_ |
| DES coefficient. | |
| Switch | lowReCorrection_ |
| Flag for low Reynolds number correction. | |
| dimensionedScalar | fwStar_ |
SpalartAllmaras IDDES turbulence model for incompressible and compressible flows.
Reference:
Shur, M. L., Spalart, P. R., Strelets, M. K., & Travin, A. K. (2008).
A hybrid RANS-LES approach with delayed-DES
and wall-modelled LES capabilities.
International journal of heat and fluid flow, 29(6), 1638-1649.
DOI:10.1016/j.ijheatfluidflow.2008.07.001
Definition at line 64 of file SpalartAllmarasIDDES.H.
| typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 131 of file SpalartAllmarasIDDES.H.
| typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 132 of file SpalartAllmarasIDDES.H.
| typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 133 of file SpalartAllmarasIDDES.H.
| SpalartAllmarasIDDES | ( | const alphaField & | alpha, |
| const rhoField & | rho, | ||
| const volVectorField & | U, | ||
| const surfaceScalarField & | alphaRhoPhi, | ||
| const surfaceScalarField & | phi, | ||
| const transportModel & | transport, | ||
| const word & | propertiesName = turbulenceModel::propertiesName, | ||
| const word & | type = typeName ) |
Construct from components.
Definition at line 134 of file SpalartAllmarasIDDES.C.
References Cdt1_, Cdt2_, Cl_, Ct_, fe_, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, IDDESDelta_, phi, rho, and U.
|
virtualdefault |
Destructor.
|
protectedvirtual |
Return the length scale.
Reimplemented from SpalartAllmarasDES< BasicTurbulenceModel >.
Definition at line 87 of file SpalartAllmarasIDDES.C.
References alpha, SpalartAllmarasDES< BasicTurbulenceModel >::CDES_, delta, Foam::dimLength, Foam::exp(), fe_, Foam::lerp(), Foam::mag(), Foam::max(), Foam::min(), Foam::pos0(), Foam::pow(), SpalartAllmarasDES< BasicTurbulenceModel >::psi(), psi, and Foam::sqr().

| TypeName | ( | "SpalartAllmarasIDDES< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Re-read model coefficients if they have changed.
Reimplemented from SpalartAllmarasDES< BasicTurbulenceModel >.
Definition at line 216 of file SpalartAllmarasIDDES.C.
References Cdt1_, Cdt2_, Cl_, Ct_, and SpalartAllmarasDES< BasicTurbulenceModel >::read().

|
virtual |
Return the shielding function.
Definition at line 233 of file SpalartAllmarasIDDES.C.
References alpha, Foam::exp(), Foam::fvc::grad(), Foam::mag(), Foam::max(), Foam::min(), Foam::pow(), and Foam::sqr().

|
protected |
Definition at line 103 of file SpalartAllmarasIDDES.H.
Referenced by read(), and SpalartAllmarasIDDES().
|
protected |
Definition at line 104 of file SpalartAllmarasIDDES.H.
Referenced by read(), and SpalartAllmarasIDDES().
|
protected |
Definition at line 105 of file SpalartAllmarasIDDES.H.
Referenced by read(), and SpalartAllmarasIDDES().
|
protected |
Definition at line 106 of file SpalartAllmarasIDDES.H.
Referenced by read(), and SpalartAllmarasIDDES().
|
protected |
Definition at line 107 of file SpalartAllmarasIDDES.H.
Referenced by dTilda(), and SpalartAllmarasIDDES().
|
protected |
IDDES delta.
Definition at line 113 of file SpalartAllmarasIDDES.H.
Referenced by SpalartAllmarasIDDES().