populationBalanceModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2017-2019 OpenFOAM Foundation
9  Copyright (C) 2020-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "populationBalanceModel.H"
30 #include "coalescenceModel.H"
31 #include "breakupModel.H"
32 #include "binaryBreakupModel.H"
33 #include "driftModel.H"
34 #include "nucleationModel.H"
35 #include "phaseSystem.H"
36 #include "surfaceTensionModel.H"
37 #include "fvmDdt.H"
38 #include "fvcDdt.H"
39 #include "fvmSup.H"
40 #include "fvcSup.H"
41 #include "fvcDiv.H"
42 #include "phaseCompressibleTurbulenceModel.H"
43 
44 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
45 
46 void Foam::diameterModels::populationBalanceModel::registerVelocityGroups()
47 {
48  forAll(fluid_.phases(), phasei)
49  {
50  if (isA<velocityGroup>(fluid_.phases()[phasei].dPtr()()))
51  {
52  const velocityGroup& velGroup =
53  refCast<const velocityGroup>(fluid_.phases()[phasei].dPtr()());
54 
55  if (velGroup.popBalName() == this->name())
56  {
57  velocityGroups_.resize(velocityGroups_.size() + 1);
58 
59  velocityGroups_.set
60  (
61  velocityGroups_.size() - 1,
62  &const_cast<velocityGroup&>(velGroup)
63  );
64 
65  forAll(velGroup.sizeGroups(), i)
66  {
67  this->registerSizeGroups
68  (
69  const_cast<sizeGroup&>(velGroup.sizeGroups()[i])
70  );
71  }
72  }
73  }
74  }
75 }
76 
77 
78 void Foam::diameterModels::populationBalanceModel::registerSizeGroups
79 (
80  sizeGroup& group
81 )
82 {
83  if
84  (
85  sizeGroups_.size() != 0
86  &&
87  group.x().value() <= sizeGroups_.last().x().value()
88  )
89  {
91  << "Size groups must be entered according to their representative"
92  << " size"
93  << exit(FatalError);
94  }
95 
96  sizeGroups_.resize(sizeGroups_.size() + 1);
97  sizeGroups_.set(sizeGroups_.size() - 1, &group);
98 
99  // Grid generation over property space
100  if (sizeGroups_.size() == 1)
101  {
102  // Set the first sizeGroup boundary to the representative value of
103  // the first sizeGroup
104  v_.append
105  (
107  (
108  "v",
109  sizeGroups_.last().x()
110  )
111  );
112 
113  // Set the last sizeGroup boundary to the representative size of the
114  // last sizeGroup
115  v_.append
116  (
118  (
119  "v",
120  sizeGroups_.last().x()
121  )
122  );
123  }
124  else
125  {
126  // Overwrite the next-to-last sizeGroup boundary to lie halfway between
127  // the last two sizeGroups
128  v_.last() =
129  0.5
130  *(
131  sizeGroups_[sizeGroups_.size()-2].x()
132  + sizeGroups_.last().x()
133  );
134 
135  // Set the last sizeGroup boundary to the representative size of the
136  // last sizeGroup
137  v_.append
138  (
140  (
141  "v",
142  sizeGroups_.last().x()
143  )
144  );
145  }
146 
147  delta_.append(new PtrList<dimensionedScalar>());
148 
149  Su_.append
150  (
152  (
153  "Su",
155  mesh_,
157  )
158  );
159 
160  SuSp_.append
161  (
163  (
164  "SuSp",
166  mesh_,
168  )
169  );
170 }
171 
172 
173 void Foam::diameterModels::populationBalanceModel::createPhasePairs()
174 {
175  forAll(velocityGroups_, i)
176  {
177  const phaseModel& phasei = velocityGroups_[i].phase();
178 
179  forAll(velocityGroups_, j)
180  {
181  const phaseModel& phasej = velocityGroups_[j].phase();
182 
183  if (&phasei != &phasej)
184  {
185  const phasePairKey key
186  (
187  phasei.name(),
188  phasej.name(),
189  false
190  );
191 
192  if (!phasePairs_.found(key))
193  {
194  phasePairs_.insert
195  (
196  key,
197  autoPtr<phasePair>
198  (
199  new phasePair
200  (
201  phasei,
202  phasej
203  )
204  )
205  );
206  }
207  }
208  }
209  }
210 }
211 
212 
213 void Foam::diameterModels::populationBalanceModel::correct()
214 {
215  calcDeltas();
216 
217  forAll(velocityGroups_, v)
218  {
219  velocityGroups_[v].preSolve();
220  }
221 
222  forAll(coalescence_, model)
223  {
224  coalescence_[model].correct();
225  }
226 
227  forAll(breakup_, model)
228  {
229  breakup_[model].correct();
230 
231  breakup_[model].dsdPtr()().correct();
232  }
233 
234  forAll(binaryBreakup_, model)
235  {
236  binaryBreakup_[model].correct();
237  }
238 
239  forAll(drift_, model)
240  {
241  drift_[model].correct();
242  }
243 
244  forAll(nucleation_, model)
245  {
246  nucleation_[model].correct();
247  }
248 }
249 
250 
251 void Foam::diameterModels::populationBalanceModel::
252 birthByCoalescence
253 (
254  const label j,
255  const label k
256 )
257 {
258  const sizeGroup& fj = sizeGroups_[j];
259  const sizeGroup& fk = sizeGroups_[k];
260 
261  dimensionedScalar Gamma;
262  dimensionedScalar v = fj.x() + fk.x();
263 
264  for (label i = j; i < sizeGroups_.size(); i++)
265  {
266  // Calculate fraction for intra-interval events
267  Gamma = gamma(i, v);
268 
269  if (Gamma.value() == 0) continue;
270 
271  const sizeGroup& fi = sizeGroups_[i];
272 
273  // Avoid double counting of events
274  if (j == k)
275  {
276  Sui_ =
277  0.5*fi.x()*coalescenceRate_()*Gamma
278  *fj*fj.phase()/fj.x()
279  *fk*fk.phase()/fk.x();
280  }
281  else
282  {
283  Sui_ =
284  fi.x()*coalescenceRate_()*Gamma
285  *fj*fj.phase()/fj.x()
286  *fk*fk.phase()/fk.x();
287  }
288 
289  Su_[i] += Sui_;
290 
291  const phasePairKey pairij
292  (
293  fi.phase().name(),
294  fj.phase().name()
295  );
296 
297  if (pDmdt_.found(pairij))
298  {
299  const scalar dmdtSign
300  (
301  Pair<word>::compare(pDmdt_.find(pairij).key(), pairij)
302  );
303 
304  pDmdt_[pairij]->ref() += dmdtSign*fj.x()/v*Sui_*fi.phase().rho();
305  }
306 
307  const phasePairKey pairik
308  (
309  fi.phase().name(),
310  fk.phase().name()
311  );
312 
313  if (pDmdt_.found(pairik))
314  {
315  const scalar dmdtSign
316  (
317  Pair<word>::compare(pDmdt_.find(pairik).key(), pairik)
318  );
319 
320  pDmdt_[pairik]->ref() += dmdtSign*fk.x()/v*Sui_*fi.phase().rho();
321  }
322  }
323 }
324 
325 
326 void Foam::diameterModels::populationBalanceModel::
327 deathByCoalescence
328 (
329  const label i,
330  const label j
331 )
332 {
333  const sizeGroup& fi = sizeGroups_[i];
334  const sizeGroup& fj = sizeGroups_[j];
335 
336  SuSp_[i] += coalescenceRate_()*fi.phase()*fj*fj.phase()/fj.x();
337 
338  if (i != j)
339  {
340  SuSp_[j] += coalescenceRate_()*fj.phase()*fi*fi.phase()/fi.x();
341  }
342 }
343 
344 
345 void Foam::diameterModels::populationBalanceModel::
346 birthByBreakup
347 (
348  const label k,
349  const label model
350 )
351 {
352  const sizeGroup& fk = sizeGroups_[k];
353 
354  for (label i = 0; i <= k; i++)
355  {
356  const sizeGroup& fi = sizeGroups_[i];
357 
358  Sui_ =
359  fi.x()*breakupRate_()*breakup_[model].dsdPtr()().nik(i, k)
360  *fk*fk.phase()/fk.x();
361 
362  Su_[i] += Sui_;
363 
364  const phasePairKey pair
365  (
366  fi.phase().name(),
367  fk.phase().name()
368  );
369 
370  if (pDmdt_.found(pair))
371  {
372  const scalar dmdtSign
373  (
374  Pair<word>::compare(pDmdt_.find(pair).key(), pair)
375  );
376 
377  pDmdt_[pair]->ref() += dmdtSign*Sui_*fi.phase().rho();
378  }
379  }
380 }
381 
382 
383 void Foam::diameterModels::populationBalanceModel::deathByBreakup(const label i)
384 {
385  const sizeGroup& fi = sizeGroups_[i];
386 
387  SuSp_[i] += breakupRate_()*fi.phase();
388 }
389 
390 
391 void Foam::diameterModels::populationBalanceModel::calcDeltas()
392 {
393  forAll(sizeGroups_, i)
394  {
395  if (delta_[i].empty())
396  {
397  for (label j = 0; j <= sizeGroups_.size() - 1; j++)
398  {
399  delta_[i].append
400  (
402  (
403  "delta",
404  dimVolume,
405  v_[i+1].value() - v_[i].value()
406  )
407  );
408 
409  if
410  (
411  v_[i].value() < 0.5*sizeGroups_[j].x().value()
412  &&
413  0.5*sizeGroups_[j].x().value() < v_[i+1].value()
414  )
415  {
416  delta_[i][j] = mag(0.5*sizeGroups_[j].x() - v_[i]);
417  }
418  else if (0.5*sizeGroups_[j].x().value() <= v_[i].value())
419  {
420  delta_[i][j].value() = 0;
421  }
422  }
423  }
424  }
425 }
426 
427 
428 void Foam::diameterModels::populationBalanceModel::
429 birthByBinaryBreakup
430 (
431  const label i,
432  const label j
433 )
434 {
435  const sizeGroup& fj = sizeGroups_[j];
436  const sizeGroup& fi = sizeGroups_[i];
437 
438  Sui_ = fi.x()*binaryBreakupRate_()*delta_[i][j]*fj*fj.phase()/fj.x();
439 
440  Su_[i] += Sui_;
441 
442  const phasePairKey pairij
443  (
444  fi.phase().name(),
445  fj.phase().name()
446  );
447 
448  if (pDmdt_.found(pairij))
449  {
450  const scalar dmdtSign
451  (
452  Pair<word>::compare(pDmdt_.find(pairij).key(), pairij)
453  );
454 
455  pDmdt_[pairij]->ref() += dmdtSign*Sui_*fi.phase().rho();
456  }
457 
458  dimensionedScalar Gamma;
459  dimensionedScalar v = fj.x() - fi.x();
460 
461  for (label k = 0; k <= j; k++)
462  {
463  // Calculate fraction for intra-interval events
464  Gamma = gamma(k, v);
465 
466  if (Gamma.value() == 0) continue;
467 
468  const sizeGroup& fk = sizeGroups_[k];
469 
470  volScalarField& Suk = Sui_;
471 
472  Suk =
473  sizeGroups_[k].x()*binaryBreakupRate_()*delta_[i][j]*Gamma
474  *fj*fj.phase()/fj.x();
475 
476  Su_[k] += Suk;
477 
478  const phasePairKey pairkj
479  (
480  fk.phase().name(),
481  fj.phase().name()
482  );
483 
484  if (pDmdt_.found(pairkj))
485  {
486  const scalar dmdtSign
487  (
489  (
490  pDmdt_.find(pairkj).key(),
491  pairkj
492  )
493  );
494 
495  pDmdt_[pairkj]->ref() += dmdtSign*Suk*fi.phase().rho();
496  }
497  }
498 }
499 
500 
501 void Foam::diameterModels::populationBalanceModel::
502 deathByBinaryBreakup
503 (
504  const label j,
505  const label i
506 )
507 {
508  const volScalarField& alphai = sizeGroups_[i].phase();
509 
510  SuSp_[i] += alphai*binaryBreakupRate_()*delta_[j][i];
511 }
512 
513 
514 void Foam::diameterModels::populationBalanceModel::drift(const label i)
515 {
516  const sizeGroup& fp = sizeGroups_[i];
517 
518  if (i == 0)
519  {
520  rx_() = pos(driftRate_())*sizeGroups_[i+1].x()/sizeGroups_[i].x();
521  }
522  else if (i == sizeGroups_.size() - 1)
523  {
524  rx_() = neg(driftRate_())*sizeGroups_[i-1].x()/sizeGroups_[i].x();
525  }
526  else
527  {
528  rx_() =
529  pos(driftRate_())*sizeGroups_[i+1].x()/sizeGroups_[i].x()
530  + neg(driftRate_())*sizeGroups_[i-1].x()/sizeGroups_[i].x();
531  }
532 
533  SuSp_[i] +=
534  (neg(1 - rx_()) + neg(rx_() - rx_()/(1 - rx_())))*driftRate_()
535  *fp.phase()/((rx_() - 1)*fp.x());
536 
537  rx_() *= 0.0;
538  rdx_() *= 0.0;
539 
540  if (i == sizeGroups_.size() - 2)
541  {
542  rx_() = pos(driftRate_())*sizeGroups_[i+1].x()/sizeGroups_[i].x();
543 
544  rdx_() =
545  pos(driftRate_())
546  *(sizeGroups_[i+1].x() - sizeGroups_[i].x())
547  /(sizeGroups_[i].x() - sizeGroups_[i-1].x());
548  }
549  else if (i < sizeGroups_.size() - 2)
550  {
551  rx_() = pos(driftRate_())*sizeGroups_[i+2].x()/sizeGroups_[i+1].x();
552 
553  rdx_() =
554  pos(driftRate_())
555  *(sizeGroups_[i+2].x() - sizeGroups_[i+1].x())
556  /(sizeGroups_[i+1].x() - sizeGroups_[i].x());
557  }
558 
559  if (i == 1)
560  {
561  rx_() += neg(driftRate_())*sizeGroups_[i-1].x()/sizeGroups_[i].x();
562 
563  rdx_() +=
564  neg(driftRate_())
565  *(sizeGroups_[i].x() - sizeGroups_[i-1].x())
566  /(sizeGroups_[i+1].x() - sizeGroups_[i].x());
567  }
568  else if (i > 1)
569  {
570  rx_() += neg(driftRate_())*sizeGroups_[i-2].x()/sizeGroups_[i-1].x();
571 
572  rdx_() +=
573  neg(driftRate_())
574  *(sizeGroups_[i-1].x() - sizeGroups_[i-2].x())
575  /(sizeGroups_[i].x() - sizeGroups_[i-1].x());
576  }
577 
578  if (i != sizeGroups_.size() - 1)
579  {
580  const sizeGroup& fe = sizeGroups_[i+1];
581  volScalarField& Sue = Sui_;
582 
583  Sue =
584  pos(driftRate_())*driftRate_()*rdx_()
585  *fp*fp.phase()/fp.x()
586  /(rx_() - 1);
587 
588  Su_[i+1] += Sue;
589 
590  const phasePairKey pairij
591  (
592  fp.phase().name(),
593  fe.phase().name()
594  );
595 
596  if (pDmdt_.found(pairij))
597  {
598  const scalar dmdtSign
599  (
600  Pair<word>::compare(pDmdt_.find(pairij).key(), pairij)
601  );
602 
603  pDmdt_[pairij]->ref() -= dmdtSign*Sue*fp.phase().rho();
604  }
605  }
606 
607  if (i != 0)
608  {
609  const sizeGroup& fw = sizeGroups_[i-1];
610  volScalarField& Suw = Sui_;
611 
612  Suw =
613  neg(driftRate_())*driftRate_()*rdx_()
614  *fp*fp.phase()/fp.x()
615  /(rx_() - 1);
616 
617  Su_[i-1] += Suw;
618 
619  const phasePairKey pairih
620  (
621  fp.phase().name(),
622  fw.phase().name()
623  );
624 
625  if (pDmdt_.found(pairih))
626  {
627  const scalar dmdtSign
628  (
629  Pair<word>::compare(pDmdt_.find(pairih).key(), pairih)
630  );
631 
632  pDmdt_[pairih]->ref() -= dmdtSign*Suw*fp.phase().rho();
633  }
634  }
635 }
636 
637 
638 void Foam::diameterModels::populationBalanceModel::nucleation(const label i)
639 {
640  Su_[i] += sizeGroups_[i].x()*nucleationRate_();
641 }
642 
643 
644 void Foam::diameterModels::populationBalanceModel::sources()
645 {
646  forAll(sizeGroups_, i)
647  {
648  Su_[i] *= 0.0;
649  SuSp_[i] *= 0.0;
650  }
651 
653  (
654  phasePairTable,
655  phasePairs(),
656  phasePairIter
657  )
658  {
659  pDmdt_(phasePairIter())->ref() *= 0.0;
660  }
661 
662  // Since the calculation of the rates is computationally expensive,
663  // they are calculated once for each sizeGroup pair and inserted into source
664  // terms as required
665  forAll(sizeGroups_, i)
666  {
667  const sizeGroup& fi = sizeGroups_[i];
668 
669  if (coalescence_.size() != 0)
670  {
671  for (label j = 0; j <= i; j++)
672  {
673  const sizeGroup& fj = sizeGroups_[j];
674 
675  if (fi.x() + fj.x() > sizeGroups_.last().x()) break;
676 
677  coalescenceRate_() *= 0.0;
678 
679  forAll(coalescence_, model)
680  {
681  coalescence_[model].addToCoalescenceRate
682  (
683  coalescenceRate_(),
684  i,
685  j
686  );
687  }
688 
689  birthByCoalescence(i, j);
690 
691  deathByCoalescence(i, j);
692  }
693  }
694 
695  if (breakup_.size() != 0)
696  {
697  forAll(breakup_, model)
698  {
699  breakup_[model].setBreakupRate(breakupRate_(), i);
700 
701  birthByBreakup(i, model);
702 
703  deathByBreakup(i);
704  }
705  }
706 
707  if (binaryBreakup_.size() != 0)
708  {
709  label j = 0;
710 
711  while (delta_[j][i].value() != 0)
712  {
713  binaryBreakupRate_() *= 0.0;
714 
715  forAll(binaryBreakup_, model)
716  {
717  binaryBreakup_[model].addToBinaryBreakupRate
718  (
719  binaryBreakupRate_(),
720  j,
721  i
722  );
723  }
724 
725  birthByBinaryBreakup(j, i);
726 
727  deathByBinaryBreakup(j, i);
728 
729  j++;
730  }
731  }
732 
733  if (drift_.size() != 0)
734  {
735  driftRate_() *= 0.0;
736 
737  forAll(drift_, model)
738  {
739  drift_[model].addToDriftRate(driftRate_(), i);
740  }
741 
742  drift(i);
743  }
744 
745  if (nucleation_.size() != 0)
746  {
747  nucleationRate_() *= 0.0;
748 
749  forAll(nucleation_, model)
750  {
751  nucleation_[model].addToNucleationRate(nucleationRate_(), i);
752  }
753 
754  nucleation(i);
755  }
756  }
757 }
758 
759 
760 void Foam::diameterModels::populationBalanceModel::dmdt()
761 {
762  forAll(velocityGroups_, v)
763  {
764  velocityGroup& velGroup = velocityGroups_[v];
765 
766  velGroup.dmdtRef() *= 0.0;
767 
768  forAll(sizeGroups_, i)
769  {
770  if (&sizeGroups_[i].phase() == &velGroup.phase())
771  {
772  sizeGroup& fi = sizeGroups_[i];
773 
774  velGroup.dmdtRef() += fi.phase().rho()*(Su_[i] - SuSp_[i]*fi);
775  }
776  }
777  }
778 }
779 
780 
781 void Foam::diameterModels::populationBalanceModel::calcAlphas()
782 {
783  alphas_() *= 0.0;
784 
785  forAll(velocityGroups_, v)
786  {
787  const phaseModel& phase = velocityGroups_[v].phase();
788 
789  alphas_() += max(phase, phase.residualAlpha());
790  }
791 }
792 
793 
795 Foam::diameterModels::populationBalanceModel::calcDsm()
796 {
797  tmp<volScalarField> tInvDsm
798  (
800  (
801  "invDsm",
802  mesh_,
804  )
805  );
806 
807  volScalarField& invDsm = tInvDsm.ref();
808 
809  forAll(velocityGroups_, v)
810  {
811  const phaseModel& phase = velocityGroups_[v].phase();
812 
813  invDsm += max(phase, phase.residualAlpha())/(phase.d()*alphas_());
814  }
815 
816  return 1.0/tInvDsm;
817 }
818 
819 
820 void Foam::diameterModels::populationBalanceModel::calcVelocity()
821 {
822  U_() *= 0.0;
823 
824  forAll(velocityGroups_, v)
825  {
826  const phaseModel& phase = velocityGroups_[v].phase();
827 
828  U_() += max(phase, phase.residualAlpha())*phase.U()/alphas_();
829  }
830 }
831 
832 bool Foam::diameterModels::populationBalanceModel::updateSources()
833 {
834  const bool result = sourceUpdateCounter_ % sourceUpdateInterval() == 0;
835 
836  ++ sourceUpdateCounter_;
837 
838  return result;
839 }
840 
841 
842 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
843 
845 (
846  const phaseSystem& fluid,
847  const word& name,
849 )
850 :
852  (
853  IOobject
854  (
855  name,
856  fluid.time().constant(),
857  fluid.mesh()
858  )
859  ),
860  fluid_(fluid),
861  pDmdt_(pDmdt),
862  mesh_(fluid_.mesh()),
863  name_(name),
864  dict_
865  (
866  fluid_.subDict("populationBalanceCoeffs").subDict(name_)
867  ),
868  pimple_(mesh_.lookupObject<pimpleControl>("solutionControl")),
869  continuousPhase_
870  (
871  mesh_.lookupObject<phaseModel>
872  (
873  IOobject::groupName("alpha", dict_.get<word>("continuousPhase"))
874  )
875  ),
876  velocityGroups_(),
877  sizeGroups_(),
878  v_(),
879  delta_(),
880  Su_(),
881  SuSp_(),
882  Sui_
883  (
884  IOobject
885  (
886  "Sui",
887  mesh_.time().timeName(),
888  mesh_
889  ),
890  mesh_,
891  dimensionedScalar("zero", inv(dimTime), Zero)
892  ),
893  coalescence_
894  (
895  dict_.lookup("coalescenceModels"),
896  coalescenceModel::iNew(*this)
897  ),
898  coalescenceRate_(),
899  breakup_
900  (
901  dict_.lookup("breakupModels"),
902  breakupModel::iNew(*this)
903  ),
904  breakupRate_(),
905  binaryBreakup_
906  (
907  dict_.lookup("binaryBreakupModels"),
908  binaryBreakupModel::iNew(*this)
909  ),
910  binaryBreakupRate_(),
911  drift_
912  (
913  dict_.lookup("driftModels"),
914  driftModel::iNew(*this)
915  ),
916  driftRate_(),
917  rx_(),
918  rdx_(),
919  nucleation_
920  (
921  dict_.lookup("nucleationModels"),
922  nucleationModel::iNew(*this)
923  ),
924  nucleationRate_(),
925  alphas_(),
926  dsm_(),
927  U_(),
928  sourceUpdateCounter_
929  (
930  (mesh_.time().timeIndex()*nCorr())%sourceUpdateInterval()
931  )
932 {
933  this->registerVelocityGroups();
934 
935  this->createPhasePairs();
936 
937  if (sizeGroups_.size() < 3)
938  {
940  << "The populationBalance " << name_
941  << " requires a minimum number of three sizeGroups to be specified."
942  << exit(FatalError);
943  }
944 
945 
946  if (coalescence_.size() != 0)
947  {
948  coalescenceRate_.reset
949  (
950  new volScalarField
951  (
952  mesh_.newIOobject("coalescenceRate"),
953  mesh_,
955  )
956  );
957  }
958 
959  if (breakup_.size() != 0)
960  {
961  breakupRate_.reset
962  (
963  new volScalarField
964  (
965  mesh_.newIOobject("breakupRate"),
966  mesh_,
968  )
969  );
970  }
971 
972  if (binaryBreakup_.size() != 0)
973  {
974  binaryBreakupRate_.reset
975  (
976  new volScalarField
977  (
978  mesh_.newIOobject("binaryBreakupRate"),
979  mesh_,
981  )
982  );
983  }
984 
985  if (drift_.size() != 0)
986  {
987  driftRate_.reset
988  (
989  new volScalarField
990  (
991  mesh_.newIOobject("driftRate"),
992  mesh_,
994  )
995  );
996 
997  rx_.reset
998  (
999  new volScalarField
1000  (
1001  mesh_.newIOobject("rx"),
1002  mesh_,
1004  )
1005  );
1006 
1007  rdx_.reset
1008  (
1009  new volScalarField
1010  (
1011  mesh_.newIOobject("rdx"),
1012  mesh_,
1014  )
1015  );
1016  }
1017 
1018  if (nucleation_.size() != 0)
1019  {
1020  nucleationRate_.reset
1021  (
1022  new volScalarField
1023  (
1024  mesh_.newIOobject("nucleationRate"),
1025  mesh_,
1027  )
1028  );
1029  }
1030 
1031  if (velocityGroups_.size() > 1)
1032  {
1033  alphas_.reset
1034  (
1035  new volScalarField
1036  (
1037  IOobject
1038  (
1039  IOobject::groupName("alpha", this->name()),
1040  fluid_.time().timeName(),
1041  mesh_,
1045  ),
1046  mesh_,
1048  )
1049  );
1050 
1051  dsm_.reset
1052  (
1053  new volScalarField
1054  (
1055  IOobject
1056  (
1057  IOobject::groupName("d", this->name()),
1058  fluid_.time().timeName(),
1059  mesh_,
1063  ),
1064  mesh_,
1066  )
1067  );
1068 
1069  U_.reset
1070  (
1071  new volVectorField
1072  (
1073  IOobject
1074  (
1075  IOobject::groupName("U", this->name()),
1076  fluid_.time().timeName(),
1077  mesh_,
1081  ),
1082  mesh_,
1084  )
1085  );
1086  }
1087 }
1088 
1089 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1090 
1092 {}
1093 
1094 
1095 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1096 
1099 {
1101  return nullptr;
1102 }
1103 
1104 
1107  return os.good();
1108 }
1109 
1110 
1113 (
1114  const label i,
1115  const dimensionedScalar& v
1116 ) const
1117 {
1118  dimensionedScalar lowerBoundary(v);
1119  dimensionedScalar upperBoundary(v);
1120  const dimensionedScalar& xi = sizeGroups_[i].x();
1121 
1122  if (i == 0)
1123  {
1124  lowerBoundary = xi;
1125  }
1126  else
1127  {
1128  lowerBoundary = sizeGroups_[i-1].x();
1129  }
1130 
1131  if (i == sizeGroups_.size() - 1)
1132  {
1133  upperBoundary = xi;
1134  }
1135  else
1136  {
1137  upperBoundary = sizeGroups_[i+1].x();
1138  }
1139 
1140  if (v < lowerBoundary || v > upperBoundary)
1141  {
1142  return 0;
1143  }
1144  else if (v.value() <= xi.value())
1145  {
1146  return (v - lowerBoundary)/(xi - lowerBoundary);
1147  }
1148  else
1149  {
1150  return (upperBoundary - v)/(upperBoundary - xi);
1151  }
1152 }
1153 
1154 
1157 (
1158  const phaseModel& dispersedPhase
1159 ) const
1160 {
1161  return
1162  fluid_.lookupSubModel<reactingMultiphaseEuler::surfaceTensionModel>
1163  (
1164  phasePair(dispersedPhase, continuousPhase_)
1165  ).sigma();
1166 }
1167 
1168 
1171 {
1172  return
1173  mesh_.lookupObject<phaseCompressibleTurbulenceModel>
1174  (
1176  (
1178  continuousPhase_.name()
1179  )
1180  );
1181 }
1182 
1183 
1185 {
1186  const dictionary& solutionControls = mesh_.solverDict(name_);
1187  const bool solveOnFinalIterOnly =
1188  solutionControls.getOrDefault("solveOnFinalIterOnly", false);
1189 
1190  if (!solveOnFinalIterOnly || pimple_.finalIter())
1191  {
1192  const label nCorr = this->nCorr();
1193  const scalar tolerance = solutionControls.get<scalar>("tolerance");
1194 
1195  if (nCorr > 0)
1196  {
1197  correct();
1198  }
1199 
1200  int iCorr = 0;
1201  scalar maxInitialResidual = 1;
1202 
1203  while (++iCorr <= nCorr && maxInitialResidual > tolerance)
1204  {
1205  Info<< "populationBalance "
1206  << this->name()
1207  << ": Iteration "
1208  << iCorr
1209  << endl;
1210 
1211  if (updateSources())
1212  {
1213  sources();
1214  }
1215 
1216  dmdt();
1217 
1218  maxInitialResidual = 0;
1219 
1220  forAll(sizeGroups_, i)
1221  {
1222  sizeGroup& fi = sizeGroups_[i];
1223  const phaseModel& phase = fi.phase();
1224  const volScalarField& alpha = phase;
1225  const dimensionedScalar& residualAlpha = phase.residualAlpha();
1226  const tmp<volScalarField> trho(phase.thermo().rho());
1227  const volScalarField& rho = trho();
1228 
1229  fvScalarMatrix sizeGroupEqn
1230  (
1231  fvm::ddt(alpha, rho, fi)
1232  + fi.VelocityGroup().mvConvection()->fvmDiv
1233  (
1234  phase.alphaRhoPhi(),
1235  fi
1236  )
1237  - fvm::Sp
1238  (
1239  fvc::ddt(alpha, rho) + fvc::div(phase.alphaRhoPhi())
1240  - fi.VelocityGroup().dmdt(),
1241  fi
1242  )
1243  ==
1244  fvc::Su(Su_[i]*rho, fi)
1245  - fvm::SuSp(SuSp_[i]*rho, fi)
1246  + fvc::ddt(residualAlpha*rho, fi)
1247  - fvm::ddt(residualAlpha*rho, fi)
1248  );
1249 
1250  sizeGroupEqn.relax();
1251 
1252  maxInitialResidual = max
1253  (
1254  sizeGroupEqn.solve().initialResidual(),
1255  maxInitialResidual
1256  );
1257  }
1258  }
1259 
1260  if (nCorr > 0)
1261  {
1262  forAll(velocityGroups_, i)
1263  {
1264  velocityGroups_[i].postSolve();
1265  }
1266  }
1267 
1268  if (velocityGroups_.size() > 1)
1269  {
1270  calcAlphas();
1271  dsm_() = calcDsm();
1272  calcVelocity();
1273  }
1274 
1275  volScalarField fAlpha0
1276  (
1277  sizeGroups_.first()*sizeGroups_.first().phase()
1278  );
1279 
1280  volScalarField fAlphaN
1281  (
1282  sizeGroups_.last()*sizeGroups_.last().phase()
1283  );
1284 
1285  Info<< this->name() << " sizeGroup phase fraction first, last = "
1286  << fAlpha0.weightedAverage(this->mesh().V()).value()
1287  << ' ' << fAlphaN.weightedAverage(this->mesh().V()).value()
1288  << endl;
1289  }
1290 }
1291 
1292 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
twoPhaseSystem & fluid
const Type & value() const noexcept
Return const reference to value.
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
Base class for drift models.
Definition: driftModel.H:48
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
label phasei
Definition: pEqn.H:27
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
populationBalanceModel(const phaseSystem &fluid, const word &name, HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > &pDmdt)
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< volScalarField > trho
const tmp< volScalarField > sigmaWithContinuousPhase(const phaseModel &dispersedPhase) const
Return the surface tension coefficient between a given dispersed.
label k
Boltzmann constant.
const dimensionSet dimless
Dimensionless.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
bool writeData(Ostream &) const
Dummy write for regIOobject.
Lookup type of boundary radiation properties.
Definition: lookup.H:57
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:478
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
dimensionedScalar neg(const dimensionedScalar &ds)
autoPtr< populationBalanceModel > clone() const
Return clone.
Base class for coalescence models.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
constexpr const char *const group
Group name for atomic constants.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:40
Base class for breakup models which give a total breakup rate and a separate daughter size distributi...
Definition: breakupModel.H:51
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
Calculate the first temporal derivative.
Base class for binary breakup models which give the breakup rate between a sizeGroup pair directly...
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:30
word timeName
Definition: getTimeIndex.H:3
dimensionedScalar pos(const dimensionedScalar &ds)
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Base class for nucleation models.
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
static int compare(const Pair< word > &a, const Pair< word > &b)
Compare Pairs.
Definition: PairI.H:24
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:69
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
Calculate the matrix for the first temporal derivative.
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Calculate the field for explicit evaluation of implicit and explicit sources.
Constant dispersed-phase particle diameter model.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition: IOobject.C:456
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
Calculate the divergence of the given field.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/chef2/andy/OpenFOAM/release/v2406/OpenFOAM-v2406/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
OBJstream os(runTime.globalPath()/outputName)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:103
PIMPLE control class to supply convergence information/checks for the PIMPLE loop.
Definition: pimpleControl.H:51
Nothing to be read.
Automatically write from objectRegistry::writeObject()
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:68
messageStream Info
Information stream (stdout output on master, null elsewhere)
const phaseCompressibleTurbulenceModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
const scalar gamma
Definition: EEqn.H:9
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
IOobject newIOobject(const word &name, IOobjectOption ioOpt) const
Create an IOobject at the current time instance (timeName) with the specified options.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
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:39
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:696
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
const dimensionedScalar gamma(const label i, const dimensionedScalar &v) const
Return allocation coefficient.
Request registration (bool: true)
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:49
Do not request registration (bool: false)
Calculate the finiteVolume matrix for implicit and explicit sources.
label timeIndex
Definition: getTimeIndex.H:24
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const dimensionSet dimVelocity