k-omega-SST DES turbulence model for incompressible and compressible flows. More...
#include <kOmegaSSTDES.H>


Public Types | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
| Public Types inherited from kOmegaSSTBase< DESModel< BasicTurbulenceModel > > | |
| typedef DESModel< BasicTurbulenceModel >::alphaField | alphaField |
| typedef DESModel< BasicTurbulenceModel >::rhoField | rhoField |
| typedef DESModel< BasicTurbulenceModel >::transportModel | transportModel |
Public Member Functions | |
| TypeName ("kOmegaSSTDES") | |
| Runtime type information. | |
| 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. | |
| virtual | ~kOmegaSSTDES ()=default |
| Destructor. | |
| virtual bool | read () |
| Re-read model coefficients if they have changed. | |
| virtual tmp< volScalarField > | lengthScaleRAS () const |
| RAS length scale. | |
| virtual tmp< volScalarField > | lengthScaleLES (const volScalarField &CDES) const |
| LES length scale. | |
| virtual tmp< volScalarField > | LESRegion () const |
| Return the LES field indicator. | |
| Public Member Functions inherited from kOmegaSSTBase< DESModel< BasicTurbulenceModel > > | |
| kOmegaSSTBase (const kOmegaSSTBase &)=delete | |
| No copy construct. | |
| void | operator= (const kOmegaSSTBase &)=delete |
| No copy assignment. | |
| virtual | ~kOmegaSSTBase ()=default |
| Destructor. | |
| tmp< volScalarField > | DkEff (const volScalarField &F1) const |
| Return the effective diffusivity for k. | |
| tmp< volScalarField > | DomegaEff (const volScalarField &F1) const |
| Return the effective diffusivity for omega. | |
| virtual tmp< volScalarField > | k () const |
| Return the turbulence kinetic energy. | |
| virtual tmp< volScalarField > | omega () const |
| Return the turbulence kinetic energy dissipation rate. | |
| virtual void | correct () |
| Solve the turbulence equations and correct the turbulence viscosity. | |
Protected Attributes | |
| Switch | useSigma_ |
| Switch to activate grey-area enhanced sigma-(D)DES. | |
| dimensionedScalar | kappa_ |
| dimensionedScalar | CDESkom_ |
| DES coefficients. | |
| dimensionedScalar | CDESkeps_ |
| Protected Attributes inherited from kOmegaSSTBase< DESModel< BasicTurbulenceModel > > | |
| 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. | |
| const volScalarField & | y_ |
| Wall distance. | |
| volScalarField | k_ |
| Turbulent kinetic energy field [m^2/s^2]. | |
| volScalarField | omega_ |
| Specific dissipation rate field [1/s]. | |
| Switch | decayControl_ |
| Flag to include the decay control. | |
| dimensionedScalar | kInf_ |
| dimensionedScalar | omegaInf_ |
k-omega-SST DES turbulence model for incompressible and compressible flows.
Reference:
Strelets, M. (2001).
Detached Eddy Simulation of Massively Separated Flows.
39th AIAA Aerospace Sciences Meeting and Exhibit, Reno, NV.
Definition at line 67 of file kOmegaSSTDES.H.
| typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 171 of file kOmegaSSTDES.H.
| typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 172 of file kOmegaSSTDES.H.
| typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 173 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 162 of file kOmegaSSTDES.C.
References alpha, CDESkeps_, CDESkom_, Foam::endl(), dimensioned< Type >::getOrAddToDict(), Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, kappa_, kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::kOmegaSSTBase(), phi, rho, U, useSigma_, and WarningInFunction.

|
virtualdefault |
Destructor.
|
inlineprotectedvirtual |
Blending for CDES parameter.
Definition at line 109 of file kOmegaSSTDES.H.
Referenced by kOmegaSSTDDES< BasicTurbulenceModel >::dTilda(), dTilda(), kOmegaSSTIDDES< BasicTurbulenceModel >::dTilda(), epsilonByk(), lengthScaleLES(), LESRegion(), kOmegaSSTDDES< BasicTurbulenceModel >::S2(), and S2().

|
protectedvirtual |
Reimplemented from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >.
Definition at line 35 of file kOmegaSSTDES.C.
References correctNut(), kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::kOmegaSSTBase(), and S2().
Referenced by correctNut(), and correctNut().


|
protectedvirtual |
Definition at line 46 of file kOmegaSSTDES.C.
References correctNut(), Foam::fvc::grad(), Foam::magSqr(), and Foam::symm().

|
protected |
Definition at line 53 of file kOmegaSSTDES.C.
References DimensionedField< Type, GeoMesh >::dimensions(), kappa_, Foam::max(), Foam::min(), Foam::sqr(), Foam::tr(), and kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::y_.

|
protectedvirtual |
Return square of strain rate.
Reimplemented from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >.
Reimplemented in kOmegaSSTDDES< BasicTurbulenceModel >, kOmegaSSTDDES< EddyDiffusivity< fluidThermoCompressibleTurbulenceModel > >, and kOmegaSSTDDES< Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel > >.
Definition at line 71 of file kOmegaSSTDES.C.
References kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::alphaOmega2_, CDES(), dTilda(), F1, kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::F1(), Foam::fvc::grad(), kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::k(), kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::k_, kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::kOmegaSSTBase(), lengthScaleRAS(), Foam::mag(), kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::omega(), kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::omega_, Foam::pos(), tmp< T >::ref(), S2(), Foam::sqr(), and useSigma_.
Referenced by correctNut(), GbyNu(), GbyNu0(), and S2().


