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