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 () |
Destructor. More... | |
virtual bool | read () |
Re-read model coefficients if they have changed. More... | |
virtual void | correct () |
Solve the turbulence equations and correct the turbulence viscosity. 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 () |
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 > | epsilon () const |
Return the turbulence kinetic energy dissipation rate. More... | |
virtual tmp< volScalarField > | omega () const |
Return the turbulence kinetic energy dissipation rate. More... | |
Private Member Functions | |
kOmegaSSTDES (const kOmegaSSTDES &) | |
kOmegaSSTDES & | operator= (const kOmegaSSTDES &) |
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
Note: DES constants CDES_kom and CDES_keps have been re-calibrated to OpenFOAM numerics via decaying isotropic turbulence test case
Definition at line 65 of file kOmegaSSTDES.H.
typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 108 of file kOmegaSSTDES.H.
typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 109 of file kOmegaSSTDES.H.
typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 110 of file kOmegaSSTDES.H.
|
private |
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 73 of file kOmegaSSTDES.C.
References Foam::type().
|
inlinevirtual |
Destructor.
Definition at line 134 of file kOmegaSSTDES.H.
|
private |
|
inlineprotectedvirtual |
Blending for CDES parameter.
Definition at line 90 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 38 of file kOmegaSSTDES.C.
|
protectedvirtual |
Reimplemented from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >.
Definition at line 49 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 57 of file kOmegaSSTDES.C.
References delta, k, Foam::min(), and Foam::sqrt().
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 136 of file kOmegaSSTDES.C.
References Foam::read().
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Reimplemented from kOmegaSSTBase< DESModel< BasicTurbulenceModel > >.
Definition at line 154 of file kOmegaSSTDES.C.
References Foam::fvc::absolute(), Foam::constant::atomic::alpha, beta(), Foam::bound(), GeometricField::boundaryField(), tmp::clear(), eddyViscosity< LESModel< BasicTurbulenceModel > >::correct(), Foam::fvm::ddt(), Foam::dev(), Foam::fvm::div(), Foam::fvc::div(), divU(), F1, Foam::constant::universal::G, Foam::fvc::grad(), k, Foam::fvm::laplacian(), Foam::mag(), Foam::magSqr(), Foam::min(), nut, phi, rho, Foam::solve(), Foam::fvm::Sp(), Foam::sqrt(), Foam::fvm::SuSp(), Foam::symm(), Foam::twoSymm(), U, and GeometricField::GeometricBoundaryField::updateCoeffs().
|
virtual |
Return the LES field indicator.
Definition at line 244 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 82 of file kOmegaSSTDES.H.
|
protected |
Definition at line 83 of file kOmegaSSTDES.H.
Referenced by kOmegaSSTDES< BasicTurbulenceModel >::CDES().
|
protected |
Definition at line 84 of file kOmegaSSTDES.H.
Referenced by kOmegaSSTDES< BasicTurbulenceModel >::CDES().
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.