103 #ifndef kOmegaSSTBase_H 104 #define kOmegaSSTBase_H 115 template<
class TurbulenceModel,
class BasicTurbulenceModel>
179 return F1*(psi1 -
psi2) + psi2;
184 return blend(F1, alphaK1_, alphaK2_);
189 return blend(F1, alphaOmega1_, alphaOmega2_);
194 return blend(F1, beta1_, beta2_);
199 return blend(F1, gamma1_, gamma2_);
230 typedef typename BasicTurbulenceModel::alphaField
alphaField;
231 typedef typename BasicTurbulenceModel::rhoField
rhoField;
232 typedef typename BasicTurbulenceModel::transportModel
transportModel;
241 const alphaField&
alpha,
299 this->mesh_.time().timeName(),
tmp< volScalarField > alphaOmega(const volScalarField &F1) const
tmp< volScalarField > beta(const volScalarField &F1) const
virtual ~kOmegaSST()
Destructor.
const transportModel & transport() const
Access function to incompressible transport model.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
dimensionedScalar gamma2_
tmp< volScalarField > nu() const
Return the laminar viscosity.
dimensionedScalar betaStar_
tmp< volScalarField > alphaK(const volScalarField &F1) const
tmp< volScalarField > gamma(const volScalarField &F1) const
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual tmp< fvScalarMatrix > kSource() const
BasicTurbulenceModel::alphaField alphaField
Templated abstract base class for turbulence models.
BasicTurbulenceModel::rhoField rhoField
virtual tmp< fvScalarMatrix > omegaSource() const
dimensionedScalar gamma1_
tmp< volScalarField > DomegaEff(const volScalarField &F1) const
Return the effective diffusivity for omega.
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none.
dimensionedScalar alphaOmega2_
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Implementation of the k-omega-SST turbulence model for incompressible and compressible flows...
tmp< volScalarField > blend(const volScalarField &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
tmp< volScalarField > F1(const volScalarField &CDkOmega) const
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
virtual tmp< volScalarField > epsilonByk(const volScalarField &F1, const volScalarField &F2) const
Return epsilon/k which for standard RAS is betaStar*omega.
tmp< volScalarField > F2() const
virtual tmp< volScalarField > omega() const
Return the turbulence kinetic energy dissipation rate.
dimensionedScalar alphaOmega1_
dimensionedScalar alphaK2_
dimensionedScalar alphaK1_
const alphaField & alpha() const
Access function to phase fraction.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
virtual tmp< fvScalarMatrix > Qsas(const volScalarField &S2, const volScalarField &gamma, const volScalarField &beta) const
tmp< volScalarField > F23() const
A class for managing temporary objects.
const volScalarField & y_
Wall distance.
tmp< volScalarField > F3() const
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
virtual bool read()
Re-read model coefficients if they have changed.
BasicTurbulenceModel::transportModel transportModel
virtual void correctNut()
tmp< volScalarField > DkEff(const volScalarField &F1) const
Return the effective diffusivity for k.