42 template<
class BasicMomentumTransportModel>
54 kEpsilon<BasicMomentumTransportModel>
65 liquidTurbulencePtr_(nullptr),
71 this->groupName(
"nutEff"),
72 this->runTime_.
name(),
92 this->printCoeffs(
type);
99 template<
class BasicMomentumTransportModel>
104 alphaInversion_.readIfPresent(this->coeffDict());
115 template<
class BasicMomentumTransportModel>
122 const phaseModel& gas = refCast<const phaseModel>(this->properties());
141 exp(
min(thetal/thetag, scalar(50))),
147 nutEff_ =
omega*liquidTurbulence.
nut();
152 template<
class BasicMomentumTransportModel>
156 if (!liquidTurbulencePtr_)
160 const phaseModel& gas = refCast<const phaseModel>(this->properties());
164 liquidTurbulencePtr_ =
168 return *liquidTurbulencePtr_;
172 template<
class BasicMomentumTransportModel>
182 (this->alpha_ - scalar(0.5))/(alphaInversion_ - 0.5),
191 this->groupName(
"nuEff"),
193 + (1.0 - blend)*rhoEff()*nutEff_/this->rho_
199 template<
class BasicMomentumTransportModel>
203 const phaseModel& gas = refCast<const phaseModel>(this->properties());
214 this->groupName(
"rhoEff"),
220 template<
class BasicMomentumTransportModel>
232 max(alphaInversion_ -
alpha, scalar(0))
236 liquidTurbulence.
epsilon()/liquidTurbulence.
k(),
237 1.0/
U.time().deltaT()
243 template<
class BasicMomentumTransportModel>
248 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
251 phaseTransferCoeff*liquidTurbulence.
k()
252 -
fvm::Sp(phaseTransferCoeff, this->k_);
256 template<
class BasicMomentumTransportModel>
261 const volScalarField phaseTransferCoeff(this->phaseTransferCoeff());
264 phaseTransferCoeff*liquidTurbulence.
epsilon()
265 -
fvm::Sp(phaseTransferCoeff, this->epsilon_);
269 template<
class BasicMomentumTransportModel>
277 this->groupName(
"R"),
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,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
virtual tmp< fvScalarMatrix > epsilonSource() const
Source term for the epsilon equation.
virtual tmp< volSymmTensorField > sigma() const
Return the stress tensor [m^2/s^2].
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
tmp< volScalarField > phaseTransferCoeff() const
virtual tmp< volScalarField > rhoEff() const
Return the effective density for the stress.
virtual void correctNut()
Correct the eddy-viscosity nut.
virtual tmp< fvScalarMatrix > kSource() const
Source term for the k equation.
continuousGasKEpsilon(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.
const momentumTransportModel & liquidTurbulence() const
Return the turbulence model for the liquid phase.
virtual bool read()
Re-read model coefficients if they have changed.
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
virtual void correctNut()
Correct the eddy-viscosity nut.
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Generic dimensioned Type class.
Class to represent a interface between phases where one phase is considered dispersed within the othe...
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
virtual const word & name() const
Return the name of the liquid.
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
scalar rho(scalar p, scalar T) const
Liquid density [kg/m^3].
Abstract base class for turbulence models (RAS, LES and laminar).
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
Convenience class to handle the input of constant rotational speed. Reads an omega entry with default...
const phaseSystem & fluid() const
Return the system to which this phase belongs.
tmp< volScalarField > d() const
Return the Sauter-mean diameter.
virtual const volScalarField & rho() const =0
Return the density field.
Class to represent a system of phases and model interfacial transfers between them.
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.
virtual tmp< volScalarField > Cvm() const =0
Return the virtual mass coefficient.
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< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
word name(const bool)
Return a word representation of a bool.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
static const Identity< scalar > I
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.