100 #ifndef kOmegaSSTBase_H
101 #define kOmegaSSTBase_H
112 template<
class MomentumTransportModel,
class BasicMomentumTransportModel>
115 public MomentumTransportModel
164 return F1*(psi1 - psi2) + psi2;
174 return F1*(psi1 - psi2) + psi2;
241 typedef typename BasicMomentumTransportModel::alphaField
alphaField;
242 typedef typename BasicMomentumTransportModel::rhoField
rhoField;
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Implementation of the k-omega-SST turbulence model for incompressible and compressible flows.
BasicMomentumTransportModel::alphaField alphaField
void boundOmega()
Bound omega.
tmp< volScalarField::Internal > gamma(const volScalarField::Internal &F1) const
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 > F1(const volScalarField &CDkOmega) const
virtual tmp< volScalarField > omega() const
Return the turbulence kinetic energy dissipation rate.
dimensionedScalar alphaOmega2_
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual tmp< volScalarField > F2() const
dimensionedScalar gamma2_
virtual tmp< volScalarField > F23() const
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField::Internal &F1, const volScalarField::Internal &F2) const
Return epsilon/k which for standard RAS is betaStar*omega.
kOmegaSST(const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct from components.
tmp< volScalarField > alphaOmega(const volScalarField &F1) const
dimensionedScalar alphaK2_
tmp< volScalarField > blend(const volScalarField &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
virtual ~kOmegaSST()
Destructor.
dimensionedScalar alphaK1_
virtual tmp< fvScalarMatrix > Qsas(const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
virtual tmp< fvScalarMatrix > omegaSource() const
virtual void correctNut()
void operator=(const kOmegaSST &)=delete
Disallow default bitwise assignment.
virtual tmp< fvScalarMatrix > kSource() const
dimensionedScalar alphaOmega1_
virtual tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Return k production rate.
tmp< volScalarField > DkEff(const volScalarField &F1) const
Return the effective diffusivity for k.
dimensionedScalar gamma1_
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
tmp< volScalarField > alphaK(const volScalarField &F1) const
dimensionedScalar betaStar_
virtual bool read()
Re-read model coefficients if they have changed.
tmp< volScalarField::Internal > beta(const volScalarField::Internal &F1) const
BasicMomentumTransportModel::rhoField rhoField
virtual tmp< volScalarField > F3() const
A class for managing temporary objects.
Abstract base class for all fluid physical properties.
A class for handling words, derived from string.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const dimensionedScalar G
Newtonian constant of gravitation.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.