Standard high Reynolds-number k-omega turbulence model for incompressible and compressible flows. More...
#include <kOmega.H>


Public Types | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
| Public Types inherited from eddyViscosity< RASModel< BasicTurbulenceModel > > | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
| Public Types inherited from linearViscousStress< BasicTurbulenceModel > | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
Public Member Functions | |
| TypeName ("kOmega") | |
| Runtime type information. | |
| kOmega (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 | ~kOmega ()=default |
| Destructor. | |
| virtual bool | read () |
| Read RASProperties dictionary. | |
| tmp< volScalarField > | DkEff () const |
| Return the effective diffusivity for k. | |
| tmp< volScalarField > | DomegaEff () 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 void | correct () |
| Solve the turbulence equations and correct the turbulence viscosity. | |
| Public Member Functions inherited from eddyViscosity< RASModel< BasicTurbulenceModel > > | |
| eddyViscosity (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName) | |
| Construct from components. | |
| virtual | ~eddyViscosity ()=default |
| Destructor. | |
| virtual tmp< volScalarField > | nut () const |
| Return the turbulence viscosity. | |
| virtual tmp< volSymmTensorField > | R () const |
| Return the Reynolds stress tensor. | |
| virtual void | validate () |
| Validate the turbulence fields after construction. | |
| Public Member Functions inherited from linearViscousStress< BasicTurbulenceModel > | |
| linearViscousStress (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName) | |
| Construct from components. | |
| virtual | ~linearViscousStress ()=default |
| Destructor. | |
| virtual tmp< volSymmTensorField > | devRhoReff () const |
| Return the effective stress tensor. | |
| virtual tmp< volSymmTensorField > | devRhoReff (const volVectorField &U) const |
| Return the effective stress tensor based on a given velocity field. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
| Return the source term for the momentum equation. | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (const volScalarField &rho, volVectorField &U) const |
| Return the source term for the momentum equation. | |
Protected Member Functions | |
| virtual void | correctNut () |
Protected Attributes | |
| dimensionedScalar | Cmu_ |
| dimensionedScalar | beta_ |
| dimensionedScalar | gamma_ |
| dimensionedScalar | alphaK_ |
| dimensionedScalar | alphaOmega_ |
| volScalarField | k_ |
| volScalarField | omega_ |
| Protected Attributes inherited from eddyViscosity< RASModel< BasicTurbulenceModel > > | |
| volScalarField | nut_ |
Standard high Reynolds-number k-omega turbulence model for incompressible and compressible flows.
References:
Wilcox, D. C. (1988).
Reassessment of the scale-determining
equation for advanced turbulence models.
AIAA Journal, 26(11), 1299-1310.
DOI:10.2514/3.10041
The default model coefficients are
kOmegaCoeffs
{
Cmu 0.09; // Equivalent to betaStar
alpha 0.52;
beta 0.072;
alphak 0.5;
alphaOmega 0.5;
}
| typedef BasicTurbulenceModel::alphaField alphaField |
| typedef BasicTurbulenceModel::rhoField rhoField |
| typedef BasicTurbulenceModel::transportModel transportModel |
| kOmega | ( | 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 49 of file kOmega.C.
References alpha, alphaK_, alphaOmega_, beta_, Foam::bound(), Cmu_, eddyViscosity< RASModel< BasicTurbulenceModel > >::eddyViscosity(), gamma_, Foam::GlobalIOList< Tuple2< scalar, vector > >::typeName, k_, omega_, phi, rho, timeName, and U.

|
virtualdefault |
Destructor.
|
protectedvirtual |
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 36 of file kOmega.C.
References optionList::correct(), k_, options::New(), eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_, and omega_.
Referenced by correct().


| TypeName | ( | "kOmega< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Read RASProperties dictionary.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 157 of file kOmega.C.
References alphaK_, alphaOmega_, beta_, Cmu_, eddyViscosity< RASModel< BasicTurbulenceModel > >::eddyViscosity(), gamma_, and read().
Referenced by read().


|
inline |
|
inline |
|
inlinevirtual |
Return the turbulence kinetic energy.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
|
inlinevirtual |
|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 175 of file kOmega.C.
References Foam::fvc::absolute(), alpha, beta_, Foam::bound(), tmp< T >::clear(), Cmu_, correct(), correctNut(), Foam::fvm::ddt(), Foam::devTwoSymm(), Foam::fvc::div(), Foam::fvm::div(), divU, DkEff(), DomegaEff(), eddyViscosity< RASModel< BasicTurbulenceModel > >::eddyViscosity(), fvOptions, gamma_, Foam::fvc::grad(), k_, Foam::fvm::laplacian(), options::New(), eddyViscosity< RASModel< BasicTurbulenceModel > >::nut(), nut, eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_, omega_, phi, tmp< T >::ref(), rho, solve(), Foam::fvm::Sp(), Foam::fvm::SuSp(), and U.

|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
Definition at line 96 of file kOmega.H.
Referenced by correct(), correctNut(), and kOmega().
|
protected |
Definition at line 97 of file kOmega.H.
Referenced by correct(), correctNut(), and kOmega().