66 #ifndef mixtureKEpsilon_H
67 #define mixtureKEpsilon_H
83 template<
class BasicMomentumTransportModel>
188 typedef typename BasicMomentumTransportModel::alphaField
alphaField;
189 typedef typename BasicMomentumTransportModel::rhoField
rhoField;
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,.
Mixture k-epsilon turbulence model for two-phase gas-liquid systems.
BasicMomentumTransportModel::alphaField alphaField
virtual tmp< fvScalarMatrix > epsilonSource() const
dimensionedScalar sigmak_
tmp< volScalarField > Ct2() const
tmp< volScalarField > rhogEff() const
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
tmp< volScalarField > mixU(const volScalarField &fc, const volScalarField &fd) const
dimensionedScalar sigmaEps_
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
void operator=(const mixtureKEpsilon &)=delete
Disallow default bitwise assignment.
tmp< volScalarField > rhom() const
tmp< volScalarField > bubbleG() const
autoPtr< volScalarField > rhom_
dimensionedScalar alphap_
autoPtr< volScalarField > km_
autoPtr< volScalarField > Ct2_
virtual void correctNut()
Correct the eddy-viscosity nut.
tmp< volScalarField > DepsilonEff(const volScalarField &nutm) const
Return the effective diffusivity for epsilon.
virtual tmp< fvScalarMatrix > kSource() const
tmp< surfaceScalarField > mixFlux(const surfaceScalarField &fc, const surfaceScalarField &fd) const
tmp< volScalarField > DkEff(const volScalarField &nutm) const
Return the effective diffusivity for k.
tmp< volScalarField > rholEff() const
autoPtr< volScalarField > epsilonm_
tmp< volScalarField > mix(const volScalarField &fc, const volScalarField &fd) const
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
TypeName("mixtureKEpsilon")
Runtime type information.
virtual ~mixtureKEpsilon()
Destructor.
void boundEpsilonm(const volScalarField &Cc2)
Bound epsilonm.
virtual bool read()
Re-read model coefficients if they have changed.
mixtureKEpsilon(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
BasicMomentumTransportModel::rhoField rhoField
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Eddy viscosity turbulence model base class.
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
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))
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.