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


Public Types | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
| Public Types inherited from kOmegaSSTDES< BasicTurbulenceModel > | |
| 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 ("kOmegaSSTIDDES") | |
| Runtime type information. | |
| kOmegaSSTIDDES (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 | ~kOmegaSSTIDDES ()=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 kOmegaSSTDES< BasicTurbulenceModel > | |
| 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 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. | |
k-omega-SST IDDES turbulence model for incompressible and compressible flows.
Reference:
Gritskevich, M. S., Garbaruk, A. V.,
Schütze, J., & Menter, F. R. (2012).
Development of DDES and IDDES formulations for
the k-ω shear stress transport model.
Flow, turbulence and combustion, 88(3), 431-449.
DOI:10.1007/s10494-011-9378-4
Definition at line 64 of file kOmegaSSTIDDES.H.
| typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 128 of file kOmegaSSTIDDES.H.
| typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 129 of file kOmegaSSTIDDES.H.
| typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 130 of file kOmegaSSTIDDES.H.
| kOmegaSSTIDDES | ( | 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 129 of file kOmegaSSTIDDES.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 kOmegaSSTDES< BasicTurbulenceModel >.
Definition at line 87 of file kOmegaSSTIDDES.C.
References alpha, kOmegaSSTDES< BasicTurbulenceModel >::CDES(), Foam::dimLength, Foam::exp(), fe_, kOmegaSSTDES< BasicTurbulenceModel >::lengthScaleLES(), kOmegaSSTDES< BasicTurbulenceModel >::lengthScaleRAS(), Foam::lerp(), Foam::max(), Foam::min(), Foam::pos0(), Foam::pow(), and Foam::sqr().

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

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

|
protected |
Definition at line 101 of file kOmegaSSTIDDES.H.
Referenced by kOmegaSSTIDDES(), and read().
|
protected |
Definition at line 102 of file kOmegaSSTIDDES.H.
Referenced by kOmegaSSTIDDES(), and read().
|
protected |
Definition at line 103 of file kOmegaSSTIDDES.H.
Referenced by kOmegaSSTIDDES(), and read().
|
protected |
Definition at line 104 of file kOmegaSSTIDDES.H.
Referenced by kOmegaSSTIDDES(), and read().
|
protected |
Definition at line 105 of file kOmegaSSTIDDES.H.
Referenced by dTilda(), and kOmegaSSTIDDES().
|
protected |