Go to the documentation of this file.
31 #include "phaseSystem.H"
32 #include "phasePairKey.H"
38 namespace diameterModels
58 pairKeys_(
dict.lookup(
"pairs")),
59 numberWeighted_(
dict.getOrDefault<
Switch>(
"numberWeighted",
false)),
75 IOobject::groupName(
type() +
":W", pair.
name()),
76 popBal_.mesh().time().timeName(),
104 if (
fluid.phasePairs().found(pairKeys_[
k]))
119 /(numberWeighted_ ? fi.
x() : fi.
d());
139 popBal_.fluid().phasePairs()[pairKeys_[
k]];
149 const scalar iDmdtSign =
150 vg.
phase().
name() == pair.first() ? +1 : -1;
152 const sizeGroup& fi = popBal_.sizeGroups()[i];
159 if (!numberWeighted_)
161 dDriftRate.
ref() *= fi.
x()/fi.
d();
164 driftRate += dDriftRate;
const phaseSystem & fluid() const
Return reference to the phaseSystem.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
const phaseModel & phase() const
Return const-reference to the phase.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
virtual void correct()
Correct diameter independent expressions.
#define forAll(list, i)
Loop across all elements in list.
const UPtrList< velocityGroup > & velocityGroups() const
Return the velocityGroups belonging to this populationBalance.
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
This diameterModel is intended for use with a populationBalanceModel in order to simulate polydispers...
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
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,...
bool contains(const phaseModel &phase) const
Return true if this phasePair contains the given phase.
Class that solves the univariate population balance equation by means of a class method (also called ...
defineTypeNameAndDebug(constantDrift, 0)
Macros for easy insertion into run-time selection tables.
const word & name() const
const dimensionedScalar & d() const
Return representative diameter of the sizeGroup.
virtual word name() const
Pair name.
label k
Boltzmann constant.
Class to represent a system of phases and model interfacial transfers between them.
const dimensionedScalar & rho() const
Base class for drift models.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
const PtrList< sizeGroup > & sizeGroups() const
Return sizeGroups belonging to this velocityGroup.
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
const dimensionSet dimVolume(pow3(dimLength))
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
phaseChange(const populationBalanceModel &popBal, const dictionary &dict)
Construct from a population balance model and a dictionary.
const phaseModel & phase() const
Return the phase.
addToRunTimeSelectionTable(driftModel, constantDrift, dictionary)
virtual void addToDriftRate(volScalarField &driftRate, const label i)
Add to driftRate.