SpalartAllmarasDES DES turbulence model for incompressible and compressible flows. More...
Public Types | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
![]() | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
![]() | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
![]() | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
![]() | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
Public Member Functions | |
TypeName ("SpalartAllmarasDES") | |
Runtime type information. More... | |
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. More... | |
virtual | ~SpalartAllmarasDES ()=default |
Destructor. More... | |
virtual bool | read () |
Re-read model coefficients if they have changed. More... | |
tmp< volScalarField > | DnuTildaEff () const |
Return the effective diffusivity for nuTilda. More... | |
virtual tmp< volScalarField > | k () const |
Return SGS kinetic energy. More... | |
tmp< volScalarField > | nuTilda () const |
virtual tmp< volScalarField > | LESRegion () const |
Return the LES field indicator. More... | |
virtual void | correct () |
Correct nuTilda and related properties. More... | |
![]() | |
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. More... | |
virtual | ~DESModel () |
Destructor. More... | |
virtual tmp< volScalarField > | LESRegion () const =0 |
Return the LES field indicator. More... | |
![]() | |
DESModelBase () | |
Constructor. More... | |
virtual | ~DESModelBase ()=default |
Destructor. More... | |
ClassName ("DESModelBase") | |
virtual tmp< volScalarField > | LESRegion () const =0 |
Return the LES field indicator. More... | |
![]() | |
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. More... | |
virtual | ~LESeddyViscosity ()=default |
Destructor. More... | |
virtual bool | read () |
Read model coefficients if they have changed. More... | |
![]() | |
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. More... | |
virtual | ~eddyViscosity ()=default |
Destructor. More... | |
virtual bool | read ()=0 |
Re-read model coefficients if they have changed. More... | |
virtual tmp< volScalarField > | nut () const |
Return the turbulence viscosity. More... | |
virtual tmp< scalarField > | nut (const label patchi) const |
Return the turbulence viscosity on patch. More... | |
virtual tmp< volScalarField > | k () const=0 |
Return the turbulence kinetic energy. More... | |
virtual tmp< volSymmTensorField > | R () const |
Return the Reynolds stress tensor. More... | |
virtual void | validate () |
Validate the turbulence fields after construction. More... | |
virtual void | correct ()=0 |
Solve the turbulence equations and correct the turbulence viscosity. More... | |
![]() | |
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. More... | |
virtual | ~linearViscousStress ()=default |
Destructor. More... | |
virtual bool | read ()=0 |
Re-read model coefficients if they have changed. More... | |
virtual tmp< volSymmTensorField > | devRhoReff () const |
Return the effective stress tensor. More... | |
virtual tmp< volSymmTensorField > | devRhoReff (const volVectorField &U) const |
Return the effective stress tensor based on a given velocity field. More... | |
virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
Return the source term for the momentum equation. More... | |
virtual tmp< fvVectorMatrix > | divDevRhoReff (const volScalarField &rho, volVectorField &U) const |
Return the source term for the momentum equation. More... | |
virtual void | correct ()=0 |
Solve the turbulence equations and correct the turbulence viscosity. More... | |
Protected Attributes | |
dimensionedScalar | sigmaNut_ |
dimensionedScalar | kappa_ |
dimensionedScalar | Cb1_ |
dimensionedScalar | Cb2_ |
dimensionedScalar | Cw1_ |
dimensionedScalar | Cw2_ |
dimensionedScalar | Cw3_ |
dimensionedScalar | Cv1_ |
dimensionedScalar | Cs_ |
dimensionedScalar | CDES_ |
dimensionedScalar | ck_ |
Switch | lowReCorrection_ |
dimensionedScalar | Ct3_ |
dimensionedScalar | Ct4_ |
dimensionedScalar | fwStar_ |
volScalarField | nuTilda_ |
const volScalarField & | y_ |
Wall distance. More... | |
![]() | |
volScalarField | nut_ |
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.Kh, Travin, A. (2006). A new version of detached-eddy simulation, resistant to ambiguous grid densities. Theor. Comput. Fluid Dyn., 20, 181-195.
Definition at line 80 of file SpalartAllmarasDES.H.
typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 188 of file SpalartAllmarasDES.H.
typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 189 of file SpalartAllmarasDES.H.
typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 190 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 240 of file SpalartAllmarasDES.C.
References Foam::type().
|
virtualdefault |
Destructor.
|
protected |
Definition at line 42 of file SpalartAllmarasDES.C.
References nu.
|
protected |
Definition at line 49 of file SpalartAllmarasDES.C.
References Foam::pow3().
|
protected |
Definition at line 60 of file SpalartAllmarasDES.C.
|
protected |
Definition at line 71 of file SpalartAllmarasDES.C.
References Foam::exp(), and Foam::sqr().
|
protected |
Definition at line 81 of file SpalartAllmarasDES.C.
References Foam::mag(), Foam::skew(), and Foam::sqrt().
|
protected |
Definition at line 91 of file SpalartAllmarasDES.C.
References Foam::max(), and Foam::sqr().
|
protected |
Definition at line 112 of file SpalartAllmarasDES.C.
References DimensionedField< Type, GeoMesh >::dimensions(), Foam::max(), Foam::min(), Foam::sqr(), and Foam::tr().
|
protected |
Definition at line 142 of file SpalartAllmarasDES.C.
References g, Foam::pow(), and Foam::pow6().
|
protected |
Definition at line 156 of file SpalartAllmarasDES.C.
References Foam::dimless, e, Foam::max(), mesh, Foam::min(), IOobject::NO_READ, IOobject::NO_WRITE, psi, tmp< T >::ref(), Foam::sqr(), Foam::sqrt(), timeName, and Foam::type().
|
protectedvirtual |
Length scale.
Reimplemented in SpalartAllmarasDDES< BasicTurbulenceModel >, and SpalartAllmarasIDDES< BasicTurbulenceModel >.
Definition at line 203 of file SpalartAllmarasDES.C.
References delta, Foam::min(), psi, and tmp< T >::ref().
|
protected |
Definition at line 217 of file SpalartAllmarasDES.C.
References GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), and Time::New().
|
protectedvirtual |
Implements eddyViscosity< LESModel< BasicTurbulenceModel > >.
Definition at line 231 of file SpalartAllmarasDES.C.
TypeName | ( | "SpalartAllmarasDES< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Re-read model coefficients if they have changed.
Reimplemented from LESeddyViscosity< BasicTurbulenceModel >.
Reimplemented in SpalartAllmarasDDES< BasicTurbulenceModel >, and SpalartAllmarasIDDES< BasicTurbulenceModel >.
Definition at line 417 of file SpalartAllmarasDES.C.
References dimensioned< Type >::readIfPresent(), and Foam::sqr().
tmp< volScalarField > DnuTildaEff |
Return the effective diffusivity for nuTilda.
Definition at line 448 of file SpalartAllmarasDES.C.
References nu.
|
virtual |
Return SGS kinetic energy.
Implements eddyViscosity< LESModel< BasicTurbulenceModel > >.
Definition at line 459 of file SpalartAllmarasDES.C.
References GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), Foam::dimLength, Foam::fvc::grad(), IOobject::NO_READ, IOobject::NO_WRITE, nut, tmp< T >::ref(), Foam::sqr(), and Foam::Zero.
|
inline |
Definition at line 228 of file SpalartAllmarasDES.H.
References SpalartAllmarasDES< BasicTurbulenceModel >::nuTilda_.
|
virtual |
Return the LES field indicator.
Implements DESModel< BasicTurbulenceModel >.
Definition at line 490 of file SpalartAllmarasDES.C.
References Foam::fvc::grad(), Foam::neg(), IOobject::NO_READ, and IOobject::NO_WRITE.
|
virtual |
Correct nuTilda and related properties.
Implements eddyViscosity< LESModel< BasicTurbulenceModel > >.
Definition at line 516 of file SpalartAllmarasDES.C.
References alpha, Foam::bound(), kEpsilonLopesdaCosta< BasicTurbulenceModel >::correct(), Foam::fvm::ddt(), Foam::fvm::div(), fvOptions, Foam::fvc::grad(), Foam::fvm::laplacian(), Foam::magSqr(), Time::New(), tmp< T >::ref(), rho, solve(), Foam::fvm::Sp(), Foam::sqr(), U, and Foam::Zero.
|
protected |
Definition at line 99 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 100 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 102 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 103 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 104 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 105 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 106 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 107 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 108 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 109 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 110 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 115 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 116 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 117 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 118 of file SpalartAllmarasDES.H.
|
protected |
Definition at line 123 of file SpalartAllmarasDES.H.
Referenced by SpalartAllmarasDES< BasicTurbulenceModel >::nuTilda().
|
protected |
Wall distance.
Note: different to wall distance in parent RASModel which is for near-wall cells only
Definition at line 128 of file SpalartAllmarasDES.H.