|
protectedvirtual |
Return length scale.
Reimplemented in kOmegaSSTDDES< BasicTurbulenceModel >, kOmegaSSTDDES< EddyDiffusivity< fluidThermoCompressibleTurbulenceModel > >, kOmegaSSTDDES< Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel > >, kOmegaSSTIDDES< BasicTurbulenceModel >, kOmegaSSTIDDES< EddyDiffusivity< fluidThermoCompressibleTurbulenceModel > >, and kOmegaSSTIDDES< Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel > >.
Definition at line 108 of file kOmegaSSTDES.C.
References CDES(), lengthScaleLES(), lengthScaleRAS(), and Foam::min().
Referenced by epsilonByk(), LESRegion(), and S2().


|
protectedvirtual |
Return epsilon/k.
Reimplemented from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >.
Definition at line 119 of file kOmegaSSTDES.C.
References CDES(), dTilda(), kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::F1(), kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::k_, Foam::mag(), and Foam::sqrt().

|
protectedvirtual |
Return (G/nu)_0.
Reimplemented from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >.
Reimplemented in kOmegaSSTDDES< BasicTurbulenceModel >, kOmegaSSTDDES< EddyDiffusivity< fluidThermoCompressibleTurbulenceModel > >, and kOmegaSSTDDES< Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel > >.
Definition at line 131 of file kOmegaSSTDES.C.
References GbyNu0(), kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::kOmegaSSTBase(), S2(), and useSigma_.
Referenced by GbyNu(), and GbyNu0().


|
protectedvirtual |
Return G/nu.
Reimplemented from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >.
Reimplemented in kOmegaSSTDDES< BasicTurbulenceModel >, kOmegaSSTDDES< EddyDiffusivity< fluidThermoCompressibleTurbulenceModel > >, and kOmegaSSTDDES< Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel > >.
Definition at line 148 of file kOmegaSSTDES.C.
References kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::F2(), GbyNu0(), and S2().

| TypeName | ( | "kOmegaSSTDES< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Re-read model coefficients if they have changed.
Reimplemented from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >.
Reimplemented in kOmegaSSTDDES< BasicTurbulenceModel >, kOmegaSSTDDES< EddyDiffusivity< fluidThermoCompressibleTurbulenceModel > >, kOmegaSSTDDES< Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel > >, kOmegaSSTIDDES< BasicTurbulenceModel >, kOmegaSSTIDDES< EddyDiffusivity< fluidThermoCompressibleTurbulenceModel > >, and kOmegaSSTIDDES< Foam::EddyDiffusivity< Foam::fluidThermoCompressibleTurbulenceModel > >.
Definition at line 242 of file kOmegaSSTDES.C.
References CDESkeps_, CDESkom_, kappa_, kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::kOmegaSSTBase(), read(), and useSigma_.
Referenced by kOmegaSSTDDES< BasicTurbulenceModel >::read(), read(), and kOmegaSSTIDDES< BasicTurbulenceModel >::read().


|
virtual |
RAS length scale.
Definition at line 260 of file kOmegaSSTDES.C.
References kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::betaStar_, kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::k(), kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::k_, kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::omega(), kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::omega_, and Foam::sqrt().
Referenced by kOmegaSSTDDES< BasicTurbulenceModel >::dTilda(), dTilda(), kOmegaSSTIDDES< BasicTurbulenceModel >::dTilda(), LESRegion(), kOmegaSSTDDES< BasicTurbulenceModel >::S2(), and S2().


|
virtual |
LES length scale.
Definition at line 271 of file kOmegaSSTDES.C.
Referenced by kOmegaSSTDDES< BasicTurbulenceModel >::dTilda(), dTilda(), kOmegaSSTIDDES< BasicTurbulenceModel >::dTilda(), and kOmegaSSTDDES< BasicTurbulenceModel >::S2().


|
virtual |
Return the LES field indicator.
Definition at line 281 of file kOmegaSSTDES.C.
References kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::alphaOmega2_, CDES(), dTilda(), F1, kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::F1(), Foam::fvc::grad(), kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::k(), kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::k_, lengthScaleRAS(), Foam::mag(), Foam::neg(), tmp< T >::New(), kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::omega(), kOmegaSSTBase< DESModel< BasicTurbulenceModel > >::omega_, IOobject::scopedName(), and U.

|
protected |
Switch to activate grey-area enhanced sigma-(D)DES.
Definition at line 91 of file kOmegaSSTDES.H.
Referenced by kOmegaSSTDDES< BasicTurbulenceModel >::GbyNu0(), GbyNu0(), kOmegaSSTDDES< BasicTurbulenceModel >::kOmegaSSTDDES(), kOmegaSSTDES(), read(), kOmegaSSTDDES< BasicTurbulenceModel >::S2(), and S2().
|
protected |
Definition at line 95 of file kOmegaSSTDES.H.
Referenced by kOmegaSSTDES(), r(), and read().
|
protected |
DES coefficients.
Definition at line 100 of file kOmegaSSTDES.H.
Referenced by kOmegaSSTDES(), and read().
|
protected |
Definition at line 101 of file kOmegaSSTDES.H.
Referenced by kOmegaSSTDES(), and read().