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