Go to the documentation of this file.
33 template<
class BasePhaseSystem>
40 if (!pDmdt_.found(key))
42 return phaseSystem::dmdt(key);
47 return pDmdtSign**pDmdt_[key];
53 template<
class BasePhaseSystem>
60 BasePhaseSystem(
mesh),
64 this->
lookup(
"populationBalances"),
68 forAll(populationBalances_, i)
71 populationBalances_[i];
77 if (!this->phasePairs_.found(key))
79 this->phasePairs_.insert
86 this->phaseModels_[key.first()],
87 this->phaseModels_[key.
second()]
117 IOobject::groupName(
"pDmdt", pair.
name()),
118 this->
mesh().time().timeName(),
120 IOobject::READ_IF_PRESENT,
133 template<
class BasePhaseSystem>
141 template<
class BasePhaseSystem>
148 return BasePhaseSystem::dmdt(key) + this->pDmdt(key);
152 template<
class BasePhaseSystem>
160 const phasePair& pair = this->phasePairs_[pDmdtIter.key()];
163 this->addField(pair.
phase1(),
"dmdt", pDmdt, dmdts);
164 this->addField(pair.
phase2(),
"dmdt", - pDmdt, dmdts);
171 template<
class BasePhaseSystem>
210 IOobject::groupName(Yi[i].member(),
phase.
name())
215 IOobject::groupName(Yi[i].member(), otherPhase.
name())
219 dmdt21*eqns[otherName]->psi()
223 dmdt12*eqns[
name]->psi()
232 template<
class BasePhaseSystem>
250 template<
class BasePhaseSystem>
255 forAll(populationBalances_, i)
257 populationBalances_[i].solve();
const T & second() const noexcept
Return second element, which is also the last element.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
virtual autoPtr< phaseSystem::massTransferTable > massTransfer() const
Return the mass transfer matrices.
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 class for handling words, derived from Foam::string.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
A class for managing temporary objects.
const dimensionSet dimDensity
dimensionedScalar posPart(const dimensionedScalar &ds)
Return a pointer to a new populationBalanceModel object created on.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
#define forAll(list, i)
Loop across all elements in list.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
word name(const complex &c)
Return string representation of complex.
const cellModel * lookup(const word &modelName)
Deprecated(2017-11) equivalent to cellModel::ptr static method.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const
Return the mass transfer rate for a pair.
Class that solves the univariate population balance equation by means of a class method (also called ...
virtual void solve()
Solve all population balance equations.
Mesh data needed to do the Finite Volume discretisation.
bool ordered() const
Return the ordered flag.
const word & name() const
const word & name() const
Return the name of this phase.
dimensionedScalar negPart(const dimensionedScalar &ds)
An ordered pair of two objects of type <T> with first() and second() elements.
virtual tmp< volScalarField > pDmdt(const phasePairKey &key) const
Return the population balance mass transfer rate.
virtual word name() const
Pair name.
const phasePairTable & phasePairs() const
Return list of unordered phasePairs in this populationBalance.
const phaseModel & phase2() const
Return phase 2.
virtual bool read()
Read base phaseProperties dictionary.
PopulationBalancePhaseSystem(const fvMesh &)
Construct from fvMesh.
phaseSystem::massTransferTable & massTransfer(massTransferPtr())
const phaseModel & phase1() const
Return phase 1.
const volScalarField & psi
virtual ~PopulationBalancePhaseSystem()
Destructor.