32#include "phasePairKey.H"
38namespace diameterModels
59 numberWeighted_(
dict.getOrDefault<
Switch>(
"numberWeighted", false)),
104 if (
fluid.phasePairs().found(pairKeys_[
k]))
108 forAll(popBal_.velocityGroups(), j)
119 /(numberWeighted_ ? fi.
x() : fi.
d());
139 popBal_.fluid().phasePairs()[pairKeys_[
k]];
149 const scalar iDmdtSign =
152 const sizeGroup& fi = popBal_.sizeGroups()[i];
159 if (!numberWeighted_)
161 dDriftRate.
ref() *= fi.
x()/fi.
d();
164 driftRate += dDriftRate;
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const Mesh & mesh() const
Return mesh.
T & first() noexcept
The first element of the list, position [0].
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
static word timeName(const scalar t, const int precision=precision_)
const phaseModel & phase() const
Return the phase.
Base class for drift models.
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
Drift induced by interfacial phase change. By default phase change mass flux is distributed between s...
virtual void correct()
Correct diameter independent expressions.
virtual void addToDriftRate(volScalarField &driftRate, const label i)
Add to driftRate.
Class that solves the univariate population balance equation by means of a class method (also called ...
const phaseSystem & fluid() const
Return reference to the phaseSystem.
const fvMesh & mesh() const
Return reference to the mesh.
This class represents a single sizeGroup belonging to a velocityGroup. The main property of a sizeGro...
const dimensionedScalar & d() const
Return representative diameter of the sizeGroup.
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
const phaseModel & phase() const
Return const-reference to the phase.
This diameterModel is intended for use with a populationBalanceModel in order to simulate polydispers...
const PtrList< sizeGroup > & sizeGroups() const
Return sizeGroups belonging to this velocityGroup.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
const Time & time() const
Return the top-level database.
const dimensionedScalar & rho() const
const word & name() const
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
virtual word name() const
Pair name.
bool contains(const phaseModel &phase) const
Return true if this phasePair contains the given phase.
Class to represent a system of phases and model interfacial transfers between them.
const phasePairTable & phasePairs() const
Return the phase pairs.
Lookup type of boundary radiation properties.
A class for managing temporary objects.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
static constexpr const zero Zero
Global zero (0)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
const dimensionSet dimVolume(pow3(dimLength))
#define forAll(list, i)
Loop across all elements in list.