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 |
Public Types inherited from TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, 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... | |
kOmegaSST (const kOmegaSST &)=delete | |
Disallow default bitwise copy construction. 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... | |
void | operator= (const kOmegaSST &)=delete |
Disallow default bitwise assignment. More... | |
Public Member Functions inherited from TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel > | |
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... | |
TurbulenceModel (const TurbulenceModel &)=delete | |
Disallow default bitwise copy construction. 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... | |
virtual tmp< volScalarField > | nu () const |
Return the laminar viscosity. More... | |
virtual tmp< scalarField > | nu (const label patchi) const |
Return the laminar viscosity on patchi. More... | |
void | operator= (const TurbulenceModel &)=delete |
Disallow default bitwise assignment. 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_ |
Protected Attributes inherited from TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel > | |
const alphaField & | alpha_ |
const transportModel & | transport_ |
Additional Inherited Members | |
Static Public Member Functions inherited from TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel > | |
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 112 of file kOmegaSSTBase.H.
typedef BasicTurbulenceModel::alphaField alphaField |
Definition at line 242 of file kOmegaSSTBase.H.
typedef BasicTurbulenceModel::rhoField rhoField |
Definition at line 243 of file kOmegaSSTBase.H.
typedef BasicTurbulenceModel::transportModel transportModel |
Definition at line 244 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 213 of file kOmegaSSTBase.C.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::Qsas().
|
delete |
Disallow default bitwise copy construction.
|
inlinevirtual |
Destructor.
Reimplemented in kOmegaSST< BasicTurbulenceModel >.
Definition at line 267 of file kOmegaSSTBase.H.
|
protectedvirtual |
Reimplemented in kOmegaSSTLM< BasicTurbulenceModel >.
|
protectedvirtual |
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().
|
protectedvirtual |
|
protectedvirtual |
|
inlineprotected |
Definition at line 162 of file kOmegaSSTBase.H.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::alphaK(), kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::alphaOmega(), kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::beta(), kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::blend(), and kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().
|
inlineprotected |
Definition at line 172 of file kOmegaSSTBase.H.
|
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 192 of file kOmegaSSTBase.H.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::alphaOmega(), and kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().
|
inlineprotected |
Definition at line 200 of file kOmegaSSTBase.H.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::beta(), and kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().
|
protectedvirtual |
Reimplemented in kOmegaSSTSato< BasicTurbulenceModel >.
Definition at line 115 of file kOmegaSSTBase.C.
|
protectedvirtual |
Definition at line 131 of file kOmegaSSTBase.C.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().
|
protectedvirtual |
Return k production rate.
Definition at line 140 of file kOmegaSSTBase.C.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::correctNut(), and kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().
|
protectedvirtual |
Return epsilon/k which for standard RAS is betaStar*omega.
Definition at line 151 of file kOmegaSSTBase.C.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma(), and kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::Pk().
|
protectedvirtual |
Definition at line 162 of file kOmegaSSTBase.C.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().
|
protectedvirtual |
Definition at line 177 of file kOmegaSSTBase.C.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma().
|
protectedvirtual |
Definition at line 192 of file kOmegaSSTBase.C.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::gamma(), and kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::omegaSource().
|
virtual |
Re-read model coefficients if they have changed.
Reimplemented in kOmegaSSTLM< BasicTurbulenceModel >, kOmegaSSTSato< BasicTurbulenceModel >, kOmegaSSTSAS< BasicTurbulenceModel >, and kOmegaSSTDES< BasicTurbulenceModel >.
Definition at line 389 of file kOmegaSSTBase.C.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::~kOmegaSST().
|
inline |
Return the effective diffusivity for k.
Definition at line 277 of file kOmegaSSTBase.H.
|
inline |
Return the effective diffusivity for omega.
Definition at line 287 of file kOmegaSSTBase.H.
|
inlinevirtual |
Return the turbulence kinetic energy.
Definition at line 297 of file kOmegaSSTBase.H.
|
inlinevirtual |
Return the turbulence kinetic energy dissipation rate.
Definition at line 303 of file kOmegaSSTBase.H.
|
inlinevirtual |
Return the turbulence kinetic energy dissipation rate.
Definition at line 314 of file kOmegaSSTBase.H.
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Reimplemented in kOmegaSSTLM< BasicTurbulenceModel >, and kOmegaSSTSato< BasicTurbulenceModel >.
Definition at line 417 of file kOmegaSSTBase.C.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::omega().
|
delete |
Disallow default bitwise assignment.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::omega().
|
protected |
Definition at line 122 of file kOmegaSSTBase.H.
|
protected |
Definition at line 123 of file kOmegaSSTBase.H.
|
protected |
Definition at line 125 of file kOmegaSSTBase.H.
|
protected |
Definition at line 126 of file kOmegaSSTBase.H.
|
protected |
Definition at line 128 of file kOmegaSSTBase.H.
|
protected |
Definition at line 129 of file kOmegaSSTBase.H.
|
protected |
Definition at line 131 of file kOmegaSSTBase.H.
|
protected |
Definition at line 132 of file kOmegaSSTBase.H.
|
protected |
Definition at line 134 of file kOmegaSSTBase.H.
|
protected |
Definition at line 136 of file kOmegaSSTBase.H.
|
protected |
Definition at line 137 of file kOmegaSSTBase.H.
|
protected |
Definition at line 138 of file kOmegaSSTBase.H.
|
protected |
Definition at line 140 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 148 of file kOmegaSSTBase.H.
|
protected |
Definition at line 150 of file kOmegaSSTBase.H.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::k().
|
protected |
Definition at line 151 of file kOmegaSSTBase.H.
Referenced by kOmegaSST< LESeddyViscosity< BasicTurbulenceModel >, BasicTurbulenceModel >::omega().