30 template<
class BasePhaseSystem>
37 BasePhaseSystem(mesh),
40 this->
lookup(
"populationBalances"),
41 diameterModels::populationBalanceModel::iNew(*this, dmdtfs_)
44 forAll(populationBalances_, popBali)
46 const diameterModels::populationBalanceModel& popBal =
47 populationBalances_[popBali];
51 HashTable<const diameterModels::velocityGroup*>,
52 popBal.velocityGroupPtrs(),
56 const diameterModels::velocityGroup& velGrp1 = *iter1();
60 HashTable<const diameterModels::velocityGroup*>,
61 popBal.velocityGroupPtrs(),
65 const diameterModels::velocityGroup& velGrp2 = *iter2();
67 const phaseInterface interface
77 !dmdtfs_.
found(interface)
80 this->
template validateMassTransfer
81 <diameterModels::populationBalanceModel>(interface);
92 "populationBalance:dmdtf",
95 this->
mesh().time().timeName(),
111 template<
class BasePhaseSystem>
119 template<
class BasePhaseSystem>
123 const phaseInterfaceKey& key
126 tmp<volScalarField> tDmdtf = BasePhaseSystem::dmdtf(key);
128 if (dmdtfs_.
found(key))
130 tDmdtf.ref() += *dmdtfs_[key];
137 template<
class BasePhaseSystem>
141 PtrList<volScalarField>
dmdts(BasePhaseSystem::dmdts());
145 const phaseInterface interface(*
this, dmdtfIter.key());
147 addField(interface.phase1(),
"dmdt", *dmdtfIter(),
dmdts);
148 addField(interface.phase2(),
"dmdt", - *dmdtfIter(),
dmdts);
155 template<
class BasePhaseSystem>
159 autoPtr<phaseSystem::momentumTransferTable> eqnsPtr =
164 this->addDmdtUfs(dmdtfs_, eqns);
170 template<
class BasePhaseSystem>
174 autoPtr<phaseSystem::momentumTransferTable> eqnsPtr =
175 BasePhaseSystem::momentumTransferf();
179 this->addDmdtUfs(dmdtfs_, eqns);
185 template<
class BasePhaseSystem>
189 autoPtr<phaseSystem::heatTransferTable> eqnsPtr =
190 BasePhaseSystem::heatTransfer();
194 this->addDmdtHefs(dmdtfs_, eqns);
200 template<
class BasePhaseSystem>
204 autoPtr<phaseSystem::specieTransferTable> eqnsPtr =
215 template<
class BasePhaseSystem>
218 const PtrList<volScalarField>& rAUs,
219 const PtrList<surfaceScalarField>& rAUfs
224 forAll(populationBalances_, i)
226 populationBalances_[i].solve();
231 template<
class BasePhaseSystem>
236 forAll(populationBalances_, i)
238 populationBalances_[i].correct();
243 template<
class BasePhaseSystem>
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-10/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
HashPtrTable< fvVectorMatrix > momentumTransferTable
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
HashPtrTable< fvScalarMatrix > heatTransferTable
HashPtrTable< fvScalarMatrix > specieTransferTable
virtual ~PopulationBalancePhaseSystem()
Destructor.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimTime
stressControl lookup("compactNormalStress") >> compactNormalStress
bool found(const Key &) const
Return true if hashedEntry is found in table.
static word groupName(Name name, const word &group)
const dimensionSet dimDensity
virtual void solve(const PtrList< volScalarField > &rAUs, const PtrList< surfaceScalarField > &rAUfs)
Solve all population balance equations.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer()
Return the momentum transfer matrices for the cell-based algorithm.
void addDmdtYfs(const phaseSystem::dmdtfTable &dmdtfs, phaseSystem::specieTransferTable &eqns) const
Add specie transfer terms which result from bulk mass transfers.
virtual void correct()
Correct derived properties.
phaseSystem::specieTransferTable & specieTransfer(specieTransferPtr())
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
Return the momentum transfer matrices for the face-based algorithm.
HashPtrTable< volScalarField, phaseInterfaceKey, phaseInterfaceKey::hash > dmdtfTable
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
phaseSystem::momentumTransferTable & momentumTransfer(momentumTransferPtr())
PopulationBalancePhaseSystem(const fvMesh &)
Construct from fvMesh.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
virtual autoPtr< phaseSystem::specieTransferTable > specieTransfer() const
Return the specie transfer matrices.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
virtual tmp< volScalarField > dmdtf(const phaseInterfaceKey &key) const
Return the mass transfer rate for an interface.
A class for managing temporary objects.
virtual bool read()
Read base phaseProperties dictionary.
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.