Implementation of the k-omega-SST turbulence model for incompressible flows. More...
#include <incompressibleRASModels/kOmegaSST.H>
Implementation of the k-omega-SST turbulence model for incompressible flows.
Turbulence model described in:
Menter, F., Esch, T. "Elements of Industrial Heat Transfer Prediction" 16th Brazilian Congress of Mechanical Engineering (COBEM), Nov. 2001
Note that this implementation is written in terms of alpha diffusion coefficients rather than the more traditional sigma (alpha = 1/sigma) so that the blending can be applied to all coefficuients in a consistent manner. The paper suggests that sigma is blended but this would not be consistent with the blending of the k-epsilon and k-omega models.
Also note that the error in the last term of equation (2) relating to sigma has been corrected.
Wall-functions are applied in this implementation by using equations (14) to specify the near-wall omega as appropriate.
The blending functions (15) and (16) are not currently used because of the uncertainty in their origin, range of applicability and that is y+ becomes sufficiently small blending u_tau in this manner clearly becomes nonsense.
The default model coefficients correspond to the following:
kOmegaSSTCoeffs { alphaK1 0.85034; alphaK2 1.0; alphaOmega1 0.5; alphaOmega2 0.85616; beta1 0.075; beta2 0.0828; betaStar0.09; gamma1 0.5532; gamma2 0.4403; a1 0.31; c1 10.0; }
Definition at line 97 of file kOmegaSST.H.
Public Member Functions | |
TypeName ("kOmegaSST") | |
Runtime type information.
| |
kOmegaSST (const volVectorField &U, const surfaceScalarField &phi, transportModel &transport) | |
Construct from components.
| |
virtual | ~kOmegaSST () |
Destructor.
| |
virtual tmp< volScalarField > | nut () const |
Return the turbulence viscosity.
| |
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 specific dissipation rate.
| |
virtual tmp< volScalarField > | epsilon () const |
Return the turbulence kinetic energy dissipation rate.
| |
virtual tmp< volSymmTensorField > | R () const |
Return the Reynolds stress tensor.
| |
virtual tmp< volSymmTensorField > | devReff () const |
Return the effective stress tensor including the laminar stress.
| |
virtual tmp< fvVectorMatrix > | divDevReff (volVectorField &U) const |
Return the source term for the momentum equation.
| |
virtual void | correct () |
Solve the turbulence equations and correct the turbulence viscosity.
| |
virtual bool | read () |
Read RASProperties dictionary.
|
kOmegaSST | ( | const volVectorField & | U, |
const surfaceScalarField & | phi, | ||
transportModel & | transport | ||
) |
Construct from components.
Definition at line 91 of file kOmegaSST.C.
References Foam::bound(), F2, Foam::fvc::grad(), Foam::mag(), Foam::max(), Foam::sqrt(), and Foam::symm().
virtual ~kOmegaSST | ( | ) | [inline, virtual]
|
Destructor.
Definition at line 198 of file kOmegaSST.H.
TypeName | ( | "kOmegaSST" | ) |
Runtime type information.
virtual tmp<volScalarField> nut | ( | ) | const [inline, virtual]
|
Return the turbulence viscosity.
Implements incompressible::RASModel.
Definition at line 205 of file kOmegaSST.H.
tmp<volScalarField> DkEff | ( | const volScalarField & | F1 ) | const [inline]
|
Return the effective diffusivity for k.
Definition at line 211 of file kOmegaSST.H.
References incompressible::turbulenceModel::nu().
Referenced by kOmegaSST::correct().
tmp<volScalarField> DomegaEff | ( | const volScalarField & | F1 ) | const [inline]
|
Return the effective diffusivity for omega.
Definition at line 220 of file kOmegaSST.H.
References incompressible::turbulenceModel::nu().
Referenced by kOmegaSST::correct().
virtual tmp<volScalarField> k | ( | ) | const [inline, virtual]
|
Return the turbulence kinetic energy.
Implements incompressible::RASModel.
Definition at line 229 of file kOmegaSST.H.
virtual tmp<volScalarField> omega | ( | ) | const [inline, virtual]
|
Return the turbulence specific dissipation rate.
Definition at line 235 of file kOmegaSST.H.
virtual tmp<volScalarField> epsilon | ( | ) | const [inline, virtual]
|
Return the turbulence kinetic energy dissipation rate.
Implements incompressible::RASModel.
Definition at line 241 of file kOmegaSST.H.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), incompressible::turbulenceModel::mesh_, fvMesh::time(), and Time::timeName().
tmp< volSymmTensorField > R | ( | ) | const [virtual]
|
Return the Reynolds stress tensor.
Implements incompressible::RASModel.
Definition at line 249 of file kOmegaSST.C.
References GeometricField< Type, PatchField, GeoMesh >::boundaryField(), Foam::fvc::grad(), Foam::I, incompressible::turbulenceModel::mesh_, IOobject::NO_READ, IOobject::NO_WRITE, incompressible::turbulenceModel::runTime_, Time::timeName(), Foam::twoSymm(), and incompressible::turbulenceModel::U_.
tmp< volSymmTensorField > devReff | ( | ) | const [virtual]
|
Return the effective stress tensor including the laminar stress.
Implements incompressible::RASModel.
Definition at line 270 of file kOmegaSST.C.
References Foam::dev(), Foam::fvc::grad(), incompressible::turbulenceModel::mesh_, IOobject::NO_READ, IOobject::NO_WRITE, incompressible::RASModel::nuEff(), incompressible::turbulenceModel::runTime_, Time::timeName(), Foam::twoSymm(), and incompressible::turbulenceModel::U_.
tmp< fvVectorMatrix > divDevReff | ( | volVectorField & | U ) | const [virtual]
|
Return the source term for the momentum equation.
Implements incompressible::RASModel.
Definition at line 290 of file kOmegaSST.C.
References Foam::dev(), Foam::fvc::div(), Foam::fvc::grad(), Foam::fvm::laplacian(), incompressible::RASModel::nuEff(), and Foam::T().
void correct | ( | ) | [virtual]
|
Solve the turbulence equations and correct the turbulence viscosity.
Implements incompressible::RASModel.
Definition at line 325 of file kOmegaSST.C.
References Foam::bound(), polyMesh::changing(), GeometricField< Type, PatchField, GeoMesh >::correctBoundaryConditions(), Foam::fvm::ddt(), Foam::fvc::div(), Foam::fvm::div(), kOmegaSST::DkEff(), kOmegaSST::DomegaEff(), Foam::fvc::grad(), incompressible::RASModel::k0_, Foam::fvm::laplacian(), Foam::magSqr(), Foam::max(), incompressible::turbulenceModel::mesh_, Foam::min(), incompressible::RASModel::omega0_, incompressible::turbulenceModel::phi_, Foam::solve(), Foam::fvm::Sp(), Foam::sqrt(), Foam::fvm::SuSp(), Foam::symm(), incompressible::RASModel::turbulence_, and incompressible::turbulenceModel::U_.
bool read | ( | ) | [virtual]
|
Read RASProperties dictionary.
Implements incompressible::RASModel.
Definition at line 300 of file kOmegaSST.C.
References incompressible::RASModel::coeffDict(), and dimensioned< Type >::readIfPresent().