Go to the documentation of this file.
66 #ifndef realizableKE_H
67 #define realizableKE_H
83 template<
class BasicTurbulenceModel>
130 typedef typename BasicTurbulenceModel::alphaField
alphaField;
131 typedef typename BasicTurbulenceModel::rhoField
rhoField;
132 typedef typename BasicTurbulenceModel::transportModel
transportModel;
171 (this->
nut_/sigmak_ + this->
nu())
184 (this->
nut_/sigmaEps_ + this->
nu())
209 this->runTime_.timeName(),
Defines the attributes of an object for which implicit objectRegistry management is supported,...
virtual tmp< fvScalarMatrix > epsilonSource() const
A class for handling words, derived from Foam::string.
BasicTurbulenceModel::rhoField rhoField
A class for managing temporary objects.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual bool read()
Re-read model coefficients if they have changed.
static const word propertiesName
Default name of the turbulence properties dictionary.
tmp< volScalarField > rCmu(const volTensorField &gradU, const volScalarField &S2, const volScalarField &magS)
dimensionedScalar sigmaEps_
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
TypeName("realizableKE")
Runtime type information.
virtual tmp< fvScalarMatrix > kSource() const
Realizable k-epsilon turbulence model for incompressible and compressible flows.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
virtual tmp< volScalarField > omega() const
Return the (estimated) specific dissipation rate.
BasicTurbulenceModel::rhoField rhoField
virtual ~realizableKE()=default
Destructor.
BasicTurbulenceModel::alphaField alphaField
virtual void correctNut()
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::transportModel transportModel
realizableKE(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.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
BasicTurbulenceModel::transportModel transportModel
Eddy viscosity turbulence model base class.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
dimensionedScalar sigmak_