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-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "populationBalanceModel.H"
27 #include "phaseSystem.H"
29 #include "coalescenceModel.H"
30 #include "breakupModel.H"
31 #include "binaryBreakupModel.H"
32 #include "driftModel.H"
33 #include "nucleationModel.H"
35 #include "fvmDdt.H"
36 #include "fvmDiv.H"
37 #include "fvmSup.H"
38 #include "fvcDdt.H"
39 #include "fvcDiv.H"
40 #include "fvcSup.H"
41 #include "shapeModel.H"
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 namespace diameterModels
48 {
50 }
51 }
52 
53 
54 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
55 
56 void Foam::diameterModels::populationBalanceModel::initialiseDmdtfs()
57 {
59  (
60  HashTable<const diameterModels::velocityGroup*>,
61  velocityGroupPtrs_,
62  iter1
63  )
64  {
65  const diameterModels::velocityGroup& velGrp1 = *iter1();
66 
68  (
69  HashTable<const diameterModels::velocityGroup*>,
70  velocityGroupPtrs_,
71  iter2
72  )
73  {
74  const diameterModels::velocityGroup& velGrp2 = *iter2();
75 
76  const phaseInterface interface(velGrp1.phase(), velGrp2.phase());
77 
78  if
79  (
80  &velGrp1 != &velGrp2
81  &&
82  !dmdtfs_.found(interface)
83  )
84  {
86  <diameterModels::populationBalanceModel>(interface);
87 
88  dmdtfs_.insert
89  (
90  interface,
91  new volScalarField
92  (
93  IOobject
94  (
96  (
97  typedName("dmdtf"),
98  interface.name()
99  ),
100  mesh().time().name(),
101  mesh()
102  ),
103  mesh(),
105  )
106  );
107  }
108  }
109  }
110 }
111 
112 
113 void Foam::diameterModels::populationBalanceModel::precompute()
114 {
115  forAll(coalescenceModels_, model)
116  {
117  coalescenceModels_[model].precompute();
118  }
119 
120  forAll(breakupModels_, model)
121  {
122  breakupModels_[model].precompute();
123 
124  breakupModels_[model].dsdPtr()->precompute();
125  }
126 
127  forAll(binaryBreakupModels_, model)
128  {
129  binaryBreakupModels_[model].precompute();
130  }
131 
132  forAll(drift_, model)
133  {
134  drift_[model].precompute();
135  }
136 
137  forAll(nucleation_, model)
138  {
139  nucleation_[model].precompute();
140  }
141 }
142 
143 
144 void Foam::diameterModels::populationBalanceModel::birthByCoalescence
145 (
146  const label j,
147  const label k
148 )
149 {
150  const sizeGroup& fj = sizeGroups()[j];
151  const sizeGroup& fk = sizeGroups()[k];
152 
153  dimensionedScalar Eta;
154  dimensionedScalar v = fj.x() + fk.x();
155 
156  for (label i = j; i < sizeGroups().size(); i++)
157  {
158  Eta = eta(i, v);
159 
160  if (Eta.value() == 0) continue;
161 
162  const sizeGroup& fi = sizeGroups()[i];
163 
164  if (j == k)
165  {
166  Sui_ =
167  0.5*fi.x()/(fj.x()*fk.x())*Eta
168  *coalescenceRate_()*fj*fj.phase()*fk*fk.phase();
169  }
170  else
171  {
172  Sui_ =
173  fi.x()/(fj.x()*fk.x())*Eta
174  *coalescenceRate_()*fj*fj.phase()*fk*fk.phase();
175  }
176 
177  Su_[i] += Sui_;
178 
179  const phaseInterface interfaceij(fi.phase(), fj.phase());
180 
181  if (dmdtfs_.found(interfaceij))
182  {
183  const scalar dmdtSign =
184  interfaceij.index(fi.phase()) == 0 ? +1 : -1;
185 
186  *dmdtfs_[interfaceij] += dmdtSign*fj.x()/v*Sui_*fj.phase().rho();
187  }
188 
189  const phaseInterface interfaceik(fi.phase(), fk.phase());
190 
191  if (dmdtfs_.found(interfaceik))
192  {
193  const scalar dmdtSign =
194  interfaceik.index(fi.phase()) == 0 ? +1 : -1;
195 
196  *dmdtfs_[interfaceik] += dmdtSign*fk.x()/v*Sui_*fk.phase().rho();
197  }
198 
199  sizeGroups_[i].shapeModelPtr()->addCoalescence(Sui_, fj, fk);
200  }
201 }
202 
203 
204 void Foam::diameterModels::populationBalanceModel::deathByCoalescence
205 (
206  const label i,
207  const label j
208 )
209 {
210  const sizeGroup& fi = sizeGroups()[i];
211  const sizeGroup& fj = sizeGroups()[j];
212 
213  Sp_[i] += coalescenceRate_()*fi.phase()*fj*fj.phase()/fj.x();
214 
215  if (i != j)
216  {
217  Sp_[j] += coalescenceRate_()*fj.phase()*fi*fi.phase()/fi.x();
218  }
219 }
220 
221 
222 void Foam::diameterModels::populationBalanceModel::birthByBreakup
223 (
224  const label k,
225  const label model
226 )
227 {
228  const sizeGroup& fk = sizeGroups()[k];
229 
230  for (label i = 0; i <= k; i++)
231  {
232  const sizeGroup& fi = sizeGroups()[i];
233 
234  Sui_ =
235  fi.x()*breakupModels_[model].dsdPtr()().nik(i, k)/fk.x()
236  *breakupRate_()*fk*fk.phase();
237 
238  Su_[i] += Sui_;
239 
240  const phaseInterface interface(fi.phase(), fk.phase());
241 
242  if (dmdtfs_.found(interface))
243  {
244  const scalar dmdtSign =
245  interface.index(fi.phase()) == 0 ? +1 : -1;
246 
247  *dmdtfs_[interface] += dmdtSign*Sui_*fk.phase().rho();
248  }
249 
250  sizeGroups_[i].shapeModelPtr()->addBreakup(Sui_, fk);
251  }
252 }
253 
254 
255 void Foam::diameterModels::populationBalanceModel::deathByBreakup(const label i)
256 {
257  Sp_[i] += breakupRate_()*sizeGroups()[i].phase();
258 }
259 
260 
261 void Foam::diameterModels::populationBalanceModel::calcDeltas()
262 {
263  forAll(sizeGroups(), i)
264  {
265  if (delta_[i].empty())
266  {
267  for (label j = 0; j <= sizeGroups().size() - 1; j++)
268  {
269  const sizeGroup& fj = sizeGroups()[j];
270 
271  delta_[i].append
272  (
274  (
275  "delta",
276  dimVolume,
277  v_[i+1].value() - v_[i].value()
278  )
279  );
280 
281  if
282  (
283  v_[i].value() < 0.5*fj.x().value()
284  &&
285  0.5*fj.x().value() < v_[i+1].value()
286  )
287  {
288  delta_[i][j] = mag(0.5*fj.x() - v_[i]);
289  }
290  else if (0.5*fj.x().value() <= v_[i].value())
291  {
292  delta_[i][j].value() = 0;
293  }
294  }
295  }
296  }
297 }
298 
299 
300 void Foam::diameterModels::populationBalanceModel::birthByBinaryBreakup
301 (
302  const label i,
303  const label j
304 )
305 {
306  const sizeGroup& fi = sizeGroups()[i];
307  const sizeGroup& fj = sizeGroups()[j];
308 
309  const volScalarField Su(binaryBreakupRate_()*fj*fj.phase());
310 
311  Sui_ = fi.x()*delta_[i][j]/fj.x()*Su;
312 
313  Su_[i] += Sui_;
314 
315  sizeGroups_[i].shapeModelPtr()->addBreakup(Sui_, fj);
316 
317  const phaseInterface interfaceij(fi.phase(), fj.phase());
318 
319  if (dmdtfs_.found(interfaceij))
320  {
321  const scalar dmdtSign =
322  interfaceij.index(fi.phase()) == 0 ? +1 : -1;
323 
324  *dmdtfs_[interfaceij] += dmdtSign*Sui_*fj.phase().rho();
325  }
326 
327  dimensionedScalar Eta;
328  dimensionedScalar v = fj.x() - fi.x();
329 
330  for (label k = 0; k <= j; k++)
331  {
332  Eta = eta(k, v);
333 
334  if (Eta.value() == 0) continue;
335 
336  const sizeGroup& fk = sizeGroups()[k];
337 
338  volScalarField& Suk = Sui_;
339 
340  Suk = fk.x()*delta_[i][j]*Eta/fj.x()*Su;
341 
342  Su_[k] += Suk;
343 
344  const phaseInterface interfacekj(fk.phase(), fj.phase());
345 
346  if (dmdtfs_.found(interfacekj))
347  {
348  const scalar dmdtSign =
349  interfacekj.index(fk.phase()) == 0 ? +1 : -1;
350 
351  *dmdtfs_[interfacekj] += dmdtSign*Suk*fj.phase().rho();
352  }
353 
354  sizeGroups_[k].shapeModelPtr()->addBreakup(Suk, fj);
355  }
356 }
357 
358 
359 void Foam::diameterModels::populationBalanceModel::deathByBinaryBreakup
360 (
361  const label j,
362  const label i
363 )
364 {
365  Sp_[i] += sizeGroups()[i].phase()*binaryBreakupRate_()*delta_[j][i];
366 }
367 
368 
369 void Foam::diameterModels::populationBalanceModel::drift
370 (
371  const label i,
372  driftModel& model
373 )
374 {
375  model.addToDriftRate(driftRate_(), i);
376 
377  const sizeGroup& fp = sizeGroups()[i];
378 
379  if (i < sizeGroups().size() - 1)
380  {
381  const sizeGroup& fe = sizeGroups()[i+1];
382  volScalarField& Sue = Sui_;
383 
384  Sp_[i] += 1/(fe.x() - fp.x())*pos(driftRate_())*driftRate_();
385 
386  Sue =
387  fe.x()/(fp.x()*(fe.x() - fp.x()))*pos(driftRate_())*driftRate_()*fp;
388 
389  Su_[i+1] += Sue;
390 
391  const phaseInterface interfacepe(fp.phase(), fe.phase());
392 
393  if (dmdtfs_.found(interfacepe))
394  {
395  const scalar dmdtSign =
396  interfacepe.index(fp.phase()) == 0 ? +1 : -1;
397 
398  *dmdtfs_[interfacepe] -= dmdtSign*Sue*fp.phase().rho();
399  }
400 
401  sizeGroups_[i+1].shapeModelPtr()->addDrift(Sue, fp, model);
402  }
403 
404  if (i == sizeGroups().size() - 1)
405  {
406  Sp_[i] -= pos(driftRate_())*driftRate_()/fp.x();
407  }
408 
409  if (i > 0)
410  {
411  const sizeGroup& fw = sizeGroups()[i-1];
412  volScalarField& Suw = Sui_;
413 
414  Sp_[i] += 1/(fw.x() - fp.x())*neg(driftRate_())*driftRate_();
415 
416  Suw =
417  fw.x()/(fp.x()*(fw.x() - fp.x()))*neg(driftRate_())*driftRate_()*fp;
418 
419  Su_[i-1] += Suw;
420 
421  const phaseInterface interfacepw(fp.phase(), fw.phase());
422 
423  if (dmdtfs_.found(interfacepw))
424  {
425  const scalar dmdtSign =
426  interfacepw.index(fp.phase()) == 0 ? +1 : -1;
427 
428  *dmdtfs_[interfacepw] -= dmdtSign*Suw*fp.phase().rho();
429  }
430 
431  sizeGroups_[i-1].shapeModelPtr()->addDrift(Suw, fp, model);
432  }
433 
434  if (i == 0)
435  {
436  Sp_[i] -= neg(driftRate_())*driftRate_()/fp.x();
437  }
438 }
439 
440 
441 void Foam::diameterModels::populationBalanceModel::nucleation
442 (
443  const label i,
444  nucleationModel& model
445 )
446 {
447  const sizeGroup& fi = sizeGroups()[i];
448 
449  model.addToNucleationRate(nucleationRate_(), i);
450 
451  Sui_ = fi.x()*nucleationRate_();
452 
453  Su_[i] += Sui_;
454 
455  sizeGroups_[i].shapeModelPtr()->addNucleation(Sui_, fi, model);
456 }
457 
458 
459 void Foam::diameterModels::populationBalanceModel::sources()
460 {
461  forAll(sizeGroups(), i)
462  {
463  sizeGroups_[i].shapeModelPtr()->reset();
464  Su_[i] = Zero;
465  Sp_[i] = Zero;
466  }
467 
468  forAllIter(phaseSystem::dmdtfTable, dmdtfs_, pDmdtIter)
469  {
470  *pDmdtIter() = Zero;
471  }
472 
473  forAll(coalescencePairs_, coalescencePairi)
474  {
475  label i = coalescencePairs_[coalescencePairi].first();
476  label j = coalescencePairs_[coalescencePairi].second();
477 
478  coalescenceRate_() = Zero;
479 
480  forAll(coalescenceModels_, model)
481  {
482  coalescenceModels_[model].addToCoalescenceRate
483  (
484  coalescenceRate_(),
485  i,
486  j
487  );
488  }
489 
490  birthByCoalescence(i, j);
491 
492  deathByCoalescence(i, j);
493  }
494 
495  forAll(binaryBreakupPairs_, binaryBreakupPairi)
496  {
497  label i = binaryBreakupPairs_[binaryBreakupPairi].first();
498  label j = binaryBreakupPairs_[binaryBreakupPairi].second();
499 
500  binaryBreakupRate_() = Zero;
501 
502  forAll(binaryBreakupModels_, model)
503  {
504  binaryBreakupModels_[model].addToBinaryBreakupRate
505  (
506  binaryBreakupRate_(),
507  j,
508  i
509  );
510  }
511 
512  birthByBinaryBreakup(j, i);
513 
514  deathByBinaryBreakup(j, i);
515  }
516 
517  forAll(sizeGroups(), i)
518  {
519  forAll(breakupModels_, model)
520  {
521  breakupModels_[model].setBreakupRate(breakupRate_(), i);
522 
523  birthByBreakup(i, model);
524 
525  deathByBreakup(i);
526  }
527 
528  forAll(drift_, model)
529  {
530  driftRate_() = Zero;
531 
532  drift(i, drift_[model]);
533  }
534 
535  forAll(nucleation_, model)
536  {
537  nucleationRate_() = Zero;
538 
539  nucleation(i, nucleation_[model]);
540  }
541  }
542 }
543 
544 
545 void Foam::diameterModels::populationBalanceModel::correctDilatationError()
546 {
547  forAllIter
548  (
549  HashTable<volScalarField>,
550  dilatationErrors_,
551  iter
552  )
553  {
554  volScalarField& dilatationError = iter();
555  const word& phaseName = iter.key();
556  const phaseModel& phase = fluid_.phases()[phaseName];
557  const velocityGroup& velGroup = *velocityGroupPtrs_[phaseName];
558  const volScalarField& alpha = phase;
559  const volScalarField& rho = phase.rho();
560 
561  dilatationError =
562  fvc::ddt(alpha) + fvc::div(phase.alphaPhi())
563  - (fluid_.fvModels().source(alpha, rho) & rho)/rho;
564 
565  forAll(velGroup.sizeGroups(), i)
566  {
567  const sizeGroup& fi = velGroup.sizeGroups()[i];
568 
569  dilatationError -= Su_[fi.i()] - fvc::Sp(Sp_[fi.i()], fi);
570  }
571  }
572 }
573 
574 
575 void Foam::diameterModels::populationBalanceModel::calcAlphas()
576 {
577  alphas_() = Zero;
578 
579  forAllConstIter(HashTable<const velocityGroup*>, velocityGroupPtrs_, iter)
580  {
581  const phaseModel& phase = iter()->phase();
582 
583  alphas_() += max(phase, phase.residualAlpha());
584  }
585 }
586 
587 
589 Foam::diameterModels::populationBalanceModel::calcDsm()
590 {
591  tmp<volScalarField> tInvDsm
592  (
594  (
595  "invDsm",
596  mesh_,
598  )
599  );
600 
601  volScalarField& invDsm = tInvDsm.ref();
602 
603  forAllConstIter(HashTable<const velocityGroup*>, velocityGroupPtrs_, iter)
604  {
605  const phaseModel& phase = iter()->phase();
606 
607  invDsm += max(phase, phase.residualAlpha())/(phase.d()*alphas_());
608  }
609 
610  return 1/tInvDsm;
611 }
612 
613 
614 void Foam::diameterModels::populationBalanceModel::calcVelocity()
615 {
616  U_() = Zero;
617 
618  forAllConstIter(HashTable<const velocityGroup*>, velocityGroupPtrs_, iter)
619  {
620  const phaseModel& phase = iter()->phase();
621 
622  U_() += max(phase, phase.residualAlpha())*phase.U()/alphas_();
623  }
624 }
625 
626 bool Foam::diameterModels::populationBalanceModel::updateSources()
627 {
628  const bool result = sourceUpdateCounter_ % sourceUpdateInterval() == 0;
629 
630  ++ sourceUpdateCounter_;
631 
632  return result;
633 }
634 
635 
636 template<class EtaType, class VType>
637 EtaType Foam::diameterModels::populationBalanceModel::eta
638 (
639  const label i,
640  const VType& v
641 ) const
642 {
643  const label n = sizeGroups().size();
644 
645  static const dimensionedScalar rootVSmallV(dimVolume, rootVSmall);
646  static const dimensionedScalar z(dimless, scalar(0));
647 
648  const dimensionedScalar& x0 =
649  i > 0 ? sizeGroups()[i - 1].x() : rootVSmallV;
650  const dimensionedScalar& xi = sizeGroups()[i].x();
651  const dimensionedScalar& x1 =
652  i < n - 1 ? sizeGroups()[i + 1].x() : rootVSmallV;
653 
654  return max(min((v - x0)/(xi - x0), (x1 - v)/(x1 - xi)), z);
655 }
656 
657 
658 template<class EtaType, class VType>
659 EtaType Foam::diameterModels::populationBalanceModel::etaV
660 (
661  const label i,
662  const VType& v
663 ) const
664 {
665  const label n = sizeGroups().size();
666 
667  static const dimensionedScalar rootVSmallInvV(inv(dimVolume), rootVSmall);
668  static const dimensionedScalar rootVGreatInvV(inv(dimVolume), rootVGreat);
669  static const dimensionedScalar z(dimless, scalar(0));
670 
671  const dimensionedScalar x0 =
672  i > 0 ? 1/sizeGroups()[i - 1].x() : rootVGreatInvV;
673  const dimensionedScalar xi = 1/sizeGroups()[i].x();
674  const dimensionedScalar x1 =
675  i < n - 1 ? 1/sizeGroups()[i + 1].x() : rootVSmallInvV;
676 
677  const VType invV
678  (
679  1/min(max(v, sizeGroups().first().x()), sizeGroups().last().x())
680  );
681 
682  return max(min((invV - x0)/(xi - x0), (x1 - invV)/(x1 - xi)), z);
683 }
684 
685 
686 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
687 
689 (
690  const phaseSystem& fluid,
691  const word& name
692 )
693 :
695  (
696  IOobject
697  (
698  name,
699  fluid.time().constant(),
700  fluid.mesh()
701  )
702  ),
703  fluid_(fluid),
704  dmdtfs_(),
705  mesh_(fluid_.mesh()),
706  name_(name),
707  dict_
708  (
709  fluid_.subDict("populationBalanceCoeffs").subDict(name_)
710  ),
711  continuousPhase_
712  (
713  mesh_.lookupObject<phaseModel>
714  (
715  IOobject::groupName("alpha", dict_.lookup("continuousPhase"))
716  )
717  ),
718  sizeGroups_(),
719  v_(),
720  delta_(),
721  Su_(),
722  Sp_(),
723  Sui_
724  (
725  IOobject
726  (
727  "Sui",
728  mesh_.time().name(),
729  mesh_
730  ),
731  mesh_,
733  ),
734  coalescenceModels_
735  (
736  dict_.lookup("coalescenceModels"),
737  coalescenceModel::iNew(*this)
738  ),
739  coalescenceRate_(),
740  coalescencePairs_(),
741  breakupModels_
742  (
743  dict_.lookup("breakupModels"),
744  breakupModel::iNew(*this)
745  ),
746  breakupRate_(),
747  binaryBreakupModels_
748  (
749  dict_.lookup("binaryBreakupModels"),
750  binaryBreakupModel::iNew(*this)
751  ),
752  binaryBreakupRate_(),
753  binaryBreakupPairs_(),
754  drift_
755  (
756  dict_.lookup("driftModels"),
757  driftModel::iNew(*this)
758  ),
759  driftRate_(),
760  nucleation_
761  (
762  dict_.lookup("nucleationModels"),
763  nucleationModel::iNew(*this)
764  ),
765  nucleationRate_(),
766  alphas_(),
767  dsm_(),
768  U_(),
769  sourceUpdateCounter_(0)
770 {
771  groups::retrieve(*this, velocityGroupPtrs_, sizeGroups_);
772 
773  if (sizeGroups().size() < 3)
774  {
776  << "The populationBalance " << name_
777  << " requires a minimum number of three sizeGroups to be specified."
778  << exit(FatalError);
779  }
780 
781  // Create velocity-group dilatation errors
783  (
785  velocityGroupPtrs_,
786  iter
787  )
788  {
789  const velocityGroup& velGroup = *iter();
790 
791  dilatationErrors_.insert
792  (
793  velGroup.phase().name(),
795  (
796  IOobject
797  (
799  (
800  "dilatationError",
801  velGroup.phase().name()
802  ),
803  fluid_.time().name(),
804  mesh_
805  ),
806  mesh_,
808  )
809  );
810  }
811 
812  // Create size-group boundaries
813  v_.setSize(sizeGroups().size() + 1);
814  v_.set(0, new dimensionedScalar("v", sizeGroups()[0].x()));
815  for (label i = 1; i < sizeGroups().size(); ++ i)
816  {
817  v_.set
818  (
819  i,
821  (
822  "v",
823  (sizeGroups()[i-1].x() + sizeGroups()[i].x())/2
824  )
825  );
826  }
827  v_.set(v_.size() - 1, new dimensionedScalar("v", sizeGroups().last().x()));
828 
829  // ???
830  delta_.setSize(sizeGroups().size());
831  forAll(sizeGroups(), i)
832  {
833  delta_.set(i, new PtrList<dimensionedScalar>());
834  }
835 
836  // Create size-group source terms
837  Su_.setSize(sizeGroups().size());
838  Sp_.setSize(sizeGroups().size());
839  forAll(sizeGroups(), i)
840  {
841  Su_.set
842  (
843  i,
844  new volScalarField
845  (
846  IOobject("Su", fluid_.time().name(), mesh_),
847  mesh_,
849  )
850  );
851 
852  Sp_.set
853  (
854  i,
855  new volScalarField
856  (
857  IOobject("Sp", fluid_.time().name(), mesh_),
858  mesh_,
860  )
861  );
862  }
863 
864  this->initialiseDmdtfs();
865 
866  if (coalescenceModels_.size() != 0)
867  {
868  coalescenceRate_.set
869  (
870  new volScalarField
871  (
872  IOobject
873  (
874  "coalescenceRate",
875  mesh_.time().name(),
876  mesh_
877  ),
878  mesh_,
880  )
881  );
882 
883  forAll(sizeGroups(), i)
884  {
885  for (label j = 0; j <= i; j++)
886  {
887  coalescencePairs_.append(labelPair(i, j));
888  }
889  }
890  }
891 
892  if (breakupModels_.size() != 0)
893  {
894  breakupRate_.set
895  (
896  new volScalarField
897  (
898  IOobject
899  (
900  "breakupRate",
901  fluid_.time().name(),
902  mesh_
903  ),
904  mesh_,
906  )
907  );
908  }
909 
910  if (binaryBreakupModels_.size() != 0)
911  {
912  binaryBreakupRate_.set
913  (
914  new volScalarField
915  (
916  IOobject
917  (
918  "binaryBreakupRate",
919  fluid_.time().name(),
920  mesh_
921  ),
922  mesh_,
924  (
925  "binaryBreakupRate",
927  Zero
928  )
929  )
930  );
931 
932  calcDeltas();
933 
934  forAll(sizeGroups(), i)
935  {
936  label j = 0;
937 
938  while (delta_[j][i].value() != 0)
939  {
940  binaryBreakupPairs_.append(labelPair(i, j));
941  j++;
942  }
943  }
944  }
945 
946  if (drift_.size() != 0)
947  {
948  driftRate_.set
949  (
950  new volScalarField
951  (
952  IOobject
953  (
954  "driftRate",
955  fluid_.time().name(),
956  mesh_
957  ),
958  mesh_,
960  )
961  );
962  }
963 
964  if (nucleation_.size() != 0)
965  {
966  nucleationRate_.set
967  (
968  new volScalarField
969  (
970  IOobject
971  (
972  "nucleationRate",
973  fluid_.time().name(),
974  mesh_
975  ),
976  mesh_,
978  (
979  "nucleationRate",
981  Zero
982  )
983  )
984  );
985  }
986 
987  if (velocityGroupPtrs_.size() > 1)
988  {
989  alphas_.set
990  (
991  new volScalarField
992  (
993  IOobject
994  (
995  IOobject::groupName("alpha", this->name()),
996  fluid_.time().name(),
997  mesh_,
1000  ),
1001  mesh_,
1003  )
1004  );
1005 
1006  dsm_.set
1007  (
1008  new volScalarField
1009  (
1010  IOobject
1011  (
1012  IOobject::groupName("d", this->name()),
1013  fluid_.time().name(),
1014  mesh_,
1017  ),
1018  mesh_,
1020  )
1021  );
1022 
1023  U_.set
1024  (
1025  new volVectorField
1026  (
1027  IOobject
1028  (
1029  IOobject::groupName("U", this->name()),
1030  fluid_.time().name(),
1031  mesh_,
1034  ),
1035  mesh_,
1037  )
1038  );
1039  }
1040 }
1041 
1042 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1043 
1045 {}
1046 
1047 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1048 
1051 {
1052  notImplemented("populationBalance::clone() const");
1053  return autoPtr<populationBalanceModel>(nullptr);
1054 }
1055 
1056 
1058 {
1059  return os.good();
1060 }
1061 
1062 
1063 Foam::dimensionedScalar Foam::diameterModels::populationBalanceModel::eta
1064 (
1065  const label i,
1066  const dimensionedScalar& v
1067 ) const
1068 {
1069  return eta<dimensionedScalar>(i, v);
1070 }
1071 
1072 
1074 Foam::diameterModels::populationBalanceModel::eta
1075 (
1076  const label i,
1077  const volScalarField::Internal& v
1078 ) const
1079 {
1080  return eta<tmp<volScalarField::Internal>>(i, v);
1081 }
1082 
1083 
1084 Foam::dimensionedScalar Foam::diameterModels::populationBalanceModel::etaV
1085 (
1086  const label i,
1087  const dimensionedScalar& v
1088 ) const
1089 {
1090  return etaV<dimensionedScalar>(i, v);
1091 }
1092 
1093 
1095 Foam::diameterModels::populationBalanceModel::etaV
1096 (
1097  const label i,
1098  const volScalarField::Internal& v
1099 ) const
1100 {
1101  return etaV<tmp<volScalarField::Internal>>(i, v);
1102 }
1103 
1104 
1107 (
1108  const phaseModel& dispersedPhase
1109 ) const
1110 {
1111  return phaseInterface(dispersedPhase, continuousPhase_).sigma();
1112 }
1113 
1114 
1117 {
1118  return
1120  (
1121  continuousPhase_.name()
1122  );
1123 }
1124 
1125 
1127 {
1128  if (!solveOnFinalIterOnly() || fluid_.pimple().finalIter())
1129  {
1130  const label nCorr = this->nCorr();
1131  const scalar tolerance =
1132  mesh_.solution().solverDict(name_).lookup<scalar>("tolerance");
1133 
1134  const bool updateSrc = updateSources();
1135 
1136  if (nCorr > 0 && updateSrc)
1137  {
1138  precompute();
1139  }
1140 
1141  int iCorr = 0;
1142  scalar maxInitialResidual = 1;
1143 
1144  while (++iCorr <= nCorr && maxInitialResidual > tolerance)
1145  {
1146  Info<< "populationBalance "
1147  << this->name()
1148  << ": Iteration "
1149  << iCorr
1150  << endl;
1151 
1152  if (updateSrc)
1153  {
1154  sources();
1155  }
1156 
1157  correctDilatationError();
1158 
1159  maxInitialResidual = 0;
1160 
1161  forAll(sizeGroups(), i)
1162  {
1163  sizeGroup& fi = sizeGroups_[i];
1164  const phaseModel& phase = fi.phase();
1165  const volScalarField& alpha = phase;
1166  const volScalarField& rho = phase.rho();
1167  const volScalarField& dilatationError =
1168  dilatationErrors_[phase.name()];
1169 
1170  fvScalarMatrix sizeGroupEqn
1171  (
1172  fvm::ddt(alpha, fi)
1173  + fvm::div(phase.alphaPhi(), fi)
1174  + fvm::Sp(-(1 - small)*dilatationError, fi)
1175  + fvm::SuSp(-small*dilatationError, fi)
1176  ==
1177  fvc::Su(Su_[i], fi)
1178  - fvm::Sp(Sp_[i], fi)
1179  + fluid_.fvModels().source(alpha, rho, fi)/rho
1180  - correction
1181  (
1182  fvm::Sp
1183  (
1184  max(phase.residualAlpha() - alpha, scalar(0))
1185  /this->mesh().time().deltaT(),
1186  fi
1187  )
1188  )
1189  );
1190 
1191  sizeGroupEqn.relax();
1192 
1193  fluid_.fvConstraints().constrain(sizeGroupEqn);
1194 
1195  maxInitialResidual = max
1196  (
1197  sizeGroupEqn.solve().initialResidual(),
1198  maxInitialResidual
1199  );
1200 
1201  fluid_.fvConstraints().constrain(fi);
1202  }
1203  }
1204 
1205  volScalarField fAlpha0
1206  (
1207  sizeGroups().first()*sizeGroups().first().phase()
1208  );
1209 
1210  volScalarField fAlphaN
1211  (
1212  sizeGroups().last()*sizeGroups().last().phase()
1213  );
1214 
1215  Info<< this->name() << " sizeGroup phase fraction first, last = "
1216  << fAlpha0.weightedAverage(this->mesh().V()).value()
1217  << ' ' << fAlphaN.weightedAverage(this->mesh().V()).value()
1218  << endl;
1219  }
1220 }
1221 
1222 
1224 {
1225  if (velocityGroupPtrs_.size() > 1)
1226  {
1227  calcAlphas();
1228  dsm_() = calcDsm();
1229  calcVelocity();
1230  }
1231 }
1232 
1233 
1234 // ************************************************************************* //
label k
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &) const
Calculate and return weighted average.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
An STL-conforming hash table.
Definition: HashTable.H:127
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const Time & time() const
Return time.
Definition: IOobject.C:318
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:167
const word & name() const
Return name.
Definition: IOobject.H:310
static word groupName(Name name, const word &group)
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:330
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
const phaseModel & phase() const
Return the phase.
Base class for binary breakup models that provide a breakup rate between a size class pair directly,...
Base class for breakup models which provide a total breakup rate and a separate daughter size distrib...
Definition: breakupModel.H:56
Base class for coalescence models.
Constant dispersed-phase particle diameter model.
Base class for drift models.
Definition: driftModel.H:54
Base class for nucleation models.
static void retrieve(const populationBalanceModel &popBal, HashTable< const velocityGroup * > &velGroupPtrs, UPtrList< sizeGroup > &szGroupPtrs)
Retrieve the pointers.
Return a pointer to a new populationBalanceModel object created on.
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
bool writeData(Ostream &) const
Dummy write for regIOobject.
autoPtr< populationBalanceModel > clone() const
Return clone.
populationBalanceModel(const phaseSystem &fluid, const word &name)
Construct for a fluid.
const UPtrList< sizeGroup > & sizeGroups() const
Return the size groups belonging to this populationBalance.
const tmp< volScalarField > sigmaWithContinuousPhase(const phaseModel &dispersedPhase) const
Return the surface tension coefficient between a given dispersed.
const fvMesh & mesh() const
Return reference to the mesh.
const phaseCompressibleMomentumTransportModel & continuousTurbulence() const
Return reference to momentumTransport model of the continuous phase.
void solve()
Solve the population balance equation.
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:102
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:37
Computes the Sauter mean diameter based on a user specified size distribution, defined in terms of si...
Definition: velocityGroup.H:87
const word & name() const
Return const reference to name.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:603
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
Templated abstract base class for multiphase compressible turbulence models.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
tmp< volScalarField > sigma() const
Surface tension coefficient.
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
Definition: phaseModel.C:169
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:145
virtual const volScalarField & rho() const =0
Return the density field.
virtual tmp< surfaceScalarField > alphaPhi() const =0
Return the volumetric flux of the phase.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
void validateMassTransfer(const phaseInterface &interface) const
Check that mass transfer is supported across the given interface.
HashPtrTable< volScalarField, phaseInterfaceKey, phaseInterfaceKey::hash > dmdtfTable
Definition: phaseSystem.H:95
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:55
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define notImplemented(functionName)
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:371
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the field for explicit evaluation of implicit and explicit sources.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for implicit and explicit sources.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
defineTypeNameAndDebug(constant, 0)
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
dimensionedScalar pos(const dimensionedScalar &ds)
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
messageStream Info
const dimensionSet dimLength
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
const dimensionSet dimDensity
const dimensionSet dimVolume
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
dimensioned< scalar > mag(const dimensioned< Type > &)
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:176
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
dimensionedScalar neg(const dimensionedScalar &ds)
const dimensionSet dimVelocity
error FatalError
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.