k-omega-SST DES turbulence model for incompressible and compressible flows More...
Public Types | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
![]() | |
typedef DESModel< BasicTurbulenceModel > ::alphaField | alphaField |
typedef DESModel< BasicTurbulenceModel > ::rhoField | rhoField |
typedef DESModel< BasicTurbulenceModel > ::transportModel | transportModel |
Public Member Functions | |
TypeName ("kOmegaSSTDES") | |
Runtime type information. More... | |
kOmegaSSTDES (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 | ~kOmegaSSTDES ()=default |
Destructor. More... | |
virtual bool | read () |
Re-read model coefficients if they have changed. More... | |
virtual tmp< volScalarField > | LESRegion () const |
Return the LES field indicator. More... | |
![]() | |
kOmegaSSTBase (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 | ~kOmegaSSTBase ()=default |
Destructor. More... | |
tmp< volScalarField > | DkEff (const volScalarField &F1) const |
Return the effective diffusivity for k. More... | |
tmp< volScalarField > | DomegaEff (const volScalarField &F1) const |
Return the effective diffusivity for omega. More... | |
virtual tmp< volScalarField > | k () const |
Return the turbulence kinetic energy. More... | |
virtual tmp< volScalarField > | omega () const |
Return the turbulence kinetic energy dissipation rate. More... | |
virtual void | correct () |
Solve the turbulence equations and correct the turbulence viscosity. More... | |
Protected Attributes | |
dimensionedScalar | kappa_ |
dimensionedScalar | CDESkom_ |
dimensionedScalar | CDESkeps_ |
![]() | |
dimensionedScalar | alphaK1_ |
dimensionedScalar | alphaK2_ |
dimensionedScalar | alphaOmega1_ |
dimensionedScalar | alphaOmega2_ |
dimensionedScalar | gamma1_ |
dimensionedScalar | gamma2_ |
dimensionedScalar | beta1_ |
dimensionedScalar | beta2_ |
dimensionedScalar | betaStar_ |
dimensionedScalar | a1_ |
dimensionedScalar | b1_ |
dimensionedScalar | c1_ |
Switch | F3_ |
Flag to include the F3 term. More... | |
const volScalarField & | y_ |
Wall distance. More... | |
volScalarField | k_ |
volScalarField | omega_ |
Switch | decayControl_ |
Flag to include the decay control. More... | |
dimensionedScalar | kInf_ |
dimensionedScalar | omegaInf_ |
k-omega-SST DES turbulence model for incompressible and compressible flows
Strelets, M. (2001) Detached Eddy Simulation of Massively Separated Flows, 39th AIAA Aerospace Sciences Meeting and Exhibit, Reno, NV
Definition at line 71 of file kOmegaSSTDES.H.
typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 131 of file kOmegaSSTDES.H.
typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 132 of file kOmegaSSTDES.H.
typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 133 of file kOmegaSSTDES.H.
kOmegaSSTDES | ( | 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 100 of file kOmegaSSTDES.C.
References Foam::type().
|
virtualdefault |
Destructor.
|
inlineprotectedvirtual |
Blending for CDES parameter.
Definition at line 98 of file kOmegaSSTDES.H.
References kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::blend(), kOmegaSSTDES< BasicTurbulenceModel >::CDESkeps_, and kOmegaSSTDES< BasicTurbulenceModel >::CDESkom_.
|
protectedvirtual |
Reimplemented from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >.
Definition at line 41 of file kOmegaSSTDES.C.
|
protectedvirtual |
Reimplemented from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >.
Definition at line 52 of file kOmegaSSTDES.C.
References Foam::fvc::grad(), Foam::magSqr(), and Foam::symm().
|
protectedvirtual |
Length scale.
Reimplemented in kOmegaSSTIDDES< BasicTurbulenceModel >, and kOmegaSSTDDES< BasicTurbulenceModel >.
Definition at line 60 of file kOmegaSSTDES.C.
References delta, k, Foam::min(), and Foam::sqrt().
|
protectedvirtual |
Return epsilon/k.
Reimplemented from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >.
Definition at line 74 of file kOmegaSSTDES.C.
References F1, Foam::mag(), and Foam::sqrt().
|
protectedvirtual |
Return G/nu.
Definition at line 86 of file kOmegaSSTDES.C.
TypeName | ( | "kOmegaSSTDES< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Re-read model coefficients if they have changed.
Reimplemented from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >.
Reimplemented in kOmegaSSTIDDES< BasicTurbulenceModel >, and kOmegaSSTDDES< BasicTurbulenceModel >.
Definition at line 163 of file kOmegaSSTDES.C.
References Foam::read().
|
virtual |
Return the LES field indicator.
Definition at line 179 of file kOmegaSSTDES.C.
References F1, Foam::fvc::grad(), k, Foam::mag(), Foam::neg(), IOobject::NO_READ, IOobject::NO_WRITE, Foam::sqrt(), and U.
|
protected |
Definition at line 90 of file kOmegaSSTDES.H.
|
protected |
Definition at line 91 of file kOmegaSSTDES.H.
Referenced by kOmegaSSTDES< BasicTurbulenceModel >::CDES().
|
protected |
Definition at line 92 of file kOmegaSSTDES.H.
Referenced by kOmegaSSTDES< BasicTurbulenceModel >::CDES().