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-2020 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 "coalescenceModel.H"
28 #include "breakupModel.H"
29 #include "binaryBreakupModel.H"
30 #include "driftModel.H"
31 #include "nucleationModel.H"
32 #include "phaseSystem.H"
33 #include "surfaceTensionModel.H"
34 #include "fvmDdt.H"
35 #include "fvcDdt.H"
36 #include "fvmSup.H"
37 #include "fvcSup.H"
38 #include "fvcDiv.H"
40 #include "shapeModel.H"
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 namespace diameterModels
47 {
48  defineTypeNameAndDebug(populationBalanceModel, 0);
49 }
50 }
51 
52 
53 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
54 
55 void Foam::diameterModels::populationBalanceModel::registerVelocityGroups()
56 {
57  forAll(fluid_.phases(), phasei)
58  {
59  if (isA<velocityGroup>(fluid_.phases()[phasei].dPtr()()))
60  {
61  const velocityGroup& velGroup =
62  refCast<const velocityGroup>(fluid_.phases()[phasei].dPtr()());
63 
64  if (velGroup.popBalName() == this->name())
65  {
66  velocityGroups_.resize(velocityGroups_.size() + 1);
67 
68  velocityGroups_.set
69  (
70  velocityGroups_.size() - 1,
71  &const_cast<velocityGroup&>(velGroup)
72  );
73 
74  forAll(velGroup.sizeGroups(), i)
75  {
76  this->registerSizeGroups
77  (
78  const_cast<sizeGroup&>(velGroup.sizeGroups()[i])
79  );
80  }
81  }
82  }
83  }
84 }
85 
86 
87 void Foam::diameterModels::populationBalanceModel::registerSizeGroups
88 (
89  sizeGroup& group
90 )
91 {
92  if
93  (
94  sizeGroups_.size() != 0
95  &&
96  group.x().value() <= sizeGroups_.last().x().value()
97  )
98  {
100  << "Size groups must be entered according to their representative"
101  << " size"
102  << exit(FatalError);
103  }
104 
105  sizeGroups_.resize(sizeGroups_.size() + 1);
106  sizeGroups_.set(sizeGroups_.size() - 1, &group);
107 
108  // Grid generation over property space
109  if (sizeGroups_.size() == 1)
110  {
111  // Set the first sizeGroup boundary to the representative value of
112  // the first sizeGroup
113  v_.append
114  (
116  (
117  "v",
118  sizeGroups_.last().x()
119  )
120  );
121 
122  // Set the last sizeGroup boundary to the representative size of the
123  // last sizeGroup
124  v_.append
125  (
127  (
128  "v",
129  sizeGroups_.last().x()
130  )
131  );
132  }
133  else
134  {
135  // Overwrite the next-to-last sizeGroup boundary to lie halfway between
136  // the last two sizeGroups
137  v_.last() =
138  0.5
139  *(
140  sizeGroups_[sizeGroups_.size()-2].x()
141  + sizeGroups_.last().x()
142  );
143 
144  // Set the last sizeGroup boundary to the representative size of the
145  // last sizeGroup
146  v_.append
147  (
149  (
150  "v",
151  sizeGroups_.last().x()
152  )
153  );
154  }
155 
156  delta_.append(new PtrList<dimensionedScalar>());
157 
158  Su_.append
159  (
160  new volScalarField
161  (
162  IOobject
163  (
164  "Su",
165  fluid_.time().timeName(),
166  mesh_
167  ),
168  mesh_,
170  )
171  );
172 
173  SuSp_.append
174  (
175  new volScalarField
176  (
177  IOobject
178  (
179  "SuSp",
180  fluid_.time().timeName(),
181  mesh_
182  ),
183  mesh_,
185  )
186  );
187 }
188 
189 
190 void Foam::diameterModels::populationBalanceModel::createPhasePairs()
191 {
192  forAll(velocityGroups_, i)
193  {
194  const phaseModel& phasei = velocityGroups_[i].phase();
195 
196  forAll(velocityGroups_, j)
197  {
198  const phaseModel& phasej = velocityGroups_[j].phase();
199 
200  if (&phasei != &phasej)
201  {
202  const phasePairKey key
203  (
204  phasei.name(),
205  phasej.name(),
206  false
207  );
208 
209  if (!phasePairs_.found(key))
210  {
211  phasePairs_.insert
212  (
213  key,
214  autoPtr<phasePair>
215  (
216  new phasePair
217  (
218  phasei,
219  phasej
220  )
221  )
222  );
223  }
224  }
225  }
226  }
227 }
228 
229 
230 void Foam::diameterModels::populationBalanceModel::correct()
231 {
232  calcDeltas();
233 
234  forAll(velocityGroups_, v)
235  {
236  velocityGroups_[v].preSolve();
237  }
238 
239  forAll(coalescence_, model)
240  {
241  coalescence_[model].correct();
242  }
243 
244  forAll(breakup_, model)
245  {
246  breakup_[model].correct();
247 
248  breakup_[model].dsdPtr()().correct();
249  }
250 
251  forAll(binaryBreakup_, model)
252  {
253  binaryBreakup_[model].correct();
254  }
255 
256  forAll(drift_, model)
257  {
258  drift_[model].correct();
259  }
260 
261  forAll(nucleation_, model)
262  {
263  nucleation_[model].correct();
264  }
265 }
266 
267 
268 void Foam::diameterModels::populationBalanceModel::
269 birthByCoalescence
270 (
271  const label j,
272  const label k
273 )
274 {
275  const sizeGroup& fj = sizeGroups_[j];
276  const sizeGroup& fk = sizeGroups_[k];
277 
278  dimensionedScalar Eta;
279  dimensionedScalar v = fj.x() + fk.x();
280 
281  for (label i = j; i < sizeGroups_.size(); i++)
282  {
283  // Calculate fraction for intra-interval events
284  Eta = eta(i, v);
285 
286  if (Eta.value() == 0) continue;
287 
288  const sizeGroup& fi = sizeGroups_[i];
289 
290  // Avoid double counting of events
291  if (j == k)
292  {
293  Sui_ =
294  0.5*fi.x()*coalescenceRate_()*Eta
295  *fj*fj.phase()/fj.x()
296  *fk*fk.phase()/fk.x();
297  }
298  else
299  {
300  Sui_ =
301  fi.x()*coalescenceRate_()*Eta
302  *fj*fj.phase()/fj.x()
303  *fk*fk.phase()/fk.x();
304  }
305 
306  Su_[i] += Sui_;
307 
308  const phasePairKey pairij
309  (
310  fi.phase().name(),
311  fj.phase().name()
312  );
313 
314  if (pDmdt_.found(pairij))
315  {
316  const scalar dmdtSign
317  (
318  Pair<word>::compare(pDmdt_.find(pairij).key(), pairij)
319  );
320 
321  pDmdt_[pairij]->ref() += dmdtSign*fj.x()/v*Sui_*fi.phase().rho();
322  }
323 
324  const phasePairKey pairik
325  (
326  fi.phase().name(),
327  fk.phase().name()
328  );
329 
330  if (pDmdt_.found(pairik))
331  {
332  const scalar dmdtSign
333  (
334  Pair<word>::compare(pDmdt_.find(pairik).key(), pairik)
335  );
336 
337  pDmdt_[pairik]->ref() += dmdtSign*fk.x()/v*Sui_*fi.phase().rho();
338  }
339 
340  sizeGroups_[i].shapeModelPtr()->addCoalescence(Sui_, fj, fk);
341  }
342 }
343 
344 
345 void Foam::diameterModels::populationBalanceModel::
346 deathByCoalescence
347 (
348  const label i,
349  const label j
350 )
351 {
352  const sizeGroup& fi = sizeGroups_[i];
353  const sizeGroup& fj = sizeGroups_[j];
354 
355  SuSp_[i] += coalescenceRate_()*fi.phase()*fj*fj.phase()/fj.x();
356 
357  if (i != j)
358  {
359  SuSp_[j] += coalescenceRate_()*fj.phase()*fi*fi.phase()/fi.x();
360  }
361 }
362 
363 
364 void Foam::diameterModels::populationBalanceModel::
365 birthByBreakup
366 (
367  const label k,
368  const label model
369 )
370 {
371  const sizeGroup& fk = sizeGroups_[k];
372 
373  for (label i = 0; i <= k; i++)
374  {
375  const sizeGroup& fi = sizeGroups_[i];
376 
377  Sui_ =
378  fi.x()*breakupRate_()*breakup_[model].dsdPtr()().nik(i, k)
379  *fk*fk.phase()/fk.x();
380 
381  Su_[i] += Sui_;
382 
383  const phasePairKey pair
384  (
385  fi.phase().name(),
386  fk.phase().name()
387  );
388 
389  if (pDmdt_.found(pair))
390  {
391  const scalar dmdtSign
392  (
393  Pair<word>::compare(pDmdt_.find(pair).key(), pair)
394  );
395 
396  pDmdt_[pair]->ref() += dmdtSign*Sui_*fi.phase().rho();
397  }
398 
399  sizeGroups_[i].shapeModelPtr()->addBreakup(Sui_, fk);
400  }
401 }
402 
403 
404 void Foam::diameterModels::populationBalanceModel::deathByBreakup(const label i)
405 {
406  const sizeGroup& fi = sizeGroups_[i];
407 
408  SuSp_[i] += breakupRate_()*fi.phase();
409 }
410 
411 
412 void Foam::diameterModels::populationBalanceModel::calcDeltas()
413 {
414  forAll(sizeGroups_, i)
415  {
416  if (delta_[i].empty())
417  {
418  for (label j = 0; j <= sizeGroups_.size() - 1; j++)
419  {
420  delta_[i].append
421  (
423  (
424  "delta",
425  dimVolume,
426  v_[i+1].value() - v_[i].value()
427  )
428  );
429 
430  if
431  (
432  v_[i].value() < 0.5*sizeGroups_[j].x().value()
433  &&
434  0.5*sizeGroups_[j].x().value() < v_[i+1].value()
435  )
436  {
437  delta_[i][j] = mag(0.5*sizeGroups_[j].x() - v_[i]);
438  }
439  else if (0.5*sizeGroups_[j].x().value() <= v_[i].value())
440  {
441  delta_[i][j].value() = 0;
442  }
443  }
444  }
445  }
446 }
447 
448 
449 void Foam::diameterModels::populationBalanceModel::
450 birthByBinaryBreakup
451 (
452  const label i,
453  const label j
454 )
455 {
456  const sizeGroup& fj = sizeGroups_[j];
457  const sizeGroup& fi = sizeGroups_[i];
458 
459  Sui_ = fi.x()*binaryBreakupRate_()*delta_[i][j]*fj*fj.phase()/fj.x();
460 
461  Su_[i] += Sui_;
462 
463  sizeGroups_[i].shapeModelPtr()->addBreakup(Sui_, fj);
464 
465  const phasePairKey pairij
466  (
467  fi.phase().name(),
468  fj.phase().name()
469  );
470 
471  if (pDmdt_.found(pairij))
472  {
473  const scalar dmdtSign
474  (
475  Pair<word>::compare(pDmdt_.find(pairij).key(), pairij)
476  );
477 
478  pDmdt_[pairij]->ref() += dmdtSign*Sui_*fi.phase().rho();
479  }
480 
481  dimensionedScalar Eta;
482  dimensionedScalar v = fj.x() - fi.x();
483 
484  for (label k = 0; k <= j; k++)
485  {
486  // Calculate fraction for intra-interval events
487  Eta = eta(k, v);
488 
489  if (Eta.value() == 0) continue;
490 
491  const sizeGroup& fk = sizeGroups_[k];
492 
493  volScalarField& Suk = Sui_;
494 
495  Suk =
496  sizeGroups_[k].x()*binaryBreakupRate_()*delta_[i][j]*Eta
497  *fj*fj.phase()/fj.x();
498 
499  Su_[k] += Suk;
500 
501  const phasePairKey pairkj
502  (
503  fk.phase().name(),
504  fj.phase().name()
505  );
506 
507  if (pDmdt_.found(pairkj))
508  {
509  const scalar dmdtSign
510  (
512  (
513  pDmdt_.find(pairkj).key(),
514  pairkj
515  )
516  );
517 
518  pDmdt_[pairkj]->ref() += dmdtSign*Suk*fi.phase().rho();
519  }
520 
521  sizeGroups_[i].shapeModelPtr()->addBreakup(Suk, fj);
522  }
523 }
524 
525 
526 void Foam::diameterModels::populationBalanceModel::
527 deathByBinaryBreakup
528 (
529  const label j,
530  const label i
531 )
532 {
533  const volScalarField& alphai = sizeGroups_[i].phase();
534 
535  SuSp_[i] += alphai*binaryBreakupRate_()*delta_[j][i];
536 }
537 
538 
539 void Foam::diameterModels::populationBalanceModel::drift
540 (
541  const label i,
542  driftModel& model
543 )
544 {
545  model.addToDriftRate(driftRate_(), i);
546 
547  const sizeGroup& fp = sizeGroups_[i];
548 
549  if (i == 0)
550  {
551  rx_() = pos(driftRate_())*sizeGroups_[i+1].x()/sizeGroups_[i].x();
552  }
553  else if (i == sizeGroups_.size() - 1)
554  {
555  rx_() = neg(driftRate_())*sizeGroups_[i-1].x()/sizeGroups_[i].x();
556  }
557  else
558  {
559  rx_() =
560  pos(driftRate_())*sizeGroups_[i+1].x()/sizeGroups_[i].x()
561  + neg(driftRate_())*sizeGroups_[i-1].x()/sizeGroups_[i].x();
562  }
563 
564  SuSp_[i] +=
565  (neg(1 - rx_()) + neg(rx_() - rx_()/(1 - rx_())))*driftRate_()
566  *fp.phase()/((rx_() - 1)*fp.x());
567 
568  rx_() = Zero;
569  rdx_() = Zero;
570 
571  if (i == sizeGroups_.size() - 2)
572  {
573  rx_() = pos(driftRate_())*sizeGroups_[i+1].x()/sizeGroups_[i].x();
574 
575  rdx_() =
576  pos(driftRate_())
577  *(sizeGroups_[i+1].x() - sizeGroups_[i].x())
578  /(sizeGroups_[i].x() - sizeGroups_[i-1].x());
579  }
580  else if (i < sizeGroups_.size() - 2)
581  {
582  rx_() = pos(driftRate_())*sizeGroups_[i+2].x()/sizeGroups_[i+1].x();
583 
584  rdx_() =
585  pos(driftRate_())
586  *(sizeGroups_[i+2].x() - sizeGroups_[i+1].x())
587  /(sizeGroups_[i+1].x() - sizeGroups_[i].x());
588  }
589 
590  if (i == 1)
591  {
592  rx_() += neg(driftRate_())*sizeGroups_[i-1].x()/sizeGroups_[i].x();
593 
594  rdx_() +=
595  neg(driftRate_())
596  *(sizeGroups_[i].x() - sizeGroups_[i-1].x())
597  /(sizeGroups_[i+1].x() - sizeGroups_[i].x());
598  }
599  else if (i > 1)
600  {
601  rx_() += neg(driftRate_())*sizeGroups_[i-2].x()/sizeGroups_[i-1].x();
602 
603  rdx_() +=
604  neg(driftRate_())
605  *(sizeGroups_[i-1].x() - sizeGroups_[i-2].x())
606  /(sizeGroups_[i].x() - sizeGroups_[i-1].x());
607  }
608 
609  if (i != sizeGroups_.size() - 1)
610  {
611  const sizeGroup& fe = sizeGroups_[i+1];
612  volScalarField& Sue = Sui_;
613 
614  Sue =
615  pos(driftRate_())*driftRate_()*rdx_()
616  *fp*fp.phase()/fp.x()
617  /(rx_() - 1);
618 
619  Su_[i+1] += Sue;
620 
621  const phasePairKey pairij
622  (
623  fp.phase().name(),
624  fe.phase().name()
625  );
626 
627  if (pDmdt_.found(pairij))
628  {
629  const scalar dmdtSign
630  (
631  Pair<word>::compare(pDmdt_.find(pairij).key(), pairij)
632  );
633 
634  pDmdt_[pairij]->ref() -= dmdtSign*Sue*fp.phase().rho();
635  }
636 
637  sizeGroups_[i+1].shapeModelPtr()->addDrift(Sue, fp, model);
638  }
639 
640  if (i != 0)
641  {
642  const sizeGroup& fw = sizeGroups_[i-1];
643  volScalarField& Suw = Sui_;
644 
645  Suw =
646  neg(driftRate_())*driftRate_()*rdx_()
647  *fp*fp.phase()/fp.x()
648  /(rx_() - 1);
649 
650  Su_[i-1] += Suw;
651 
652  const phasePairKey pairih
653  (
654  fp.phase().name(),
655  fw.phase().name()
656  );
657 
658  if (pDmdt_.found(pairih))
659  {
660  const scalar dmdtSign
661  (
662  Pair<word>::compare(pDmdt_.find(pairih).key(), pairih)
663  );
664 
665  pDmdt_[pairih]->ref() -= dmdtSign*Suw*fp.phase().rho();
666  }
667 
668  sizeGroups_[i-1].shapeModelPtr()->addDrift(Suw, fp, model);
669  }
670 }
671 
672 
673 void Foam::diameterModels::populationBalanceModel::
674 nucleation
675 (
676  const label i,
677  nucleationModel& model
678 )
679 {
680  const sizeGroup& fi = sizeGroups_[i];
681 
682  model.addToNucleationRate(nucleationRate_(), i);
683 
684  Sui_ = sizeGroups_[i].x()*nucleationRate_();
685 
686  Su_[i] += Sui_;
687 
688  sizeGroups_[i].shapeModelPtr()->addNucleation(Sui_, fi, model);
689 }
690 
691 
692 void Foam::diameterModels::populationBalanceModel::sources()
693 {
694  forAll(sizeGroups_, i)
695  {
696  sizeGroups_[i].shapeModelPtr()->reset();
697  Su_[i] = Zero;
698  SuSp_[i] = Zero;
699  }
700 
702  (
703  phasePairTable,
704  phasePairs(),
705  phasePairIter
706  )
707  {
708  pDmdt_(phasePairIter())->ref() = Zero;
709  }
710 
711  // Since the calculation of the rates is computationally expensive,
712  // they are calculated once for each sizeGroup pair and inserted into source
713  // terms as required
714  forAll(sizeGroups_, i)
715  {
716  if (coalescence_.size() != 0)
717  {
718  for (label j = 0; j <= i; j++)
719  {
720  coalescenceRate_() = Zero;
721 
722  forAll(coalescence_, model)
723  {
724  coalescence_[model].addToCoalescenceRate
725  (
726  coalescenceRate_(),
727  i,
728  j
729  );
730  }
731 
732  birthByCoalescence(i, j);
733 
734  deathByCoalescence(i, j);
735  }
736  }
737 
738  if (breakup_.size() != 0)
739  {
740  forAll(breakup_, model)
741  {
742  breakup_[model].setBreakupRate(breakupRate_(), i);
743 
744  birthByBreakup(i, model);
745 
746  deathByBreakup(i);
747  }
748  }
749 
750  if (binaryBreakup_.size() != 0)
751  {
752  label j = 0;
753 
754  while (delta_[j][i].value() != 0)
755  {
756  binaryBreakupRate_() = Zero;
757 
758  forAll(binaryBreakup_, model)
759  {
760  binaryBreakup_[model].addToBinaryBreakupRate
761  (
762  binaryBreakupRate_(),
763  j,
764  i
765  );
766  }
767 
768  birthByBinaryBreakup(j, i);
769 
770  deathByBinaryBreakup(j, i);
771 
772  j++;
773  }
774  }
775 
776  if (drift_.size() != 0)
777  {
778  forAll(drift_, model)
779  {
780  driftRate_() = Zero;
781  drift(i, drift_[model]);
782  }
783  }
784 
785  if (nucleation_.size() != 0)
786  {
787  forAll(nucleation_, model)
788  {
789  nucleationRate_() = Zero;
790  nucleation(i, nucleation_[model]);
791  }
792  }
793  }
794 }
795 
796 
797 void Foam::diameterModels::populationBalanceModel::dmdt()
798 {
799  forAll(velocityGroups_, v)
800  {
801  velocityGroup& velGroup = velocityGroups_[v];
802 
803  velGroup.dmdtRef() = Zero;
804 
805  forAll(sizeGroups_, i)
806  {
807  if (&sizeGroups_[i].phase() == &velGroup.phase())
808  {
809  sizeGroup& fi = sizeGroups_[i];
810 
811  velGroup.dmdtRef() += fi.phase().rho()*(Su_[i] - SuSp_[i]*fi);
812  }
813  }
814  }
815 }
816 
817 
818 void Foam::diameterModels::populationBalanceModel::calcAlphas()
819 {
820  alphas_() = Zero;
821 
822  forAll(velocityGroups_, v)
823  {
824  const phaseModel& phase = velocityGroups_[v].phase();
825 
826  alphas_() += max(phase, phase.residualAlpha());
827  }
828 }
829 
830 
832 Foam::diameterModels::populationBalanceModel::calcDsm()
833 {
834  tmp<volScalarField> tInvDsm
835  (
837  (
838  "invDsm",
839  mesh_,
841  )
842  );
843 
844  volScalarField& invDsm = tInvDsm.ref();
845 
846  forAll(velocityGroups_, v)
847  {
848  const phaseModel& phase = velocityGroups_[v].phase();
849 
850  invDsm += max(phase, phase.residualAlpha())/(phase.d()*alphas_());
851  }
852 
853  return 1.0/tInvDsm;
854 }
855 
856 
857 void Foam::diameterModels::populationBalanceModel::calcVelocity()
858 {
859  U_() = Zero;
860 
861  forAll(velocityGroups_, v)
862  {
863  const phaseModel& phase = velocityGroups_[v].phase();
864 
865  U_() += max(phase, phase.residualAlpha())*phase.U()/alphas_();
866  }
867 }
868 
869 bool Foam::diameterModels::populationBalanceModel::updateSources()
870 {
871  const bool result = sourceUpdateCounter_ % sourceUpdateInterval() == 0;
872 
873  ++ sourceUpdateCounter_;
874 
875  return result;
876 }
877 
878 
879 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
880 
882 (
883  const phaseSystem& fluid,
884  const word& name,
885  HashPtrTable<volScalarField, phasePairKey, phasePairKey::hash>& pDmdt
886 )
887 :
888  regIOobject
889  (
890  IOobject
891  (
892  name,
893  fluid.time().constant(),
894  fluid.mesh()
895  )
896  ),
897  fluid_(fluid),
898  pDmdt_(pDmdt),
899  mesh_(fluid_.mesh()),
900  name_(name),
901  dict_
902  (
903  fluid_.subDict("populationBalanceCoeffs").subDict(name_)
904  ),
905  pimple_(mesh_.lookupObject<pimpleControl>("solutionControl")),
906  continuousPhase_
907  (
908  mesh_.lookupObject<phaseModel>
909  (
910  IOobject::groupName("alpha", dict_.lookup("continuousPhase"))
911  )
912  ),
913  velocityGroups_(),
914  sizeGroups_(),
915  v_(),
916  delta_(),
917  Su_(),
918  SuSp_(),
919  Sui_
920  (
921  IOobject
922  (
923  "Sui",
924  mesh_.time().timeName(),
925  mesh_
926  ),
927  mesh_,
929  ),
930  coalescence_
931  (
932  dict_.lookup("coalescenceModels"),
933  coalescenceModel::iNew(*this)
934  ),
935  coalescenceRate_(),
936  breakup_
937  (
938  dict_.lookup("breakupModels"),
939  breakupModel::iNew(*this)
940  ),
941  breakupRate_(),
942  binaryBreakup_
943  (
944  dict_.lookup("binaryBreakupModels"),
945  binaryBreakupModel::iNew(*this)
946  ),
947  binaryBreakupRate_(),
948  drift_
949  (
950  dict_.lookup("driftModels"),
951  driftModel::iNew(*this)
952  ),
953  driftRate_(),
954  rx_(),
955  rdx_(),
956  nucleation_
957  (
958  dict_.lookup("nucleationModels"),
959  nucleationModel::iNew(*this)
960  ),
961  nucleationRate_(),
962  alphas_(),
963  dsm_(),
964  U_(),
965  sourceUpdateCounter_(0)
966 {
967  this->registerVelocityGroups();
968 
969  this->createPhasePairs();
970 
971  if (sizeGroups_.size() < 3)
972  {
974  << "The populationBalance " << name_
975  << " requires a minimum number of three sizeGroups to be specified."
976  << exit(FatalError);
977  }
978 
979 
980  if (coalescence_.size() != 0)
981  {
982  coalescenceRate_.set
983  (
984  new volScalarField
985  (
986  IOobject
987  (
988  "coalescenceRate",
989  mesh_.time().timeName(),
990  mesh_
991  ),
992  mesh_,
994  )
995  );
996  }
997 
998  if (breakup_.size() != 0)
999  {
1000  breakupRate_.set
1001  (
1002  new volScalarField
1003  (
1004  IOobject
1005  (
1006  "breakupRate",
1007  fluid_.time().timeName(),
1008  mesh_
1009  ),
1010  mesh_,
1012  )
1013  );
1014  }
1015 
1016  if (binaryBreakup_.size() != 0)
1017  {
1018  binaryBreakupRate_.set
1019  (
1020  new volScalarField
1021  (
1022  IOobject
1023  (
1024  "binaryBreakupRate",
1025  fluid_.time().timeName(),
1026  mesh_
1027  ),
1028  mesh_,
1030  (
1031  "binaryBreakupRate",
1033  Zero
1034  )
1035  )
1036  );
1037  }
1038 
1039  if (drift_.size() != 0)
1040  {
1041  driftRate_.set
1042  (
1043  new volScalarField
1044  (
1045  IOobject
1046  (
1047  "driftRate",
1048  fluid_.time().timeName(),
1049  mesh_
1050  ),
1051  mesh_,
1053  )
1054  );
1055 
1056  rx_.set
1057  (
1058  new volScalarField
1059  (
1060  IOobject
1061  (
1062  "r",
1063  fluid_.time().timeName(),
1064  mesh_
1065  ),
1066  mesh_,
1068  )
1069  );
1070 
1071  rdx_.set
1072  (
1073  new volScalarField
1074  (
1075  IOobject
1076  (
1077  "r",
1078  fluid_.time().timeName(),
1079  mesh_
1080  ),
1081  mesh_,
1083  )
1084  );
1085  }
1086 
1087  if (nucleation_.size() != 0)
1088  {
1089  nucleationRate_.set
1090  (
1091  new volScalarField
1092  (
1093  IOobject
1094  (
1095  "nucleationRate",
1096  fluid_.time().timeName(),
1097  mesh_
1098  ),
1099  mesh_,
1101  (
1102  "nucleationRate",
1104  Zero
1105  )
1106  )
1107  );
1108  }
1109 
1110  if (velocityGroups_.size() > 1)
1111  {
1112  alphas_.set
1113  (
1114  new volScalarField
1115  (
1116  IOobject
1117  (
1118  IOobject::groupName("alpha", this->name()),
1119  fluid_.time().timeName(),
1120  mesh_,
1123  ),
1124  mesh_,
1126  )
1127  );
1128 
1129  dsm_.set
1130  (
1131  new volScalarField
1132  (
1133  IOobject
1134  (
1135  IOobject::groupName("d", this->name()),
1136  fluid_.time().timeName(),
1137  mesh_,
1140  ),
1141  mesh_,
1143  )
1144  );
1145 
1146  U_.set
1147  (
1148  new volVectorField
1149  (
1150  IOobject
1151  (
1152  IOobject::groupName("U", this->name()),
1153  fluid_.time().timeName(),
1154  mesh_,
1157  ),
1158  mesh_,
1160  )
1161  );
1162  }
1163 }
1164 
1165 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1166 
1168 {}
1169 
1170 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1171 
1174 {
1175  notImplemented("populationBalance::clone() const");
1176  return autoPtr<populationBalanceModel>(nullptr);
1177 }
1178 
1179 
1181 {
1182  return os.good();
1183 }
1184 
1185 
1188 (
1189  const label i,
1190  const dimensionedScalar& v
1191 ) const
1192 {
1193  const dimensionedScalar& x0 = sizeGroups_[0].x();
1194  const dimensionedScalar& xi = sizeGroups_[i].x();
1195  const dimensionedScalar& xm = sizeGroups_.last().x();
1196  dimensionedScalar lowerBoundary(x0);
1197  dimensionedScalar upperBoundary(xm);
1198 
1199  if (i != 0) lowerBoundary = sizeGroups_[i-1].x();
1200 
1201  if (i != sizeGroups_.size() - 1) upperBoundary = sizeGroups_[i+1].x();
1202 
1203  if ((i == 0 && v < x0) || (i == sizeGroups_.size() - 1 && v > xm))
1204  {
1205  return v/xi;
1206  }
1207  else if (v < lowerBoundary || v > upperBoundary)
1208  {
1209  return 0;
1210  }
1211  else if (v.value() == xi.value())
1212  {
1213  return 1.0;
1214  }
1215  else if (v > xi)
1216  {
1217  return (upperBoundary - v)/(upperBoundary - xi);
1218  }
1219  else
1220  {
1221  return (v - lowerBoundary)/(xi - lowerBoundary);
1222  }
1223 }
1224 
1225 
1228 (
1229  const phaseModel& dispersedPhase
1230 ) const
1231 {
1232  return
1233  fluid_.lookupSubModel<surfaceTensionModel>
1234  (
1235  phasePair(dispersedPhase, continuousPhase_)
1236  ).sigma();
1237 }
1238 
1239 
1242 {
1243  return
1244  mesh_.lookupObject<phaseCompressibleMomentumTransportModel>
1245  (
1247  (
1248  momentumTransportModel::typeName,
1249  continuousPhase_.name()
1250  )
1251  );
1252 }
1253 
1254 
1256 {
1257  const dictionary& solutionControls = mesh_.solverDict(name_);
1258  bool solveOnFinalIterOnly =
1259  solutionControls.lookupOrDefault<bool>("solveOnFinalIterOnly", false);
1260 
1261  if (!solveOnFinalIterOnly || pimple_.finalPimpleIter())
1262  {
1263  const label nCorr = this->nCorr();
1264  const scalar tolerance =
1265  solutionControls.lookup<scalar>("tolerance");
1266 
1267  if (nCorr > 0)
1268  {
1269  correct();
1270  }
1271 
1272  int iCorr = 0;
1273  scalar maxInitialResidual = 1;
1274 
1275  while (++iCorr <= nCorr && maxInitialResidual > tolerance)
1276  {
1277  Info<< "populationBalance "
1278  << this->name()
1279  << ": Iteration "
1280  << iCorr
1281  << endl;
1282 
1283  if (updateSources())
1284  {
1285  sources();
1286  dmdt();
1287  }
1288 
1289  maxInitialResidual = 0;
1290 
1291  forAll(sizeGroups_, i)
1292  {
1293  sizeGroup& fi = sizeGroups_[i];
1294  const phaseModel& phase = fi.phase();
1295  const volScalarField& alpha = phase;
1296  const dimensionedScalar& residualAlpha = phase.residualAlpha();
1297  const volScalarField& rho = phase.thermo().rho();
1298 
1299  fvScalarMatrix sizeGroupEqn
1300  (
1301  fvm::ddt(alpha, rho, fi)
1302  + fi.VelocityGroup().mvConvection()->fvmDiv
1303  (
1304  phase.alphaRhoPhi(),
1305  fi
1306  )
1307  + fvm::SuSp
1308  (
1309  fi.VelocityGroup().dmdt()
1310  - (fvc::ddt(alpha, rho) + fvc::div(phase.alphaRhoPhi())),
1311  fi
1312  )
1313  ==
1314  fvc::Su(Su_[i]*rho, fi)
1315  - fvm::SuSp(SuSp_[i]*rho, fi)
1316  + fvc::ddt(residualAlpha*rho, fi)
1317  - fvm::ddt(residualAlpha*rho, fi)
1318  );
1319 
1320  sizeGroupEqn.relax();
1321 
1322  maxInitialResidual = max
1323  (
1324  sizeGroupEqn.solve().initialResidual(),
1325  maxInitialResidual
1326  );
1327  }
1328  }
1329 
1330  if (nCorr > 0)
1331  {
1332  forAll(sizeGroups_, i)
1333  {
1334  sizeGroups_[i].correct();
1335  }
1336 
1337  forAll(velocityGroups_, i)
1338  {
1339  velocityGroups_[i].postSolve();
1340  }
1341  }
1342 
1343  if (velocityGroups_.size() > 1)
1344  {
1345  calcAlphas();
1346  dsm_() = calcDsm();
1347  calcVelocity();
1348  }
1349 
1350  volScalarField fAlpha0
1351  (
1352  sizeGroups_.first()*sizeGroups_.first().phase()
1353  );
1354 
1355  volScalarField fAlphaN
1356  (
1357  sizeGroups_.last()*sizeGroups_.last().phase()
1358  );
1359 
1360  Info<< this->name() << " sizeGroup phase fraction first, last = "
1361  << fAlpha0.weightedAverage(this->mesh().V()).value()
1362  << ' ' << fAlphaN.weightedAverage(this->mesh().V()).value()
1363  << endl;
1364  }
1365 }
1366 
1367 // ************************************************************************* //
Templated abstract base class for multiphase compressible turbulence models.
const char *const group
Group name for atomic constants.
virtual void correct()
Correct diameter independent expressions.
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
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
PhaseCompressibleMomentumTransportModel< phaseModel > phaseCompressibleMomentumTransportModel
Typedef for phaseCompressibleMomentumTransportModel.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
label phasei
Definition: pEqn.H:27
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
populationBalanceModel(const phaseSystem &fluid, const word &name, HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > &pDmdt)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
const dimensionedScalar & sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
static int compare(const Pair< word > &a, const Pair< word > &b)
Compare Pairs.
Definition: Pair.H:142
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
label k
Boltzmann constant.
Generic dimensioned Type class.
bool writeData(Ostream &) const
Dummy write for regIOobject.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
dimensionedScalar neg(const dimensionedScalar &ds)
const dimensionedScalar eta(const label i, const dimensionedScalar &v) const
Return allocation coefficient.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
autoPtr< populationBalanceModel > clone() const
Return clone.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Calculate the first temporal derivative.
dimensionedScalar pos(const dimensionedScalar &ds)
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
static word groupName(Name name, const word &group)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Calculate the matrix for the first temporal derivative.
const phaseCompressibleMomentumTransportModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
Calculate the field for explicit evaluation of implicit and explicit sources.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Calculate the divergence of the given field.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-8/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()
Definition: pEqn.H:68
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Internal & ref()
Return a reference to the dimensioned internal field.
#define notImplemented(functionName)
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:356
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
int nCorr
Definition: createControls.H:3
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const tmp< volScalarField > sigmaWithContinuousPhase(const phaseModel &dispersedPhase) const
Return the surface tension coefficient between a given dispersed.
A class for managing temporary objects.
Definition: PtrList.H:53
void solve()
Solve the population balance equation.
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
Calculate the matrix for implicit and explicit sources.
Namespace for OpenFOAM.
const dimensionSet dimVelocity