41 template<
class BasicMomentumTransportModel>
53 kEpsilon<BasicMomentumTransportModel>
64 gasTurbulencePtr_(nullptr),
66 alphaInversion_(
"alphaInversion", this->coeffDict(), 0.3),
67 Cp_(
"Cp", this->coeffDict(), 0.25),
68 C4_(
"C4", this->coeffDict(), this->C2_.value()),
69 Cmub_(
"Cmub", this->coeffDict(), 0.6)
75 template<
class BasicMomentumTransportModel>
80 alphaInversion_.readIfPresent(this->coeffDict());
81 Cp_.readIfPresent(this->coeffDict());
82 C4_.readIfPresent(this->coeffDict());
83 Cmub_.readIfPresent(this->coeffDict());
94 template<
class BasicMomentumTransportModel>
98 if (!gasTurbulencePtr_)
103 refCast<const phaseModel>(this->properties());
115 momentumTransportModel::typeName,
121 return *gasTurbulencePtr_;
125 template<
class BasicMomentumTransportModel>
129 this->gasTurbulence();
136 this->Cmu_*
sqr(this->k_)/this->epsilon_
137 + Cmub_*gas.
d()*gasTurbulence.
alpha()
138 *(
mag(this->U_ - gasTurbulence.
U()));
140 this->nut_.correctBoundaryConditions();
145 template<
class BasicMomentumTransportModel>
149 this->gasTurbulence();
177 template<
class BasicMomentumTransportModel>
189 max(alphaInversion_ -
alpha, scalar(0))
191 *
min(gasTurbulence.
epsilon()/gasTurbulence.
k(), 1.0/
U.time().deltaT())
196 template<
class BasicMomentumTransportModel>
203 this->gasTurbulence();
205 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
209 + phaseTransferCoeff*gasTurbulence.
k()
210 -
fvm::Sp(phaseTransferCoeff, this->k_);
214 template<
class BasicMomentumTransportModel>
222 this->gasTurbulence();
224 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
227 alpha*
rho*this->C4_*this->epsilon_*bubbleG()/this->k_
228 + phaseTransferCoeff*gasTurbulence.
epsilon()
229 -
fvm::Sp(phaseTransferCoeff, this->epsilon_);
233 template<
class BasicMomentumTransportModel>
Generic GeometricField class.
static word groupName(Name name, const word &group)
Continuous-phase k-epsilon model including bubble-generated turbulence.
virtual tmp< fvScalarMatrix > epsilonSource() const
Source term for the epsilon equation.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
LaheyKEpsilon(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.
tmp< volScalarField > bubbleG() const
tmp< volScalarField > phaseTransferCoeff() const
virtual void correctNut()
Correct the eddy-viscosity nut.
virtual tmp< fvScalarMatrix > kSource() const
Source term for the k equation.
virtual bool read()
Read model coefficients if they have changed.
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Class to represent a interface between phases where one phase is considered dispersed within the othe...
Base class for Lagrangian drag models.
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Abstract base class for turbulence models (RAS, LES and laminar).
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
const volVectorField & U() const
Access function to velocity field.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
Templated abstract base class for multiphase compressible turbulence models.
const alphaField & alpha() const
Access function to phase fraction.
tmp< volScalarField > d() const
Return the Sauter-mean diameter.
const word & name() const
Return the name of this phase.
Class to represent a system of phases.
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument in a two-phase system.
const ModelType & lookupInterfacialModel(const phaseInterface &interface) const
Return a sub model for an interface.
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))
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.