Go to the documentation of this file.
32 #include "phaseCompressibleTurbulenceModel.H"
38 namespace diameterModels
40 namespace coalescenceModels
65 C1_(dimensionedScalar::getOrDefault(
"C1",
dict,
dimless, 0.356)),
68 dimensionedScalar::getOrDefault
78 dimensionedScalar::getOrDefault
86 turbulence_(
dict.lookup(
"turbulence")),
87 buoyancy_(
dict.lookup(
"buoyancy")),
88 laminarShear_(
dict.lookup(
"laminarShear"))
102 const phaseModel& continuousPhase = popBal_.continuousPhase();
103 const sizeGroup& fi = popBal_.sizeGroups()[i];
104 const sizeGroup& fj = popBal_.sizeGroups()[j];
116 pow3(rij)*continuousPhase.
rho()
117 /(16.0*popBal_.sigmaWithContinuousPhase(fi.
phase()))
120 *
cbrt(popBal_.continuousTurbulence().epsilon())/
pow(rij, 2.0/3.0)
129 *
cbrt(popBal_.continuousTurbulence().epsilon())
132 *collisionEfficiency;
146 2.14*popBal_.sigmaWithContinuousPhase(fi.
phase())
147 /(continuousPhase.
rho()*fi.
d()) + 0.505*
mag(
g)*fi.
d()
151 2.14*popBal_.sigmaWithContinuousPhase(fi.
phase())
152 /(continuousPhase.
rho()*fj.
d()) + 0.505*
mag(
g)*fj.
d()
156 *collisionEfficiency;
162 <<
"Laminar shear collision contribution not implemented for "
163 << this->
type() <<
" coalescence model."
Base class for coalescence models.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
const phaseModel & phase() const
Return const-reference to the phase.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
virtual tmp< volScalarField > rho() const =0
Return the density field.
dimensionedScalar exp(const dimensionedScalar &ds)
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
PrinceBlanch(const populationBalanceModel &popBal, const dictionary &dict)
addToRunTimeSelectionTable(coalescenceModel, constantCoalescence, dictionary)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
This class represents a single sizeGroup belonging to a velocityGroup. The main property of a sizeGro...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Class that solves the univariate population balance equation by means of a class method (also called ...
dimensionedScalar log(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
const uniformDimensionedVectorField & g
errorManipArg< error, int > exit(error &err, const int errNo=1)
const dimensionedScalar & d() const
Return representative diameter of the sizeGroup.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr scalar pi(M_PI)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const dimensionedScalar e
Elementary charge.
defineTypeNameAndDebug(constantCoalescence, 0)
dimensionedScalar cbrt(const dimensionedScalar &ds)