SpalartAllmarasDES DES turbulence model for incompressible and compressible flows. More...
#include <SpalartAllmarasDES.H>


Public Types | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
Public Member Functions | |
| 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 bool | read () |
| Re-read model coefficients if they have changed. | |
| virtual tmp< volScalarField > | LESRegion () const |
| Return the LES field indicator. | |
Protected Member Functions | |
| 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 tmp< volScalarField > | dTilda (const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const |
| Return the length scale. | |
| virtual void | correctNut () |
Protected Attributes | |
| 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_ |
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
Reference:
Spalart, P. R., Jou, W. H., Strelets, M., & Allmaras, S. R. (1997).
Comments on the feasibility of LES for wings, and on a hybrid
RANS/LES approach.
Advances in DNS/LES, 1, 4-8.
Including the low Reynolds number correction documented in
Spalart, P. R., Deck, S., Shur, M. L., Squires,
K. D., Strelets, M. K., & Travin, A. (2006).
A new version of detached-eddy simulation,
resistant to ambiguous grid densities.
Theoretical and computational fluid dynamics, 20(3), 181-195.
DOI:10.1007/s00162-006-0015-0
Definition at line 78 of file SpalartAllmarasDES.H.
| typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 163 of file SpalartAllmarasDES.H.
| typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 164 of file SpalartAllmarasDES.H.
| typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 165 of file SpalartAllmarasDES.H.
| 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.
Definition at line 162 of file SpalartAllmarasDES.C.
References alpha, CDES_, Foam::endl(), fwStar_, dimensioned< Type >::getOrAddToDict(), Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, lowReCorrection_, phi, rho, U, useSigma_, and WarningInFunction.

|
virtualdefault |
Destructor.
|
protectedvirtual |
Return the low Reynolds number correction function.
Definition at line 36 of file SpalartAllmarasDES.C.
References Foam::dimless, e, fwStar_, lowReCorrection_, Foam::max(), mesh, Foam::min(), tmp< T >::New(), IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, psi(), IOobject::scopedName(), Foam::sqr(), Foam::sqrt(), timeName, Foam::type(), and Foam::Zero.
Referenced by SpalartAllmarasIDDES< BasicTurbulenceModel >::dTilda(), and psi().


|
protectedvirtual |
Return the LES length scale.
Definition at line 81 of file SpalartAllmarasDES.C.
References CDES_, delta, and psi.
Referenced by SpalartAllmarasDDES< BasicTurbulenceModel >::dTilda(), dTilda(), SpalartAllmarasDDES< BasicTurbulenceModel >::Stilda(), and Stilda().

|
protectedvirtual |
Return the production term.
Reimplemented in SpalartAllmarasDDES< BasicTurbulenceModel >, SpalartAllmarasDDES< EddyDiffusivity< fluidThermoCompressibleTurbulenceModel > >, and SpalartAllmarasDDES< Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel > >.
Definition at line 92 of file SpalartAllmarasDES.C.
References dTilda(), lengthScaleLES(), Foam::max(), Foam::pos(), Foam::sqr(), Stilda(), and useSigma_.
Referenced by Stilda().


|
protectedvirtual |
Return the length scale.
Reimplemented in SpalartAllmarasDDES< BasicTurbulenceModel >, SpalartAllmarasDDES< EddyDiffusivity< fluidThermoCompressibleTurbulenceModel > >, SpalartAllmarasDDES< Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel > >, SpalartAllmarasIDDES< BasicTurbulenceModel >, SpalartAllmarasIDDES< EddyDiffusivity< fluidThermoCompressibleTurbulenceModel > >, and SpalartAllmarasIDDES< Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel > >.
Definition at line 132 of file SpalartAllmarasDES.C.
References lengthScaleLES(), Foam::min(), and tmp< T >::ref().
Referenced by LESRegion(), and Stilda().


|
protectedvirtual |
Definition at line 149 of file SpalartAllmarasDES.C.
References correctNut().
Referenced by correctNut().


| TypeName | ( | "SpalartAllmarasDES< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Re-read model coefficients if they have changed.
Reimplemented in SpalartAllmarasDDES< BasicTurbulenceModel >, SpalartAllmarasDDES< EddyDiffusivity< fluidThermoCompressibleTurbulenceModel > >, SpalartAllmarasDDES< Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel > >, SpalartAllmarasIDDES< BasicTurbulenceModel >, SpalartAllmarasIDDES< EddyDiffusivity< fluidThermoCompressibleTurbulenceModel > >, and SpalartAllmarasIDDES< Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel > >.
Definition at line 242 of file SpalartAllmarasDES.C.
References CDES_, fwStar_, lowReCorrection_, read(), and useSigma_.
Referenced by SpalartAllmarasDDES< BasicTurbulenceModel >::read(), read(), and SpalartAllmarasIDDES< BasicTurbulenceModel >::read().


|
virtual |
Return the LES field indicator.
Definition at line 260 of file SpalartAllmarasDES.C.
References dTilda(), Foam::fvc::grad(), Foam::neg(), tmp< T >::New(), IOobjectOption::NO_READ, IOobjectOption::NO_WRITE, and IOobject::scopedName().

|
protected |
Switch to activate grey-area enhanced sigma-(D)DES.
Definition at line 102 of file SpalartAllmarasDES.H.
Referenced by read(), SpalartAllmarasDDES< BasicTurbulenceModel >::SpalartAllmarasDDES(), SpalartAllmarasDES(), SpalartAllmarasDDES< BasicTurbulenceModel >::Stilda(), and Stilda().
|
protected |
DES coefficient.
Definition at line 109 of file SpalartAllmarasDES.H.
Referenced by SpalartAllmarasIDDES< BasicTurbulenceModel >::dTilda(), lengthScaleLES(), read(), and SpalartAllmarasDES().
|
protected |
Flag for low Reynolds number correction.
Definition at line 114 of file SpalartAllmarasDES.H.
Referenced by psi(), read(), and SpalartAllmarasDES().
|
protected |
Definition at line 115 of file SpalartAllmarasDES.H.
Referenced by psi(), read(), and SpalartAllmarasDES().