populationBalanceModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2017-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "populationBalanceModel.H"
27 #include "phaseSystem.H"
29 #include "coalescenceModel.H"
30 #include "breakupModel.H"
31 #include "binaryBreakupModel.H"
32 #include "fvmDdt.H"
33 #include "fvmDiv.H"
34 #include "fvmSup.H"
35 #include "fvcDdt.H"
36 #include "fvcDiv.H"
37 #include "shapeModel.H"
38 #include "distribution.H"
39 #include "fvSpecificSource.H"
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 namespace diameterModels
47 {
49 }
50 }
51 
52 
53 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
54 
55 const Foam::dictionary&
56 Foam::diameterModels::populationBalanceModel::coeffDict() const
57 {
58  return fluid_.subDict("populationBalanceCoeffs").subDict(name_);
59 }
60 
61 
62 void
63 Foam::diameterModels::populationBalanceModel::precomputeCoalescenceAndBreakup()
64 {
65  forAll(coalescenceModels_, model)
66  {
67  coalescenceModels_[model].precompute();
68  }
69 
70  forAll(breakupModels_, model)
71  {
72  breakupModels_[model].precompute();
73 
74  breakupModels_[model].dsdPtr()->precompute();
75  }
76 
77  forAll(binaryBreakupModels_, model)
78  {
79  binaryBreakupModels_[model].precompute();
80  }
81 }
82 
83 
84 void Foam::diameterModels::populationBalanceModel::birthByCoalescence
85 (
86  const label j,
87  const label k
88 )
89 {
90  const sizeGroup& fj = sizeGroups()[j];
91  const sizeGroup& fk = sizeGroups()[k];
92 
93  const dimensionedScalar v = fj.x() + fk.x();
94 
95  for (label i = j; i < sizeGroups().size(); i++)
96  {
97  const dimensionedScalar Eta = eta(i, v);
98 
99  if (Eta.value() == 0) continue;
100 
101  const sizeGroup& fi = sizeGroups()[i];
102 
103  tmp<volScalarField::Internal> tSui;
104  if (j == k)
105  {
106  tSui =
107  0.5*fi.x()/(fj.x()*fk.x())*Eta
108  *coalescenceRate_()*fj*fj.phase()*fk*fk.phase();
109  }
110  else
111  {
112  tSui =
113  fi.x()/(fj.x()*fk.x())*Eta
114  *coalescenceRate_()*fj*fj.phase()*fk*fk.phase();
115  }
116  const volScalarField::Internal& Sui = tSui();
117 
118  Su_[i] += Sui;
119 
120  const phaseInterface interfaceij(fi.phase(), fj.phase());
121 
122  if (dmdtfs_.found(interfaceij))
123  {
124  const scalar dmdtSign =
125  interfaceij.index(fi.phase()) == 0 ? +1 : -1;
126 
127  *dmdtfs_[interfaceij] += dmdtSign*fj.x()/v*Sui*fj.phase().rho();
128  }
129 
130  const phaseInterface interfaceik(fi.phase(), fk.phase());
131 
132  if (dmdtfs_.found(interfaceik))
133  {
134  const scalar dmdtSign =
135  interfaceik.index(fi.phase()) == 0 ? +1 : -1;
136 
137  *dmdtfs_[interfaceik] += dmdtSign*fk.x()/v*Sui*fk.phase().rho();
138  }
139 
140  sizeGroups_[i].shape().addCoalescence(Sui, fj, fk);
141  }
142 }
143 
144 
145 void Foam::diameterModels::populationBalanceModel::deathByCoalescence
146 (
147  const label i,
148  const label j
149 )
150 {
151  const sizeGroup& fi = sizeGroups()[i];
152  const sizeGroup& fj = sizeGroups()[j];
153 
154  Sp_[i] -= coalescenceRate_()*fi.phase()*fj*fj.phase()/fj.x();
155 
156  if (i != j)
157  {
158  Sp_[j] -= coalescenceRate_()*fj.phase()*fi*fi.phase()/fi.x();
159  }
160 }
161 
162 
163 void Foam::diameterModels::populationBalanceModel::birthByBreakup
164 (
165  const label k,
166  const label model
167 )
168 {
169  const sizeGroup& fk = sizeGroups()[k];
170 
171  for (label i = 0; i <= k; i++)
172  {
173  const sizeGroup& fi = sizeGroups()[i];
174 
175  tmp<volScalarField::Internal> tSui =
176  fi.x()*breakupModels_[model].dsdPtr()().nik(i, k)/fk.x()
177  *breakupRate_()*fk*fk.phase();
178  const volScalarField::Internal& Sui = tSui();
179 
180  Su_[i] += Sui;
181 
182  const phaseInterface interface(fi.phase(), fk.phase());
183 
184  if (dmdtfs_.found(interface))
185  {
186  const scalar dmdtSign =
187  interface.index(fi.phase()) == 0 ? +1 : -1;
188 
189  *dmdtfs_[interface] += dmdtSign*Sui*fk.phase().rho();
190  }
191 
192  sizeGroups_[i].shape().addBreakup(Sui, fk);
193  }
194 }
195 
196 
197 void Foam::diameterModels::populationBalanceModel::deathByBreakup(const label i)
198 {
199  Sp_[i] -= breakupRate_()*sizeGroups()[i].phase();
200 }
201 
202 
203 void Foam::diameterModels::populationBalanceModel::birthByBinaryBreakup
204 (
205  const label i,
206  const label j
207 )
208 {
209  const sizeGroup& fi = sizeGroups()[i];
210  const sizeGroup& fj = sizeGroups()[j];
211 
212  const volScalarField::Internal Su(binaryBreakupRate_()*fj*fj.phase());
213 
214  {
215  tmp<volScalarField::Internal> tSui = fi.x()*delta_[i][j]/fj.x()*Su;
216  const volScalarField::Internal& Sui = tSui();
217 
218  Su_[i] += Sui;
219 
220  sizeGroups_[i].shape().addBreakup(Sui, fj);
221 
222  const phaseInterface interfaceij(fi.phase(), fj.phase());
223 
224  if (dmdtfs_.found(interfaceij))
225  {
226  const scalar dmdtSign =
227  interfaceij.index(fi.phase()) == 0 ? +1 : -1;
228 
229  *dmdtfs_[interfaceij] += dmdtSign*Sui*fj.phase().rho();
230  }
231  }
232 
233  const dimensionedScalar v = fj.x() - fi.x();
234 
235  for (label k = 0; k <= j; k ++)
236  {
237  const dimensionedScalar Eta = eta(k, v);
238 
239  if (Eta.value() == 0) continue;
240 
241  const sizeGroup& fk = sizeGroups()[k];
242 
243  tmp<volScalarField::Internal> tSuk = fk.x()*delta_[i][j]*Eta/fj.x()*Su;
244  const volScalarField::Internal& Suk = tSuk();
245 
246  Su_[k] += Suk;
247 
248  const phaseInterface interfacekj(fk.phase(), fj.phase());
249 
250  if (dmdtfs_.found(interfacekj))
251  {
252  const scalar dmdtSign =
253  interfacekj.index(fk.phase()) == 0 ? +1 : -1;
254 
255  *dmdtfs_[interfacekj] += dmdtSign*Suk*fj.phase().rho();
256  }
257 
258  sizeGroups_[k].shape().addBreakup(Suk, fj);
259  }
260 }
261 
262 
263 void Foam::diameterModels::populationBalanceModel::deathByBinaryBreakup
264 (
265  const label j,
266  const label i
267 )
268 {
269  Sp_[i] -= sizeGroups()[i].phase()*binaryBreakupRate_()*delta_[j][i];
270 }
271 
272 
273 void
274 Foam::diameterModels::populationBalanceModel::computeCoalescenceAndBreakup()
275 {
276  forAll(sizeGroups(), i)
277  {
278  sizeGroups_[i].shape().reset();
279  }
280 
281  forAll(sizeGroups(), i)
282  {
283  Su_[i] = Zero;
284  Sp_[i] = Zero;
285  }
286 
287  forAllIter(dmdtfTable, dmdtfs_, dmdtfIter)
288  {
289  *dmdtfIter() = Zero;
290  }
291 
292  forAll(coalescencePairs_, coalescencePairi)
293  {
294  label i = coalescencePairs_[coalescencePairi].first();
295  label j = coalescencePairs_[coalescencePairi].second();
296 
297  coalescenceRate_() = Zero;
298 
299  forAll(coalescenceModels_, model)
300  {
301  coalescenceModels_[model].addToCoalescenceRate
302  (
303  coalescenceRate_(),
304  i,
305  j
306  );
307  }
308 
309  birthByCoalescence(i, j);
310 
311  deathByCoalescence(i, j);
312  }
313 
314  forAll(sizeGroups(), i)
315  {
316  forAll(breakupModels_, model)
317  {
318  breakupModels_[model].setBreakupRate(breakupRate_(), i);
319 
320  birthByBreakup(i, model);
321 
322  deathByBreakup(i);
323  }
324  }
325 
326  forAll(binaryBreakupPairs_, binaryBreakupPairi)
327  {
328  label i = binaryBreakupPairs_[binaryBreakupPairi].first();
329  label j = binaryBreakupPairs_[binaryBreakupPairi].second();
330 
331  binaryBreakupRate_() = Zero;
332 
333  forAll(binaryBreakupModels_, model)
334  {
335  binaryBreakupModels_[model].addToBinaryBreakupRate
336  (
337  binaryBreakupRate_(),
338  j,
339  i
340  );
341  }
342 
343  birthByBinaryBreakup(j, i);
344 
345  deathByBinaryBreakup(j, i);
346  }
347 }
348 
349 
350 void Foam::diameterModels::populationBalanceModel::precomputeExpansion()
351 {
352  forAllConstIter(HashTable<const velocityGroup*>, velocityGroupPtrs_, iter)
353  {
354  const velocityGroup& velGrp = *iter();
355  const phaseModel& phase = velGrp.phase();
356  const volScalarField& alpha = phase;
357  const volScalarField& rho = phase.rho();
358 
359  expansionRates_.set
360  (
361  phase.index(),
362  - (fvc::ddt(alpha, rho)()() + fvc::div(phase.alphaRhoPhi())()())/rho()
363  + fvc::ddt(alpha)()() + fvc::div(phase.alphaPhi())()()
364  );
365  }
366 }
367 
368 
370 Foam::diameterModels::populationBalanceModel::expansionSus
371 (
372  const label i,
373  const UPtrList<const volScalarField>& flds
374 ) const
375 {
376  const sizeGroup& fi = sizeGroups()[i];
377 
378  auto fiFld = [&](const label deltai)
379  {
380  return
381  flds.empty()
382  ? tmp<volScalarField::Internal>(sizeGroups()[i + deltai])
383  : sizeGroups()[i + deltai]()*flds[i + deltai]();
384  };
385 
386  Pair<tmp<volScalarField::Internal>> tSus;
387 
388  if (i != 0)
389  {
390  const sizeGroup& fiMinus1 = sizeGroups()[i - 1];
391  const phaseModel& phaseMinus1 = fiMinus1.phase();
392 
393  tSus.first() =
394  posPart(expansionRates_[phaseMinus1.index()])
395  *fi.x()/(fi.x() - fiMinus1.x())
396  *fiFld(-1);
397  }
398 
399  if (i != sizeGroups().size() - 1)
400  {
401  const sizeGroup& fiPlus1 = sizeGroups()[i + 1];
402  const phaseModel& phasePlus1 = fiPlus1.phase();
403 
404  tSus.second() =
405  - negPart(expansionRates_[phasePlus1.index()])
406  *fi.x()/(fiPlus1.x() - fi.x())
407  *fiFld(+1);
408  }
409 
410  return tSus;
411 }
412 
413 
414 void Foam::diameterModels::populationBalanceModel::computeExpansion()
415 {
416  forAllIter(dmdtfTable, expansionDmdtfs_, dmdtfIter)
417  {
418  *dmdtfIter() = Zero;
419  }
420 
421  forAllConstIter(HashTable<const velocityGroup*>, velocityGroupPtrs_, iter)
422  {
423  const sizeGroup& fi0 = iter()->sizeGroups().last();
424 
425  if (fi0.i() == sizeGroups().size() - 1) continue;
426 
427  const sizeGroup& fi1 = sizeGroups()[fi0.i() + 1];
428 
429  Pair<tmp<volScalarField::Internal>> tSus0 = expansionSus(fi0.i());
430  Pair<tmp<volScalarField::Internal>> tSus1 = expansionSus(fi1.i());
431 
432  const phaseInterface interface01(fi0.phase(), fi1.phase());
433  const scalar sign = interface01.index(fi0.phase()) == 0 ? -1 : +1;
434 
435  *expansionDmdtfs_[interface01] +=
436  sign
437  *(
438  - tSus0.second()*fi0.phase().rho()
439  + tSus1.first()*fi1.phase().rho()
440  );
441  }
442 }
443 
444 
445 void Foam::diameterModels::populationBalanceModel::precomputeModelSources()
446 {
447  // nothing to do
448 }
449 
450 
452 Foam::diameterModels::populationBalanceModel::modelSourceRhoSus
453 (
454  const label i,
455  const UPtrList<const volScalarField>& flds
456 ) const
457 {
458  const sizeGroup& fi = sizeGroups()[i];
459  const velocityGroup& velGrp = fi.group();
460 
461  const volScalarField& fldi = flds.empty() ? fi : flds[fi.i()];
462 
463  Pair<tmp<volScalarField::Internal>> tRhoSus;
464 
465  forAll(tRhoSus, Sui)
466  {
467  const sizeGroup& fiPopBalEnd =
468  Sui == 0 ? sizeGroups().first() : sizeGroups().last();
469 
470  if (fi.i() == fiPopBalEnd.i()) continue;
471 
472  const sizeGroup& fiVelGrpEnd =
473  Sui == 0 ? velGrp.sizeGroups().first() : velGrp.sizeGroups().last();
474 
475  if (fi.i() != fiVelGrpEnd.i()) continue;
476 
477  const sizeGroup& fiOther = sizeGroups()[i + (Sui == 0 ? -1 : +1)];
478 
479  forAll(fluid_.fvModels(), modeli)
480  {
481  if (!isA<fvSpecificSource>(fluid_.fvModels()[modeli])) continue;
482 
483  const fvSpecificSource& source =
484  refCast<const fvSpecificSource>(fluid_.fvModels()[modeli]);
485 
486  if
487  (
488  source.addsSupToField
489  (
490  fiOther.phase().volScalarField::name()
491  )
492  && isA<growthFvScalarFieldSource>
493  (
494  fldi.sources()[source.name()]
495  )
496  )
497  {
498  const growthFvScalarFieldSource& growthSource =
499  refCast<const growthFvScalarFieldSource>
500  (
501  fldi.sources()[source.name()]
502  );
503 
504  const volScalarField::Internal S(source.S(fiOther.name()));
505 
506  Pair<tmp<volScalarField::Internal>> sourceCoeffs =
507  growthSource.sourceCoeffs(source);
508 
509  tRhoSus[Sui] =
510  Sui == 0
511  ? posPart(S)*sourceCoeffs.first()
512  : negPart(S)*sourceCoeffs.second();
513  }
514  }
515  }
516 
517  return tRhoSus;
518 }
519 
520 
521 void Foam::diameterModels::populationBalanceModel::computeModelSources()
522 {
523  forAllIter(dmdtfTable, modelSourceDmdtfs_, dmdtfIter)
524  {
525  *dmdtfIter() = Zero;
526  }
527 
528  forAllConstIter(HashTable<const velocityGroup*>, velocityGroupPtrs_, iter)
529  {
530  const sizeGroup& fi0 = iter()->sizeGroups().last();
531 
532  if (fi0.i() == sizeGroups().size() - 1) continue;
533 
534  const sizeGroup& fi1 = sizeGroups()[fi0.i() + 1];
535 
536  Pair<tmp<volScalarField::Internal>> tRhoSus0 =
537  modelSourceRhoSus(fi0.i());
538  Pair<tmp<volScalarField::Internal>> tRhoSus1 =
539  modelSourceRhoSus(fi1.i());
540 
541  const phaseInterface interface01(fi0.phase(), fi1.phase());
542  const scalar sign = interface01.index(fi0.phase()) == 0 ? -1 : +1;
543 
544  if (tRhoSus0.second().valid())
545  {
546  *modelSourceDmdtfs_[interface01] -= sign*tRhoSus0.second();
547  }
548 
549  if (tRhoSus1.first().valid())
550  {
551  *modelSourceDmdtfs_[interface01] -= sign*tRhoSus1.first();
552  }
553  }
554 }
555 
556 
557 void Foam::diameterModels::populationBalanceModel::computeDilatationErrors()
558 {
559  PtrList<volScalarField::Internal> modelSourceDmdts(fluid_.phases().size());
560  forAllConstIter(dmdtfTable, modelSourceDmdtfs_, dmdtfIter)
561  {
562  const phaseInterface interface(fluid_, dmdtfIter.key());
563 
564  addField(interface.phase1(), "dmdt", *dmdtfIter(), modelSourceDmdts);
565  addField(interface.phase2(), "dmdt", - *dmdtfIter(), modelSourceDmdts);
566  }
567 
568  forAllConstIter(HashTable<const velocityGroup*>, velocityGroupPtrs_, iter)
569  {
570  const velocityGroup& velGrp = *iter();
571  const phaseModel& phase = velGrp.phase();
572  const volScalarField& alpha = phase;
573  const volScalarField& rho = phase.rho();
574 
575  dilatationErrors_.set
576  (
577  phase.index(),
578  fvc::ddt(alpha)()() + fvc::div(phase.alphaPhi())()()
579  - (fluid_.fvModels().source(alpha, rho) & rho)()()/rho()
580  );
581 
582  forAll(velGrp.sizeGroups(), i)
583  {
584  const sizeGroup& fi = velGrp.sizeGroups()[i];
585 
586  dilatationErrors_[phase.index()] -=
587  Su_[fi.i()] + expansionSu(fi.i())
588  + (Sp_[fi.i()] + expansionSp(fi.i()))*fi;
589  }
590 
591  if (modelSourceDmdts.set(phase.index()))
592  {
593  dilatationErrors_[phase.index()] -=
594  modelSourceDmdts[phase.index()]/rho;
595  }
596  }
597 }
598 
599 
600 bool Foam::diameterModels::populationBalanceModel::updateSources()
601 {
602  const bool result = sourceUpdateCounter_ % sourceUpdateInterval() == 0;
603 
604  ++ sourceUpdateCounter_;
605 
606  return result;
607 }
608 
609 
611 Foam::diameterModels::populationBalanceModel::etaCoeffs0(const label i) const
612 {
613  static const dimensionedScalar z(dimless, scalar(0));
614 
615  const dimensionedScalar& xi = sizeGroups()[i].x();
616 
617  if (i == 0) return Pair<dimensionedScalar>(z, 1/xi);
618 
619  const dimensionedScalar& x0 = sizeGroups()[i - 1].x();
620 
621  return Pair<dimensionedScalar>(- x0/(xi - x0), 1/(xi - x0));
622 }
623 
624 
626 Foam::diameterModels::populationBalanceModel::etaCoeffs1(const label i) const
627 {
628  static const dimensionedScalar z(dimless, scalar(0));
629 
630  const label n = sizeGroups().size();
631 
632  const dimensionedScalar& xi = sizeGroups()[i].x();
633 
634  if (i == n - 1) return Pair<dimensionedScalar>(z, 1/xi);
635 
636  const dimensionedScalar& x1 = sizeGroups()[i + 1].x();
637 
638  return Pair<dimensionedScalar>(x1/(x1 - xi), - 1/(x1 - xi));
639 }
640 
641 
643 Foam::diameterModels::populationBalanceModel::etaVCoeffs0(const label i) const
644 {
645  static const dimensionedScalar o(scalar(1)), zV(dimVolume, scalar(0));
646 
647  if (i == 0) return Pair<dimensionedScalar>(o, zV);
648 
649  const dimensionedScalar& x0 = 1/sizeGroups()[i - 1].x();
650  const dimensionedScalar& xi = 1/sizeGroups()[i].x();
651 
652  return Pair<dimensionedScalar>(- x0/(xi - x0), 1/(xi - x0));
653 }
654 
655 
657 Foam::diameterModels::populationBalanceModel::etaVCoeffs1(const label i) const
658 {
659  static const dimensionedScalar o(scalar(1)), zV(dimVolume, scalar(0));
660 
661  const label n = sizeGroups().size();
662 
663  if (i == n - 1) return Pair<dimensionedScalar>(o, zV);
664 
665  const dimensionedScalar& xi = 1/sizeGroups()[i].x();
666  const dimensionedScalar& x1 = 1/sizeGroups()[i + 1].x();
667 
668  return Pair<dimensionedScalar>(x1/(x1 - xi), - 1/(x1 - xi));
669 }
670 
671 
672 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
673 
675 (
676  const phaseSystem& fluid,
677  const word& name
678 )
679 :
681  (
682  IOobject
683  (
684  name,
685  fluid.time().constant(),
686  fluid.mesh()
687  )
688  ),
689  fluid_(fluid),
690  mesh_(fluid_.mesh()),
691  name_(name),
692  continuousPhase_
693  (
694  mesh_.lookupObject<phaseModel>
695  (
696  IOobject::groupName("alpha", coeffDict().lookup("continuousPhase"))
697  )
698  ),
699  sizeGroups_(),
700  v_(),
701  delta_(),
702  Su_(),
703  Sp_(),
704  dmdtfs_(),
705  expansionDmdtfs_(),
706  modelSourceDmdtfs_(),
707  expansionRates_(fluid_.phases().size()),
708  dilatationErrors_(fluid_.phases().size()),
709  coalescenceModels_
710  (
711  coeffDict().lookup("coalescenceModels"),
712  coalescenceModel::iNew(*this)
713  ),
714  coalescenceRate_(),
715  coalescencePairs_(),
716  breakupModels_
717  (
718  coeffDict().lookup("breakupModels"),
719  breakupModel::iNew(*this)
720  ),
721  breakupRate_(),
722  binaryBreakupModels_
723  (
724  coeffDict().lookup("binaryBreakupModels"),
725  binaryBreakupModel::iNew(*this)
726  ),
727  binaryBreakupRate_(),
728  binaryBreakupPairs_(),
729  alphas_(),
730  dsm_(),
731  U_(),
732  sourceUpdateCounter_(0)
733 {
734  groups::retrieve(*this, velocityGroupPtrs_, sizeGroups_);
735 
736  if (sizeGroups().size() < 3)
737  {
739  << "The populationBalance " << name_
740  << " requires a minimum number of three sizeGroups to be specified."
741  << exit(FatalError);
742  }
743 
744  // Create size-group boundaries
745  v_.setSize(sizeGroups().size() + 1);
746  v_.set(0, new dimensionedScalar("v", sizeGroups()[0].x()));
747  for (label i = 1; i < sizeGroups().size(); ++ i)
748  {
749  v_.set
750  (
751  i,
753  (
754  "v",
755  (sizeGroups()[i-1].x() + sizeGroups()[i].x())/2
756  )
757  );
758  }
759  v_.set(v_.size() - 1, new dimensionedScalar("v", sizeGroups().last().x()));
760 
761  // Create section widths if needed
762  if (binaryBreakupModels_.size() != 0)
763  {
764  delta_.setSize(sizeGroups().size());
765 
766  forAll(sizeGroups(), i)
767  {
768  delta_.set(i, new PtrList<dimensionedScalar>());
769  }
770 
771  forAll(sizeGroups(), i)
772  {
773  if (!delta_[i].empty()) continue;
774 
775  for (label j = 0; j <= sizeGroups().size() - 1; j++)
776  {
777  const sizeGroup& fj = sizeGroups()[j];
778 
779  delta_[i].append
780  (
782  (
783  "delta",
784  dimVolume,
785  v_[i+1].value() - v_[i].value()
786  )
787  );
788 
789  if
790  (
791  v_[i].value() < 0.5*fj.x().value()
792  && 0.5*fj.x().value() < v_[i+1].value()
793  )
794  {
795  delta_[i][j] = mag(0.5*fj.x() - v_[i]);
796  }
797  else if (0.5*fj.x().value() <= v_[i].value())
798  {
799  delta_[i][j].value() = 0;
800  }
801  }
802  }
803  }
804 
805  // Create size-group source terms
806  Su_.setSize(sizeGroups().size());
807  Sp_.setSize(sizeGroups().size());
808  forAll(sizeGroups(), i)
809  {
810  Su_.set
811  (
812  i,
814  (
815  IOobject
816  (
817  IOobject::groupName("Su" + Foam::name(i), this->name()),
818  fluid_.time().name(),
819  mesh_
820  ),
821  mesh_,
823  )
824  );
825 
826  Sp_.set
827  (
828  i,
830  (
831  IOobject
832  (
833  IOobject::groupName("Sp" + Foam::name(i), this->name()),
834  fluid_.time().name(),
835  mesh_
836  ),
837  mesh_,
839  )
840  );
841  }
842 
843  // Create interfacial mass transfer rates
845  (
847  velocityGroupPtrs_,
848  iter1
849  )
850  {
851  const diameterModels::velocityGroup& velGrp1 = *iter1();
852 
854  (
856  velocityGroupPtrs_,
857  iter2
858  )
859  {
860  const diameterModels::velocityGroup& velGrp2 = *iter2();
861 
862  const phaseInterface interface(velGrp1.phase(), velGrp2.phase());
863 
864  if (&velGrp1 != &velGrp2 && !dmdtfs_.found(interface))
865  {
866  dmdtfs_.insert
867  (
868  interface,
870  (
871  IOobject
872  (
874  (
875  typedName("dmdtf"),
876  interface.name()
877  ),
878  mesh().time().name(),
879  mesh()
880  ),
881  mesh(),
883  )
884  );
885 
886  expansionDmdtfs_.insert
887  (
888  interface,
890  (
891  IOobject
892  (
894  (
895  typedName("expansionDmdtf"),
896  interface.name()
897  ),
898  mesh().time().name(),
899  mesh()
900  ),
901  mesh(),
903  )
904  );
905 
906  modelSourceDmdtfs_.insert
907  (
908  interface,
910  (
911  IOobject
912  (
914  (
915  typedName("modelSourceDmdtf"),
916  interface.name()
917  ),
918  mesh().time().name(),
919  mesh()
920  ),
921  mesh(),
923  )
924  );
925  }
926  }
927  }
928 
929  if (coalescenceModels_.size() != 0)
930  {
931  coalescenceRate_.set
932  (
934  (
935  IOobject
936  (
937  IOobject::groupName("coalescenceRate", this->name()),
938  mesh_.time().name(),
939  mesh_
940  ),
941  mesh_,
943  )
944  );
945 
946  forAll(sizeGroups(), i)
947  {
948  for (label j = 0; j <= i; j++)
949  {
950  coalescencePairs_.append(labelPair(i, j));
951  }
952  }
953  }
954 
955  if (breakupModels_.size() != 0)
956  {
957  breakupRate_.set
958  (
960  (
961  IOobject
962  (
963  IOobject::groupName("breakupRate", this->name()),
964  fluid_.time().name(),
965  mesh_
966  ),
967  mesh_,
969  )
970  );
971  }
972 
973  if (binaryBreakupModels_.size() != 0)
974  {
975  binaryBreakupRate_.set
976  (
977  new volScalarField
978  (
979  IOobject
980  (
981  IOobject::groupName("binaryBreakupRate", this->name()),
982  fluid_.time().name(),
983  mesh_
984  ),
985  mesh_,
987  (
988  "binaryBreakupRate",
990  Zero
991  )
992  )
993  );
994 
995  forAll(sizeGroups(), i)
996  {
997  label j = 0;
998 
999  while (delta_[j][i].value() != 0)
1000  {
1001  binaryBreakupPairs_.append(labelPair(i, j));
1002  j++;
1003  }
1004  }
1005  }
1006 
1007  correct();
1008 }
1009 
1010 
1011 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1012 
1014 {}
1015 
1016 
1017 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1018 
1021 {
1022  notImplemented("populationBalance::clone() const");
1023  return autoPtr<populationBalanceModel>(nullptr);
1024 }
1025 
1026 
1028 {
1029  return os.good();
1030 }
1031 
1032 
1034 (
1035  const label i,
1036  const dimensionedScalar& v
1037 ) const
1038 {
1039  const Pair<dimensionedScalar> coeffs0 = etaCoeffs0(i);
1040  const Pair<dimensionedScalar> coeffs1 = etaCoeffs1(i);
1041 
1042  return
1043  max
1044  (
1045  min
1046  (
1047  coeffs0.first() + coeffs0.second()*v,
1048  coeffs1.first() + coeffs1.second()*v
1049  ),
1050  dimensionedScalar(dimless, scalar(0))
1051  );
1052 }
1053 
1054 
1057 (
1058  const label i,
1059  const volScalarField::Internal& v
1060 ) const
1061 {
1062  const Pair<dimensionedScalar> coeffs0 = etaCoeffs0(i);
1063  const Pair<dimensionedScalar> coeffs1 = etaCoeffs1(i);
1064 
1065  return
1066  max
1067  (
1068  min
1069  (
1070  coeffs0.first() + coeffs0.second()*v,
1071  coeffs1.first() + coeffs1.second()*v
1072  ),
1073  dimensionedScalar(dimless, scalar(0))
1074  );
1075 }
1076 
1077 
1079 (
1080  const label i,
1081  const dimensionedScalar& v
1082 ) const
1083 {
1084  const Pair<dimensionedScalar> coeffs0 = etaVCoeffs0(i);
1085  const Pair<dimensionedScalar> coeffs1 = etaVCoeffs1(i);
1086 
1087  return
1088  max
1089  (
1090  min
1091  (
1092  coeffs0.first() + coeffs0.second()/v,
1093  coeffs1.first() + coeffs1.second()/v
1094  ),
1095  dimensionedScalar(dimless, scalar(0))
1096  );
1097 }
1098 
1099 
1102 (
1103  const label i,
1104  const volScalarField::Internal& v
1105 ) const
1106 {
1107  const Pair<dimensionedScalar> coeffs0 = etaVCoeffs0(i);
1108  const Pair<dimensionedScalar> coeffs1 = etaVCoeffs1(i);
1109 
1110  return
1111  max
1112  (
1113  min
1114  (
1115  coeffs0.first() + coeffs0.second()/v,
1116  coeffs1.first() + coeffs1.second()/v
1117  ),
1118  dimensionedScalar(dimless, scalar(0))
1119  );
1120 }
1121 
1122 
1124 (
1125  const labelPair is,
1126  const distribution& d
1127 ) const
1128 {
1129  // Check that the distribution is generating volumes
1130  if (d.sampleQ() != 3)
1131  {
1133  << "The volumetric allocation coefficient should be evaluated "
1134  << "with a volumetrically sampled distribution (i.e., sampleQ "
1135  << "should equal 3)"
1136  << exit(FatalError);
1137  }
1138 
1139  // Get the four diameters that bound the sections of the size-group range's
1140  // basis function. Diameters #1 and #2 are the representative diameters of
1141  // the first and last size groups in the range. Diameters #0 and #3 are the
1142  // diameters of the adjacent size groups, or if at an end (or ends), the
1143  // upper or lower bounding diameter (or diameters) of the distribution.
1144  //
1145  // ^ .
1146  // | o - o - - - - o .
1147  // F_first | . . .\ .
1148  // | . . . \ .
1149  // +----+---+---------+--o------------------------------+---->
1150  // 0 1 2 3 .
1151  // ^ . .
1152  // | . o - - - - o .
1153  // F_middle | . /. .\ .
1154  // | . / . . \ .
1155  // +----+-------------o--+---------+--o-----------------+---->
1156  // . 0 1 2 3 .
1157  // ^ . .
1158  // | . o - - - - o - - - o
1159  // F_last | . /. . .
1160  // | . / . . .
1161  // +----+--------------------------o--+---------+-------+---->
1162  // . 0 1 2 3
1163  // . .
1164  // dMin dMax
1165  //
1166  const bool isFirst = is.first() == 0;
1167  const bool isLast = is.second() == sizeGroups().size() - 1;
1168  const scalarField dSphs
1169  (
1170  scalarList
1171  ({
1172  isFirst
1173  ? d.min()*(1 - small)
1174  : sizeGroups()[is.first() - 1].dSph().value(),
1175  sizeGroups()[is.first()].dSph().value(),
1176  sizeGroups()[is.second()].dSph().value(),
1177  isLast
1178  ? d.max()*(1 + small)
1179  : sizeGroups()[is.second() + 1].dSph().value()
1180  })
1181  );
1182 
1183  // Integrate the distribution, and the distribution divided by volume,
1184  // across the size-group range's basis function
1185  const scalarField integralPDFs
1186  (
1187  d.integralPDFxPow(dSphs, 0, true)
1188  );
1189  const scalarField integralPDFByVs
1190  (
1191  d.integralPDFxPow(dSphs, -3, true)
1193  );
1194 
1195  // Compute the integral of the size-group range's basis function multiplied
1196  // by the PDF.
1197  //
1198  // Between diameters #1 and #2 the basis function is a constant value of
1199  // one, so the integral is just the same as that of the distribution.
1200  //
1201  // Between diameters #0 and #1, and between #2 and #3, the basis function
1202  // is given by 'C0 + C1/v', where C0 and C1 are thecoefficients given by
1203  // etaVCoeffs0 and etaVCoeffs1. So, the integral is 'Integral(C0*PDF +
1204  // C1/v*PDF)'. C0 and C1 are constants, so this can be rearranged to
1205  // 'C0*Integral(PDF) + C1*Integral(PDF/v)'. The integrals in this
1206  // expression are those calculated above.
1207  const Pair<dimensionedScalar> etaVCoeffs0 = this->etaVCoeffs0(is.first());
1208  const Pair<dimensionedScalar> etaVCoeffs1 = this->etaVCoeffs1(is.second());
1209  return
1210  etaVCoeffs0.first().value()
1211  *(integralPDFs[1] - integralPDFs[0])
1212  + etaVCoeffs0.second().value()
1213  *(integralPDFByVs[1] - integralPDFByVs[0])
1214  + integralPDFs[2] - integralPDFs[1]
1215  + etaVCoeffs1.first().value()
1216  *(integralPDFs[3] - integralPDFs[2])
1217  + etaVCoeffs1.second().value()
1218  *(integralPDFByVs[3] - integralPDFByVs[2]);
1219 }
1220 
1221 
1223 (
1224  const label i,
1225  const distribution& d
1226 ) const
1227 {
1228  const sizeGroup& fi = sizeGroups()[i];
1229 
1230  const PtrList<sizeGroup>& vgSizeGroups = fi.group().sizeGroups();
1231 
1232  // Compute the integral for the size-group's basis function and for the
1233  // velocity-group's basis function. The allocation coefficient for this
1234  // size-group is then the ratio of these two integrals. The size-group
1235  // integral is "scaled up" to be a proportion for the phase, rather than
1236  // for the population balance as a whole.
1237 
1238  const dimensionedScalar sgIntegralPDFetaV =
1239  etaV(labelPair(i, i), d);
1240  const dimensionedScalar vgIntegralPDFetaV =
1241  etaV(labelPair(vgSizeGroups.first().i(), vgSizeGroups.last().i()), d);
1242 
1243  const dimensionedScalar sgSmall(dimless, rootVSmall/vgSizeGroups.size());
1244  const dimensionedScalar vgSmall(dimless, rootVSmall);
1245 
1246  return max(sgIntegralPDFetaV, sgSmall)/max(vgIntegralPDFetaV, vgSmall);
1247 }
1248 
1249 
1252 (
1253  const phaseModel& dispersedPhase
1254 ) const
1255 {
1256  return phaseInterface(dispersedPhase, continuousPhase_).sigma();
1257 }
1258 
1259 
1262 {
1263  return
1265  (
1266  continuousPhase_.name()
1267  );
1268 }
1269 
1270 
1273 {
1274  return Sp_[i];
1275 }
1276 
1277 
1280 (
1281  const label i,
1282  const UPtrList<const volScalarField>& flds
1283 ) const
1284 {
1285  Pair<tmp<volScalarField::Internal>> tSus = expansionSus(i, flds);
1286 
1287  return
1288  !tSus.first().valid() ? tSus.second()
1289  : !tSus.second().valid() ? tSus.first()
1290  : tSus.first() + tSus.second();
1291 }
1292 
1293 
1296 {
1297  const sizeGroup& fi = sizeGroups()[i];
1298  const phaseModel& phase = fi.phase();
1299 
1301 
1302  if (i == 0)
1303  {
1304  tSp = negPart(expansionRates_[phase.index()]);
1305  }
1306  else
1307  {
1308  const sizeGroup& fiMinus1 = sizeGroups()[i - 1];
1309 
1310  tSp =
1311  negPart(expansionRates_[phase.index()])
1312  *fi.x()/(fi.x() - fiMinus1.x());
1313  }
1314 
1315  if (i != sizeGroups().size() - 1)
1316  {
1317  const sizeGroup& fiPlus1 = sizeGroups()[i + 1];
1318 
1319  tSp.ref() -=
1320  posPart(expansionRates_[phase.index()])
1321  *fi.x()/(fiPlus1.x() - fi.x());
1322  }
1323  else
1324  {
1325  tSp.ref() += posPart(expansionRates_[phase.index()]);
1326  }
1327 
1328  return tSp;
1329 }
1330 
1331 
1334 (
1335  const label i,
1336  const UPtrList<const volScalarField>& flds
1337 ) const
1338 {
1339  Pair<tmp<volScalarField::Internal>> tRhoSus = modelSourceRhoSus(i, flds);
1340 
1341  const volScalarField::Internal rho = sizeGroups()[i].phase().rho();
1342 
1343  const dimensionedScalar zeroSu
1344  (
1345  (flds.empty() ? dimless : flds[i].dimensions())/dimTime,
1346  scalar(0)
1347  );
1348 
1349  return
1350  tRhoSus.first().valid() && tRhoSus.second().valid()
1351  ? (tRhoSus.first() + tRhoSus.second())/rho
1352  : tRhoSus.first().valid() ? tRhoSus.first()/rho
1353  : tRhoSus.second().valid() ? tRhoSus.second()/rho
1354  : volScalarField::Internal::New(zeroSu.name(), mesh(), zeroSu);
1355 }
1356 
1357 
1359 {
1360  if (!solveOnFinalIterOnly() || fluid_.pimple().finalIter())
1361  {
1362  const label nCorr =
1363  solverDict().lookupBackwardsCompatible<label>
1364  (
1365  {"nCorrectors", "nCorr"}
1366  );
1367 
1368  const scalar tolerance = solverDict().lookup<scalar>("tolerance");
1369 
1370  const bool updateSrc = updateSources();
1371 
1372  if (nCorr > 0 && updateSrc)
1373  {
1374  precomputeCoalescenceAndBreakup();
1375  }
1376  precomputeExpansion();
1377  precomputeModelSources();
1378 
1379  int iCorr = 0;
1380  scalar maxInitialResidual = 1;
1381  while (++iCorr <= nCorr && maxInitialResidual > tolerance)
1382  {
1383  Info<< "populationBalance " << this->name()
1384  << ": Iteration " << iCorr << endl;
1385 
1386  if (updateSrc)
1387  {
1388  computeCoalescenceAndBreakup();
1389  }
1390  computeExpansion();
1391  computeModelSources();
1392 
1393  computeDilatationErrors();
1394 
1395  maxInitialResidual = 0;
1396 
1397  forAll(sizeGroups(), i)
1398  {
1399  sizeGroup& fi = sizeGroups_[i];
1400  const phaseModel& phase = fi.phase();
1401  const volScalarField& alpha = phase;
1402  const volScalarField& rho = phase.rho();
1403 
1404  fvScalarMatrix fiEqn
1405  (
1406  fvm::ddt(alpha, fi)
1407  + fvm::div(phase.alphaPhi(), fi)
1408  + fvm::Sp(-(1 - small)*dilatationErrors_[phase.index()], fi)
1409  + fvm::SuSp(-small*dilatationErrors_[phase.index()], fi)
1410  ==
1411  Su_[i] + fvm::Sp(Sp_[i], fi)
1412  + expansionSu(i) + fvm::Sp(expansionSp(i), fi)
1413  + modelSourceSu(i)
1414  + fluid_.fvModels().source(alpha, rho, fi)/rho
1415  - correction
1416  (
1417  fvm::Sp
1418  (
1419  max(phase.residualAlpha() - alpha, scalar(0))
1420  /mesh().time().deltaT(),
1421  fi
1422  )
1423  )
1424  );
1425 
1426  fiEqn.relax();
1427 
1428  fluid_.fvConstraints().constrain(fiEqn);
1429 
1430  maxInitialResidual = max
1431  (
1432  fiEqn.solve().initialResidual(),
1433  maxInitialResidual
1434  );
1435 
1436  fluid_.fvConstraints().constrain(fi);
1437  }
1438  }
1439 
1440  const volScalarField alphaF0
1441  (
1442  sizeGroups().first().phase()*sizeGroups().first()
1443  );
1444 
1445  const volScalarField alphaFNm1
1446  (
1447  sizeGroups().last().phase()*sizeGroups().last()
1448  );
1449 
1450  Info<< this->name() << " sizeGroup phase fraction first, last = "
1451  << alphaF0.weightedAverage(mesh().V()).value()
1452  << ' ' << alphaFNm1.weightedAverage(mesh().V()).value()
1453  << endl;
1454  }
1455 }
1456 
1457 
1459 {
1460  if (velocityGroupPtrs_.size() <= 1) return;
1461 
1462  // Calculate the total void fraction
1463  if (alphas_.empty())
1464  {
1465  alphas_.set
1466  (
1467  new volScalarField
1468  (
1469  IOobject
1470  (
1471  IOobject::groupName("alpha", this->name()),
1472  fluid_.time().name(),
1473  mesh_,
1476  ),
1477  mesh_,
1479  )
1480  );
1481  }
1482  else
1483  {
1484  alphas_() = Zero;
1485  }
1486 
1487  forAllConstIter(HashTable<const velocityGroup*>, velocityGroupPtrs_, iter)
1488  {
1489  const phaseModel& phase = iter()->phase();
1490 
1491  alphas_() += max(phase, phase.residualAlpha());
1492  }
1493 
1494  // Calculate the Sauter mean diameter
1495  tmp<volScalarField> tInvDsm
1496  (
1498  (
1499  "invDsm",
1500  mesh_,
1502  )
1503  );
1504  volScalarField& invDsm = tInvDsm.ref();
1505 
1506  forAllConstIter(HashTable<const velocityGroup*>, velocityGroupPtrs_, iter)
1507  {
1508  const phaseModel& phase = iter()->phase();
1509 
1510  invDsm += max(phase, phase.residualAlpha())/(phase.d()*alphas_());
1511  }
1512 
1513  if (dsm_.empty())
1514  {
1515  dsm_.set
1516  (
1517  new volScalarField
1518  (
1519  IOobject
1520  (
1521  IOobject::groupName("d", this->name()),
1522  fluid_.time().name(),
1523  mesh_,
1526  ),
1527  1/tInvDsm
1528  )
1529  );
1530  }
1531  else
1532  {
1533  dsm_() = 1/tInvDsm;
1534  }
1535 
1536  // Calculate the average velocity
1537  if (U_.empty())
1538  {
1539  U_.set
1540  (
1541  new volVectorField
1542  (
1543  IOobject
1544  (
1545  IOobject::groupName("U", this->name()),
1546  fluid_.time().name(),
1547  mesh_,
1550  ),
1551  mesh_,
1553  )
1554  );
1555  }
1556  else
1557  {
1558  U_() = Zero;
1559  }
1560 
1561  forAllConstIter(HashTable<const velocityGroup*>, velocityGroupPtrs_, iter)
1562  {
1563  const phaseModel& phase = iter()->phase();
1564 
1565  U_() += max(phase, phase.residualAlpha())*phase.U()/alphas_();
1566  }
1567 }
1568 
1569 
1570 // ************************************************************************* //
label k
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:458
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh, PrimitiveField2 > &) const
Calculate and return weighted average.
static tmp< DimensionedField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< Type > &)
Return a temporary field constructed from name, mesh,.
Generic GeometricField class.
DimensionedField< Type, GeoMesh, PrimitiveField > Internal
Type of the internal field from which this GeometricField is derived.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
An STL-conforming hash table.
Definition: HashTable.H:127
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const Time & time() const
Return time.
Definition: IOobject.C:315
const word & name() const
Return name.
Definition: IOobject.H:307
static word groupName(Name name, const word &group)
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
An ordered pair of two objects of type <Type> with first() and second() elements.
Definition: Pair.H:66
const Type & second() const
Return second.
Definition: PairI.H:121
const Type & first() const
Return first.
Definition: PairI.H:107
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:66
T & first()
Return reference to the first element of the list.
Definition: UPtrListI.H:43
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
bool empty() const
Return true if the UPtrList is empty (ie, size() is zero)
Definition: UPtrListI.H:36
T & last()
Return reference to the last element of the list.
Definition: UPtrListI.H:57
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
const phaseModel & phase() const
Return the phase.
Base class for binary breakup models that provide a breakup rate between a size class pair directly,...
Base class for breakup models which provide a total breakup rate and a separate daughter size distrib...
Definition: breakupModel.H:56
Base class for coalescence models.
Constant dispersed-phase particle diameter model.
static void retrieve(const populationBalanceModel &popBal, HashTable< const velocityGroup * > &velGroupPtrs, UPtrList< sizeGroup > &szGroupPtrs)
Retrieve the pointers.
Return a pointer to a new populationBalanceModel object created on.
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
tmp< volScalarField::Internal > Sp(const label i) const
Return the implicit coalescence and breakup source term.
bool writeData(Ostream &) const
Dummy write for regIOobject.
autoPtr< populationBalanceModel > clone() const
Return clone.
populationBalanceModel(const phaseSystem &fluid, const word &name)
Construct for a fluid.
const UPtrList< sizeGroup > & sizeGroups() const
Return the size groups belonging to this populationBalance.
tmp< volScalarField::Internal > modelSourceSu(const label i, const UPtrList< const volScalarField > &flds=UPtrList< const volScalarField >()) const
Return the explicit model source source term.
tmp< volScalarField::Internal > expansionSp(const label i) const
Return the implicit expansion source term.
dimensionedScalar eta(const label i, const dimensionedScalar &v) const
Return the number allocation coefficient for a single volume.
tmp< volScalarField::Internal > expansionSu(const label i, const UPtrList< const volScalarField > &flds=UPtrList< const volScalarField >()) const
Return the explicit expansion source term.
dimensionedScalar etaV(const label i, const dimensionedScalar &v) const
Return the volume allocation coefficient for a single volume.
const tmp< volScalarField > sigmaWithContinuousPhase(const phaseModel &dispersedPhase) const
Return the surface tension coefficient between a given dispersed.
const fvMesh & mesh() const
Return reference to the mesh.
const phaseCompressibleMomentumTransportModel & continuousTurbulence() const
Return reference to momentumTransport model of the continuous phase.
void solve()
Solve the population balance equation.
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:96
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition: sizeGroupI.H:58
const velocityGroup & group() const
Return const-reference to the velocityGroup.
Definition: sizeGroupI.H:44
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:37
Computes the Sauter mean diameter based on a user specified size distribution, defined in terms of si...
Definition: velocityGroup.H:87
const PtrList< sizeGroup > & sizeGroups() const
Return sizeGroups belonging to this velocityGroup.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:849
const Type & value() const
Return const reference to value.
const word & name() const
Return const reference to name.
Base class for statistical distributions.
Definition: distribution.H:76
virtual scalar max() const =0
Return the maximum value.
virtual tmp< scalarField > integralPDFxPow(const scalarField &x, const label e, const bool consistent=false) const =0
Return the integral of the PDF multiplied by an integer power of x.
label sampleQ() const
Access the sample size exponent.
Definition: distribution.C:134
virtual scalar min() const =0
Return the minimum value.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:603
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
Templated abstract base class for multiphase compressible turbulence models.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
virtual word name() const
Name.
tmp< volScalarField > sigma() const
Surface tension coefficient.
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
Definition: phaseModel.C:169
tmp< volScalarField > d() const
Return the Sauter-mean diameter.
Definition: phaseModel.C:181
label index() const
Return the index of the phase.
Definition: phaseModel.C:157
virtual tmp< volVectorField > U() const =0
Return the velocity.
virtual const volScalarField & rho() const =0
Return the density field.
virtual tmp< surfaceScalarField > alphaPhi() const =0
Return the volumetric flux of the phase.
Class to represent a system of phases.
Definition: phaseSystem.H:74
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:55
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define notImplemented(functionName)
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:371
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for implicit and explicit sources.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
defineTypeNameAndDebug(constant, 0)
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
dimensionedScalar sign(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
void addField(const label phasei, const word &name, tmp< GeoField > field, PtrList< GeoField > &fieldList)
Definition: phaseSystemI.H:265
const dimensionSet dimless
messageStream Info
const dimensionSet dimLength
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
dimensionedScalar negPart(const dimensionedScalar &ds)
const dimensionSet dimDensity
const dimensionSet dimVolume
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:181
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimVelocity
error FatalError
void inv(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
dimensionedScalar posPart(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.