Go to the documentation of this file.
31 #include "phaseCompressibleTurbulenceModel.H"
37 namespace diameterModels
39 namespace coalescenceModels
45 CoulaloglouTavlaridesCoalescence,
79 const phaseModel& continuousPhase = popBal_.continuousPhase();
80 const sizeGroup& fi = popBal_.sizeGroups()[i];
81 const sizeGroup& fj = popBal_.sizeGroups()[j];
84 C1_*(
pow(fi.
x(), 2.0/3.0) +
pow(fj.
x(), 2.0/3.0))
86 *
cbrt(popBal_.continuousTurbulence().epsilon())/(1 + popBal_.alphas())
89 - C2_*continuousPhase.
mu()*continuousPhase.
rho()
90 *popBal_.continuousTurbulence().epsilon()
91 /
sqr(popBal_.sigmaWithContinuousPhase(fi.
phase()))
92 /
pow3(1 + popBal_.alphas())
virtual tmp< volScalarField > mu() const
Base class for coalescence models.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
const phaseModel & phase() const
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
dimensionedScalar exp(const dimensionedScalar &ds)
CoulaloglouTavlaridesCoalescence(const populationBalanceModel &popBal, const dictionary &dict)
dimensionedScalar pow4(const dimensionedScalar &ds)
const dimensionSet dimArea(sqr(dimLength))
addToRunTimeSelectionTable(coalescenceModel, constantCoalescence, dictionary)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
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 ...
Macros for easy insertion into run-time selection tables.
Generic dimensioned Type class.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
defineTypeNameAndDebug(constantCoalescence, 0)
const dimensionedScalar & rho() const
dimensionedScalar cbrt(const dimensionedScalar &ds)
Generic GeometricField class.
const dimensionedScalar & x() const
const dimensionSet dimless