Implementation of the k-omega-SST turbulence model for incompressible and compressible flows. More...
Public Types | |
typedef BasicTurbulenceModel::alphaField | alphaField |
typedef BasicTurbulenceModel::rhoField | rhoField |
typedef BasicTurbulenceModel::transportModel | transportModel |
![]() | |
typedef Alpha | alphaField |
typedef Rho | rhoField |
typedef TransportModel | transportModel |
Public Member Functions | |
kOmegaSST (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 | ~kOmegaSST () |
Destructor. More... | |
virtual bool | read () |
Re-read model coefficients if they have changed. 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... | |
virtual void | correct () |
Solve the turbulence equations and correct the turbulence viscosity. More... | |
![]() | |
declareRunTimeNewSelectionTable (autoPtr, TurbulenceModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName),(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)) | |
TurbulenceModel (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName) | |
Construct. More... | |
virtual | ~TurbulenceModel () |
Destructor. More... | |
const alphaField & | alpha () const |
Access function to phase fraction. More... | |
const transportModel & | transport () const |
Access function to incompressible transport model. More... | |
tmp< volScalarField > | nu () const |
Return the laminar viscosity. More... | |
tmp< scalarField > | nu (const label patchi) const |
Return the laminar viscosity on patchi. More... | |
Protected Attributes | |
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_ |
const volScalarField & | y_ |
Wall distance. More... | |
volScalarField | k_ |
volScalarField | omega_ |
![]() | |
const alphaField & | alpha_ |
const transportModel & | transport_ |
Additional Inherited Members | |
![]() | |
static autoPtr< TurbulenceModel > | New (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName) |
Return a reference to the selected turbulence model. More... | |
Implementation of the k-omega-SST turbulence model for incompressible and compressible flows.
Turbulence model described in:
Menter, F. R. & Esch, T. (2001). Elements of Industrial Heat Transfer Prediction. 16th Brazilian Congress of Mechanical Engineering (COBEM).
with updated coefficients from
Menter, F. R., Kuntz, M., and Langtry, R. (2003). Ten Years of Industrial Experience with the SST Turbulence Model. Turbulence, Heat and Mass Transfer 4, ed: K. Hanjalic, Y. Nagano, & M. Tummers, Begell House, Inc., 625 - 632.
but with the consistent production terms from the 2001 paper as form in the 2003 paper is a typo, see
http://turbmodels.larc.nasa.gov/sst.html
and the addition of the optional F3 term for rough walls from
Hellsten, A. (1998). "Some Improvements in Menter’s k-omega-SST turbulence model" 29th AIAA Fluid Dynamics Conference, AIAA-98-2554.
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 if y+ becomes sufficiently small blending u_tau in this manner clearly becomes nonsense.
The default model coefficients are
kOmegaSSTCoeffs { alphaK1 0.85; alphaK2 1.0; alphaOmega1 0.5; alphaOmega2 0.856; beta1 0.075; beta2 0.0828; betaStar 0.09; gamma1 5/9; gamma2 0.44; a1 0.31; b1 1.0; c1 10.0; F3 no; }
Definition at line 115 of file kOmegaSSTBase.H.
typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 229 of file kOmegaSSTBase.H.
typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 230 of file kOmegaSSTBase.H.
typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 231 of file kOmegaSSTBase.H.
kOmegaSST | ( | 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.
Definition at line 201 of file kOmegaSSTBase.C.
References Foam::bound().
|
inlinevirtual |
Destructor.
Reimplemented in kOmegaSST< BasicTurbulenceModel >.
Definition at line 251 of file kOmegaSSTBase.H.
|
protected |
|
protected |
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().
|
protected |
|
protected |
|
inlineprotected |
Definition at line 172 of file kOmegaSSTBase.H.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::alphaK(), kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::alphaOmega(), kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::beta(), and kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().
|
inlineprotected |
Definition at line 181 of file kOmegaSSTBase.H.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::DkEff().
|
inlineprotected |
Definition at line 186 of file kOmegaSSTBase.H.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::DomegaEff().
|
inlineprotected |
Definition at line 191 of file kOmegaSSTBase.H.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().
|
inlineprotected |
Definition at line 196 of file kOmegaSSTBase.H.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().
|
protected |
Definition at line 115 of file kOmegaSSTBase.C.
References optionList::correct(), Foam::max(), options::New(), and Foam::sqrt().
|
protectedvirtual |
Reimplemented in kOmegaSSTSato< BasicTurbulenceModel >.
Definition at line 131 of file kOmegaSSTBase.C.
References kOmegaSST< TurbulenceModel, BasicTurbulenceModel >::epsilonByk(), Foam::fvc::grad(), Foam::magSqr(), and Foam::symm().
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().
|
protectedvirtual |
Return epsilon/k which for standard RAS is betaStar*omega.
Reimplemented in kOmegaSSTDES< BasicTurbulenceModel >.
Definition at line 139 of file kOmegaSSTBase.C.
Referenced by kOmegaSST< TurbulenceModel, BasicTurbulenceModel >::correctNut(), and kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().
|
protectedvirtual |
Definition at line 150 of file kOmegaSSTBase.C.
References Foam::dimTime, and Foam::dimVolume.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().
|
protectedvirtual |
Definition at line 165 of file kOmegaSSTBase.C.
References Foam::dimTime, Foam::dimVolume, and kOmegaSST< TurbulenceModel, BasicTurbulenceModel >::Qsas().
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().
|
protectedvirtual |
Reimplemented in kOmegaSSTSAS< BasicTurbulenceModel >.
Definition at line 180 of file kOmegaSSTBase.C.
References Foam::dimTime, and Foam::dimVolume.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma(), and kOmegaSST< TurbulenceModel, BasicTurbulenceModel >::omegaSource().
|
virtual |
Re-read model coefficients if they have changed.
Reimplemented in kOmegaSSTSato< BasicTurbulenceModel >, kOmegaSSTSAS< BasicTurbulenceModel >, and kOmegaSSTDES< BasicTurbulenceModel >.
Definition at line 377 of file kOmegaSSTBase.C.
References Foam::read().
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::~kOmegaSST().
|
inline |
Return the effective diffusivity for k.
Definition at line 261 of file kOmegaSSTBase.H.
|
inline |
Return the effective diffusivity for omega.
Definition at line 270 of file kOmegaSSTBase.H.
|
inlinevirtual |
Return the turbulence kinetic energy.
Definition at line 283 of file kOmegaSSTBase.H.
|
inlinevirtual |
Return the turbulence kinetic energy dissipation rate.
Definition at line 289 of file kOmegaSSTBase.H.
|
inlinevirtual |
Return the turbulence kinetic energy dissipation rate.
Definition at line 308 of file kOmegaSSTBase.H.
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Reimplemented in kOmegaSSTSato< BasicTurbulenceModel >.
Definition at line 405 of file kOmegaSSTBase.C.
References Foam::fvc::absolute(), beta(), Foam::bound(), fvMatrix< Type >::boundaryManipulate(), tmp< T >::clear(), optionList::constrain(), correct, optionList::correct(), Foam::fvm::ddt(), Foam::dev(), Foam::fvm::div(), Foam::fvc::div(), divU(), F1, fvOptions, Foam::constant::universal::G, Foam::fvc::grad(), Foam::fvm::laplacian(), Foam::magSqr(), Foam::max(), Foam::min(), options::New(), nut, phi, tmp< T >::ref(), fvMatrix< Type >::relax(), Foam::solve(), Foam::fvm::Sp(), Foam::sqrt(), Foam::fvm::SuSp(), Foam::symm(), and Foam::twoSymm().
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::omega().
|
protected |
Definition at line 132 of file kOmegaSSTBase.H.
|
protected |
Definition at line 133 of file kOmegaSSTBase.H.
|
protected |
Definition at line 135 of file kOmegaSSTBase.H.
|
protected |
Definition at line 136 of file kOmegaSSTBase.H.
|
protected |
Definition at line 138 of file kOmegaSSTBase.H.
|
protected |
Definition at line 139 of file kOmegaSSTBase.H.
|
protected |
Definition at line 141 of file kOmegaSSTBase.H.
|
protected |
Definition at line 142 of file kOmegaSSTBase.H.
|
protected |
Definition at line 144 of file kOmegaSSTBase.H.
|
protected |
Definition at line 146 of file kOmegaSSTBase.H.
|
protected |
Definition at line 147 of file kOmegaSSTBase.H.
|
protected |
Definition at line 148 of file kOmegaSSTBase.H.
|
protected |
Definition at line 150 of file kOmegaSSTBase.H.
|
protected |
Wall distance.
Note: different to wall distance in parent RASModel which is for near-wall cells only
Definition at line 158 of file kOmegaSSTBase.H.
|
protected |
Definition at line 160 of file kOmegaSSTBase.H.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::k().
|
protected |
Definition at line 161 of file kOmegaSSTBase.H.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::omega().