All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2022 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"
28 #include "coalescenceModel.H"
29 #include "breakupModel.H"
30 #include "binaryBreakupModel.H"
31 #include "driftModel.H"
32 #include "nucleationModel.H"
33 #include "surfaceTensionModel.H"
34 #include "fvmDdt.H"
35 #include "fvmDiv.H"
36 #include "fvmSup.H"
37 #include "fvcDdt.H"
38 #include "fvcDiv.H"
39 #include "fvcSup.H"
40 #include "shapeModel.H"
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 namespace diameterModels
47 {
48  defineTypeNameAndDebug(populationBalanceModel, 0);
49 }
50 }
51 
52 
53 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
54 
55 void Foam::diameterModels::populationBalanceModel::registerVelocityGroups()
56 {
57  forAll(fluid_.phases(), phasei)
58  {
59  if (isA<velocityGroup>(fluid_.phases()[phasei].dPtr()()))
60  {
61  const velocityGroup& velGroup =
62  refCast<const velocityGroup>(fluid_.phases()[phasei].dPtr()());
63 
64  if (velGroup.popBalName() == this->name())
65  {
66  velocityGroupPtrs_.insert(velGroup.phase().name(), &velGroup);
67 
68  dilatationErrors_.insert
69  (
70  velGroup.phase().name(),
72  (
73  IOobject
74  (
76  (
77  "dilatationError",
78  velGroup.phase().name()
79  ),
80  fluid_.time().timeName(),
81  mesh_
82  ),
83  mesh_,
84  dimensionedScalar(inv(dimTime), 0)
85  )
86  );
87 
88  forAll(velGroup.sizeGroups(), i)
89  {
90  this->registerSizeGroups
91  (
92  const_cast<sizeGroup&>(velGroup.sizeGroups()[i])
93  );
94  }
95  }
96  }
97  }
98 }
99 
100 
101 void Foam::diameterModels::populationBalanceModel::registerSizeGroups
102 (
103  sizeGroup& group
104 )
105 {
106  if
107  (
108  sizeGroups().size() != 0
109  &&
110  group.x().value() <= sizeGroups().last().x().value()
111  )
112  {
114  << "Size groups must be entered according to their representative"
115  << " size"
116  << exit(FatalError);
117  }
118 
119  sizeGroups_.resize(sizeGroups().size() + 1);
120  sizeGroups_.set(sizeGroups().size() - 1, &group);
121 
122  if (sizeGroups().size() == 1)
123  {
124  v_.append
125  (
127  (
128  "v",
129  sizeGroups().last().x()
130  )
131  );
132 
133  v_.append
134  (
136  (
137  "v",
138  sizeGroups().last().x()
139  )
140  );
141  }
142  else
143  {
144  v_.last() =
145  0.5
146  *(
147  sizeGroups()[sizeGroups().size()-2].x()
148  + sizeGroups().last().x()
149  );
150 
151  v_.append
152  (
154  (
155  "v",
156  sizeGroups().last().x()
157  )
158  );
159  }
160 
161  delta_.append(new PtrList<dimensionedScalar>());
162 
163  Su_.append
164  (
165  new volScalarField
166  (
167  IOobject
168  (
169  "Su",
170  fluid_.time().timeName(),
171  mesh_
172  ),
173  mesh_,
174  dimensionedScalar(inv(dimTime), 0)
175  )
176  );
177 
178  Sp_.append
179  (
180  new volScalarField
181  (
182  IOobject
183  (
184  "Sp",
185  fluid_.time().timeName(),
186  mesh_
187  ),
188  mesh_,
189  dimensionedScalar(inv(dimTime), 0)
190  )
191  );
192 }
193 
194 
195 void Foam::diameterModels::populationBalanceModel::precompute()
196 {
197  forAll(coalescenceModels_, model)
198  {
199  coalescenceModels_[model].precompute();
200  }
201 
202  forAll(breakupModels_, model)
203  {
204  breakupModels_[model].precompute();
205 
206  breakupModels_[model].dsdPtr()->precompute();
207  }
208 
209  forAll(binaryBreakupModels_, model)
210  {
211  binaryBreakupModels_[model].precompute();
212  }
213 
214  forAll(drift_, model)
215  {
216  drift_[model].precompute();
217  }
218 
219  forAll(nucleation_, model)
220  {
221  nucleation_[model].precompute();
222  }
223 }
224 
225 
226 void Foam::diameterModels::populationBalanceModel::birthByCoalescence
227 (
228  const label j,
229  const label k
230 )
231 {
232  const sizeGroup& fj = sizeGroups()[j];
233  const sizeGroup& fk = sizeGroups()[k];
234 
235  dimensionedScalar Eta;
236  dimensionedScalar v = fj.x() + fk.x();
237 
238  for (label i = j; i < sizeGroups().size(); i++)
239  {
240  Eta = eta(i, v);
241 
242  if (Eta.value() == 0) continue;
243 
244  const sizeGroup& fi = sizeGroups()[i];
245 
246  if (j == k)
247  {
248  Sui_ =
249  0.5*fi.x()/(fj.x()*fk.x())*Eta
250  *coalescenceRate_()*fj*fj.phase()*fk*fk.phase();
251  }
252  else
253  {
254  Sui_ =
255  fi.x()/(fj.x()*fk.x())*Eta
256  *coalescenceRate_()*fj*fj.phase()*fk*fk.phase();
257  }
258 
259  Su_[i] += Sui_;
260 
261  const phaseInterface interfaceij(fi.phase(), fj.phase());
262 
263  if (pDmdt_.found(interfaceij))
264  {
265  const scalar dmdtSign =
266  interfaceij.index(fi.phase()) == 0 ? +1 : -1;
267 
268  *pDmdt_[interfaceij] += dmdtSign*fj.x()/v*Sui_*fj.phase().rho();
269  }
270 
271  const phaseInterface interfaceik(fi.phase(), fk.phase());
272 
273  if (pDmdt_.found(interfaceik))
274  {
275  const scalar dmdtSign =
276  interfaceik.index(fi.phase()) == 0 ? +1 : -1;
277 
278  *pDmdt_[interfaceik] += dmdtSign*fk.x()/v*Sui_*fk.phase().rho();
279  }
280 
281  sizeGroups_[i].shapeModelPtr()->addCoalescence(Sui_, fj, fk);
282  }
283 }
284 
285 
286 void Foam::diameterModels::populationBalanceModel::deathByCoalescence
287 (
288  const label i,
289  const label j
290 )
291 {
292  const sizeGroup& fi = sizeGroups()[i];
293  const sizeGroup& fj = sizeGroups()[j];
294 
295  Sp_[i] += coalescenceRate_()*fi.phase()*fj*fj.phase()/fj.x();
296 
297  if (i != j)
298  {
299  Sp_[j] += coalescenceRate_()*fj.phase()*fi*fi.phase()/fi.x();
300  }
301 }
302 
303 
304 void Foam::diameterModels::populationBalanceModel::birthByBreakup
305 (
306  const label k,
307  const label model
308 )
309 {
310  const sizeGroup& fk = sizeGroups()[k];
311 
312  for (label i = 0; i <= k; i++)
313  {
314  const sizeGroup& fi = sizeGroups()[i];
315 
316  Sui_ =
317  fi.x()*breakupModels_[model].dsdPtr()().nik(i, k)/fk.x()
318  *breakupRate_()*fk*fk.phase();
319 
320  Su_[i] += Sui_;
321 
322  const phaseInterface interface(fi.phase(), fk.phase());
323 
324  if (pDmdt_.found(interface))
325  {
326  const scalar dmdtSign =
327  interface.index(fi.phase()) == 0 ? +1 : -1;
328 
329  *pDmdt_[interface] += dmdtSign*Sui_*fk.phase().rho();
330  }
331 
332  sizeGroups_[i].shapeModelPtr()->addBreakup(Sui_, fk);
333  }
334 }
335 
336 
337 void Foam::diameterModels::populationBalanceModel::deathByBreakup(const label i)
338 {
339  Sp_[i] += breakupRate_()*sizeGroups()[i].phase();
340 }
341 
342 
343 void Foam::diameterModels::populationBalanceModel::calcDeltas()
344 {
345  forAll(sizeGroups(), i)
346  {
347  if (delta_[i].empty())
348  {
349  for (label j = 0; j <= sizeGroups().size() - 1; j++)
350  {
351  const sizeGroup& fj = sizeGroups()[j];
352 
353  delta_[i].append
354  (
356  (
357  "delta",
358  dimVolume,
359  v_[i+1].value() - v_[i].value()
360  )
361  );
362 
363  if
364  (
365  v_[i].value() < 0.5*fj.x().value()
366  &&
367  0.5*fj.x().value() < v_[i+1].value()
368  )
369  {
370  delta_[i][j] = mag(0.5*fj.x() - v_[i]);
371  }
372  else if (0.5*fj.x().value() <= v_[i].value())
373  {
374  delta_[i][j].value() = 0;
375  }
376  }
377  }
378  }
379 }
380 
381 
382 void Foam::diameterModels::populationBalanceModel::birthByBinaryBreakup
383 (
384  const label i,
385  const label j
386 )
387 {
388  const sizeGroup& fi = sizeGroups()[i];
389  const sizeGroup& fj = sizeGroups()[j];
390 
391  const volScalarField Su(binaryBreakupRate_()*fj*fj.phase());
392 
393  Sui_ = fi.x()*delta_[i][j]/fj.x()*Su;
394 
395  Su_[i] += Sui_;
396 
397  sizeGroups_[i].shapeModelPtr()->addBreakup(Sui_, fj);
398 
399  const phaseInterface interfaceij(fi.phase(), fj.phase());
400 
401  if (pDmdt_.found(interfaceij))
402  {
403  const scalar dmdtSign =
404  interfaceij.index(fi.phase()) == 0 ? +1 : -1;
405 
406  *pDmdt_[interfaceij] += dmdtSign*Sui_*fj.phase().rho();
407  }
408 
409  dimensionedScalar Eta;
410  dimensionedScalar v = fj.x() - fi.x();
411 
412  for (label k = 0; k <= j; k++)
413  {
414  Eta = eta(k, v);
415 
416  if (Eta.value() == 0) continue;
417 
418  const sizeGroup& fk = sizeGroups()[k];
419 
420  volScalarField& Suk = Sui_;
421 
422  Suk = fk.x()*delta_[i][j]*Eta/fj.x()*Su;
423 
424  Su_[k] += Suk;
425 
426  const phaseInterface interfacekj(fk.phase(), fj.phase());
427 
428  if (pDmdt_.found(interfacekj))
429  {
430  const scalar dmdtSign =
431  interfacekj.index(fk.phase()) == 0 ? +1 : -1;
432 
433  *pDmdt_[interfacekj] += dmdtSign*Suk*fj.phase().rho();
434  }
435 
436  sizeGroups_[k].shapeModelPtr()->addBreakup(Suk, fj);
437  }
438 }
439 
440 
441 void Foam::diameterModels::populationBalanceModel::deathByBinaryBreakup
442 (
443  const label j,
444  const label i
445 )
446 {
447  Sp_[i] += sizeGroups()[i].phase()*binaryBreakupRate_()*delta_[j][i];
448 }
449 
450 
451 void Foam::diameterModels::populationBalanceModel::drift
452 (
453  const label i,
454  driftModel& model
455 )
456 {
457  model.addToDriftRate(driftRate_(), i);
458 
459  const sizeGroup& fp = sizeGroups()[i];
460 
461  if (i < sizeGroups().size() - 1)
462  {
463  const sizeGroup& fe = sizeGroups()[i+1];
464  volScalarField& Sue = Sui_;
465 
466  Sp_[i] += 1/(fe.x() - fp.x())*pos(driftRate_())*driftRate_();
467 
468  Sue =
469  fe.x()/(fp.x()*(fe.x() - fp.x()))*pos(driftRate_())*driftRate_()*fp;
470 
471  Su_[i+1] += Sue;
472 
473  const phaseInterface interfacepe(fp.phase(), fe.phase());
474 
475  if (pDmdt_.found(interfacepe))
476  {
477  const scalar dmdtSign =
478  interfacepe.index(fp.phase()) == 0 ? +1 : -1;
479 
480  *pDmdt_[interfacepe] -= dmdtSign*Sue*fp.phase().rho();
481  }
482 
483  sizeGroups_[i+1].shapeModelPtr()->addDrift(Sue, fp, model);
484  }
485 
486  if (i == sizeGroups().size() - 1)
487  {
488  Sp_[i] -= pos(driftRate_())*driftRate_()/fp.x();
489  }
490 
491  if (i > 0)
492  {
493  const sizeGroup& fw = sizeGroups()[i-1];
494  volScalarField& Suw = Sui_;
495 
496  Sp_[i] += 1/(fw.x() - fp.x())*neg(driftRate_())*driftRate_();
497 
498  Suw =
499  fw.x()/(fp.x()*(fw.x() - fp.x()))*neg(driftRate_())*driftRate_()*fp;
500 
501  Su_[i-1] += Suw;
502 
503  const phaseInterface interfacepw(fp.phase(), fw.phase());
504 
505  if (pDmdt_.found(interfacepw))
506  {
507  const scalar dmdtSign =
508  interfacepw.index(fp.phase()) == 0 ? +1 : -1;
509 
510  *pDmdt_[interfacepw] -= dmdtSign*Suw*fp.phase().rho();
511  }
512 
513  sizeGroups_[i-1].shapeModelPtr()->addDrift(Suw, fp, model);
514  }
515 
516  if (i == 0)
517  {
518  Sp_[i] -= neg(driftRate_())*driftRate_()/fp.x();
519  }
520 }
521 
522 
523 void Foam::diameterModels::populationBalanceModel::nucleation
524 (
525  const label i,
526  nucleationModel& model
527 )
528 {
529  const sizeGroup& fi = sizeGroups()[i];
530 
531  model.addToNucleationRate(nucleationRate_(), i);
532 
533  Sui_ = fi.x()*nucleationRate_();
534 
535  Su_[i] += Sui_;
536 
537  sizeGroups_[i].shapeModelPtr()->addNucleation(Sui_, fi, model);
538 }
539 
540 
541 void Foam::diameterModels::populationBalanceModel::sources()
542 {
543  forAll(sizeGroups(), i)
544  {
545  sizeGroups_[i].shapeModelPtr()->reset();
546  Su_[i] = Zero;
547  Sp_[i] = Zero;
548  }
549 
550  forAllIter(phaseSystem::dmdtfTable, pDmdt_, pDmdtIter)
551  {
552  *pDmdtIter() = Zero;
553  }
554 
555  forAll(coalescencePairs_, coalescencePairi)
556  {
557  label i = coalescencePairs_[coalescencePairi].first();
558  label j = coalescencePairs_[coalescencePairi].second();
559 
560  coalescenceRate_() = Zero;
561 
562  forAll(coalescenceModels_, model)
563  {
564  coalescenceModels_[model].addToCoalescenceRate
565  (
566  coalescenceRate_(),
567  i,
568  j
569  );
570  }
571 
572  birthByCoalescence(i, j);
573 
574  deathByCoalescence(i, j);
575  }
576 
577  forAll(binaryBreakupPairs_, binaryBreakupPairi)
578  {
579  label i = binaryBreakupPairs_[binaryBreakupPairi].first();
580  label j = binaryBreakupPairs_[binaryBreakupPairi].second();
581 
582  binaryBreakupRate_() = Zero;
583 
584  forAll(binaryBreakupModels_, model)
585  {
586  binaryBreakupModels_[model].addToBinaryBreakupRate
587  (
588  binaryBreakupRate_(),
589  j,
590  i
591  );
592  }
593 
594  birthByBinaryBreakup(j, i);
595 
596  deathByBinaryBreakup(j, i);
597  }
598 
599  forAll(sizeGroups(), i)
600  {
601  forAll(breakupModels_, model)
602  {
603  breakupModels_[model].setBreakupRate(breakupRate_(), i);
604 
605  birthByBreakup(i, model);
606 
607  deathByBreakup(i);
608  }
609 
610  forAll(drift_, model)
611  {
612  driftRate_() = Zero;
613 
614  drift(i, drift_[model]);
615  }
616 
617  forAll(nucleation_, model)
618  {
619  nucleationRate_() = Zero;
620 
621  nucleation(i, nucleation_[model]);
622  }
623  }
624 }
625 
626 
627 void Foam::diameterModels::populationBalanceModel::correctDilatationError()
628 {
629  forAllIter
630  (
631  HashTable<volScalarField>,
632  dilatationErrors_,
633  iter
634  )
635  {
636  volScalarField& dilatationError = iter();
637  const word& phaseName = iter.key();
638  const phaseModel& phase = fluid_.phases()[phaseName];
639  const velocityGroup& velGroup = *velocityGroupPtrs_[phaseName];
640  const volScalarField& alpha = phase;
641  const volScalarField& rho = phase.thermo().rho();
642 
643  dilatationError =
644  fvc::ddt(alpha) + fvc::div(phase.alphaPhi())
645  - (fluid_.fvModels().source(alpha, rho) & rho)/rho;
646 
647  forAll(velGroup.sizeGroups(), i)
648  {
649  const sizeGroup& fi = velGroup.sizeGroups()[i];
650 
651  dilatationError -= Su_[fi.i()] - fvc::Sp(Sp_[fi.i()], fi);
652  }
653  }
654 }
655 
656 
657 void Foam::diameterModels::populationBalanceModel::calcAlphas()
658 {
659  alphas_() = Zero;
660 
661  forAllConstIter(HashTable<const velocityGroup*>, velocityGroupPtrs_, iter)
662  {
663  const phaseModel& phase = iter()->phase();
664 
665  alphas_() += max(phase, phase.residualAlpha());
666  }
667 }
668 
669 
671 Foam::diameterModels::populationBalanceModel::calcDsm()
672 {
673  tmp<volScalarField> tInvDsm
674  (
676  (
677  "invDsm",
678  mesh_,
680  )
681  );
682 
683  volScalarField& invDsm = tInvDsm.ref();
684 
685  forAllConstIter(HashTable<const velocityGroup*>, velocityGroupPtrs_, iter)
686  {
687  const phaseModel& phase = iter()->phase();
688 
689  invDsm += max(phase, phase.residualAlpha())/(phase.d()*alphas_());
690  }
691 
692  return 1/tInvDsm;
693 }
694 
695 
696 void Foam::diameterModels::populationBalanceModel::calcVelocity()
697 {
698  U_() = Zero;
699 
700  forAllConstIter(HashTable<const velocityGroup*>, velocityGroupPtrs_, iter)
701  {
702  const phaseModel& phase = iter()->phase();
703 
704  U_() += max(phase, phase.residualAlpha())*phase.U()/alphas_();
705  }
706 }
707 
708 bool Foam::diameterModels::populationBalanceModel::updateSources()
709 {
710  const bool result = sourceUpdateCounter_ % sourceUpdateInterval() == 0;
711 
712  ++ sourceUpdateCounter_;
713 
714  return result;
715 }
716 
717 
718 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
719 
721 (
722  const phaseSystem& fluid,
723  const word& name,
725 )
726 :
727  regIOobject
728  (
729  IOobject
730  (
731  name,
732  fluid.time().constant(),
733  fluid.mesh()
734  )
735  ),
736  fluid_(fluid),
737  pDmdt_(pDmdt),
738  mesh_(fluid_.mesh()),
739  name_(name),
740  dict_
741  (
742  fluid_.subDict("populationBalanceCoeffs").subDict(name_)
743  ),
744  pimple_(mesh_.lookupObject<pimpleControl>("solutionControl")),
745  continuousPhase_
746  (
747  mesh_.lookupObject<phaseModel>
748  (
749  IOobject::groupName("alpha", dict_.lookup("continuousPhase"))
750  )
751  ),
752  sizeGroups_(),
753  v_(),
754  delta_(),
755  Su_(),
756  Sp_(),
757  Sui_
758  (
759  IOobject
760  (
761  "Sui",
762  mesh_.time().timeName(),
763  mesh_
764  ),
765  mesh_,
766  dimensionedScalar(inv(dimTime), Zero)
767  ),
768  coalescenceModels_
769  (
770  dict_.lookup("coalescenceModels"),
771  coalescenceModel::iNew(*this)
772  ),
773  coalescenceRate_(),
774  coalescencePairs_(),
775  breakupModels_
776  (
777  dict_.lookup("breakupModels"),
778  breakupModel::iNew(*this)
779  ),
780  breakupRate_(),
781  binaryBreakupModels_
782  (
783  dict_.lookup("binaryBreakupModels"),
784  binaryBreakupModel::iNew(*this)
785  ),
786  binaryBreakupRate_(),
787  binaryBreakupPairs_(),
788  drift_
789  (
790  dict_.lookup("driftModels"),
791  driftModel::iNew(*this)
792  ),
793  driftRate_(),
794  nucleation_
795  (
796  dict_.lookup("nucleationModels"),
797  nucleationModel::iNew(*this)
798  ),
799  nucleationRate_(),
800  alphas_(),
801  dsm_(),
802  U_(),
803  sourceUpdateCounter_(0)
804 {
805  this->registerVelocityGroups();
806 
807  if (sizeGroups().size() < 3)
808  {
810  << "The populationBalance " << name_
811  << " requires a minimum number of three sizeGroups to be specified."
812  << exit(FatalError);
813  }
814 
815 
816  if (coalescenceModels_.size() != 0)
817  {
818  coalescenceRate_.set
819  (
820  new volScalarField
821  (
822  IOobject
823  (
824  "coalescenceRate",
825  mesh_.time().timeName(),
826  mesh_
827  ),
828  mesh_,
830  )
831  );
832 
833  forAll(sizeGroups(), i)
834  {
835  for (label j = 0; j <= i; j++)
836  {
837  coalescencePairs_.append(labelPair(i, j));
838  }
839  }
840  }
841 
842  if (breakupModels_.size() != 0)
843  {
844  breakupRate_.set
845  (
846  new volScalarField
847  (
848  IOobject
849  (
850  "breakupRate",
851  fluid_.time().timeName(),
852  mesh_
853  ),
854  mesh_,
855  dimensionedScalar(inv(dimTime), Zero)
856  )
857  );
858  }
859 
860  if (binaryBreakupModels_.size() != 0)
861  {
862  binaryBreakupRate_.set
863  (
864  new volScalarField
865  (
866  IOobject
867  (
868  "binaryBreakupRate",
869  fluid_.time().timeName(),
870  mesh_
871  ),
872  mesh_,
874  (
875  "binaryBreakupRate",
876  inv(dimVolume*dimTime),
877  Zero
878  )
879  )
880  );
881 
882  calcDeltas();
883 
884  forAll(sizeGroups(), i)
885  {
886  label j = 0;
887 
888  while (delta_[j][i].value() != 0)
889  {
890  binaryBreakupPairs_.append(labelPair(i, j));
891  j++;
892  }
893  }
894  }
895 
896  if (drift_.size() != 0)
897  {
898  driftRate_.set
899  (
900  new volScalarField
901  (
902  IOobject
903  (
904  "driftRate",
905  fluid_.time().timeName(),
906  mesh_
907  ),
908  mesh_,
910  )
911  );
912  }
913 
914  if (nucleation_.size() != 0)
915  {
916  nucleationRate_.set
917  (
918  new volScalarField
919  (
920  IOobject
921  (
922  "nucleationRate",
923  fluid_.time().timeName(),
924  mesh_
925  ),
926  mesh_,
928  (
929  "nucleationRate",
930  inv(dimTime*dimVolume),
931  Zero
932  )
933  )
934  );
935  }
936 
937  if (velocityGroupPtrs_.size() > 1)
938  {
939  alphas_.set
940  (
941  new volScalarField
942  (
943  IOobject
944  (
945  IOobject::groupName("alpha", this->name()),
946  fluid_.time().timeName(),
947  mesh_,
950  ),
951  mesh_,
953  )
954  );
955 
956  dsm_.set
957  (
958  new volScalarField
959  (
960  IOobject
961  (
962  IOobject::groupName("d", this->name()),
963  fluid_.time().timeName(),
964  mesh_,
967  ),
968  mesh_,
970  )
971  );
972 
973  U_.set
974  (
975  new volVectorField
976  (
977  IOobject
978  (
979  IOobject::groupName("U", this->name()),
980  fluid_.time().timeName(),
981  mesh_,
984  ),
985  mesh_,
987  )
988  );
989  }
990 }
991 
992 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
993 
995 {}
996 
997 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
998 
1001 {
1002  notImplemented("populationBalance::clone() const");
1003  return autoPtr<populationBalanceModel>(nullptr);
1004 }
1005 
1006 
1008 {
1009  return os.good();
1010 }
1011 
1012 
1015 (
1016  const label i,
1017  const dimensionedScalar& v
1018 ) const
1019 {
1020  const dimensionedScalar& x0 = sizeGroups()[0].x();
1021  const dimensionedScalar& xi = sizeGroups()[i].x();
1022  const dimensionedScalar& xm = sizeGroups().last().x();
1023  dimensionedScalar lowerBoundary(x0);
1024  dimensionedScalar upperBoundary(xm);
1025 
1026  if (i != 0) lowerBoundary = sizeGroups()[i-1].x();
1027 
1028  if (i != sizeGroups().size() - 1) upperBoundary = sizeGroups()[i+1].x();
1029 
1030  if ((i == 0 && v < x0) || (i == sizeGroups().size() - 1 && v > xm))
1031  {
1032  return v/xi;
1033  }
1034  else if (v < lowerBoundary || v > upperBoundary)
1035  {
1036  return 0;
1037  }
1038  else if (v.value() == xi.value())
1039  {
1040  return 1;
1041  }
1042  else if (v > xi)
1043  {
1044  return (upperBoundary - v)/(upperBoundary - xi);
1045  }
1046  else
1047  {
1048  return (v - lowerBoundary)/(xi - lowerBoundary);
1049  }
1050 }
1051 
1052 
1055 (
1056  const phaseModel& dispersedPhase
1057 ) const
1058 {
1059  return phaseInterface(dispersedPhase, continuousPhase_).sigma();
1060 }
1061 
1062 
1065 {
1066  return
1067  mesh_.lookupObject<phaseCompressible::momentumTransportModel>
1068  (
1070  (
1071  momentumTransportModel::typeName,
1072  continuousPhase_.name()
1073  )
1074  );
1075 }
1076 
1077 
1079 {
1080  if (!solveOnFinalIterOnly() || pimple_.finalPimpleIter())
1081  {
1082  const label nCorr = this->nCorr();
1083  const scalar tolerance =
1084  mesh_.solution().solverDict(name_).lookup<scalar>("tolerance");
1085 
1086  if (nCorr > 0)
1087  {
1088  precompute();
1089  }
1090 
1091  int iCorr = 0;
1092  scalar maxInitialResidual = 1;
1093 
1094  while (++iCorr <= nCorr && maxInitialResidual > tolerance)
1095  {
1096  Info<< "populationBalance "
1097  << this->name()
1098  << ": Iteration "
1099  << iCorr
1100  << endl;
1101 
1102  if (updateSources())
1103  {
1104  sources();
1105  }
1106 
1107  correctDilatationError();
1108 
1109  maxInitialResidual = 0;
1110 
1111  forAll(sizeGroups(), i)
1112  {
1113  sizeGroup& fi = sizeGroups_[i];
1114  const phaseModel& phase = fi.phase();
1115  const volScalarField& alpha = phase;
1116  const volScalarField& rho = phase.thermo().rho();
1117  const volScalarField& dilatationError =
1118  dilatationErrors_[phase.name()];
1119 
1120  fvScalarMatrix sizeGroupEqn
1121  (
1122  fvm::ddt(alpha, fi)
1123  + fvm::div(phase.alphaPhi(), fi)
1124  - fvm::Sp(dilatationError, fi)
1125  ==
1126  fvc::Su(Su_[i], fi)
1127  - fvm::Sp(Sp_[i], fi)
1128  + fluid_.fvModels().source(alpha, rho, fi)/rho
1129  - correction
1130  (
1131  fvm::Sp
1132  (
1133  max(phase.residualAlpha() - alpha, scalar(0))
1134  /this->mesh().time().deltaT(),
1135  fi
1136  )
1137  )
1138  );
1139 
1140  sizeGroupEqn.relax();
1141  fluid_.fvConstraints().constrain(sizeGroupEqn);
1142 
1143  maxInitialResidual = max
1144  (
1145  sizeGroupEqn.solve().initialResidual(),
1146  maxInitialResidual
1147  );
1148 
1149  fluid_.fvConstraints().constrain(fi);
1150  }
1151  }
1152 
1153  volScalarField fAlpha0
1154  (
1155  sizeGroups().first()*sizeGroups().first().phase()
1156  );
1157 
1158  volScalarField fAlphaN
1159  (
1160  sizeGroups().last()*sizeGroups().last().phase()
1161  );
1162 
1163  Info<< this->name() << " sizeGroup phase fraction first, last = "
1164  << fAlpha0.weightedAverage(this->mesh().V()).value()
1165  << ' ' << fAlphaN.weightedAverage(this->mesh().V()).value()
1166  << endl;
1167  }
1168 }
1169 
1170 
1172 {
1173  if (velocityGroupPtrs_.size() > 1)
1174  {
1175  calcAlphas();
1176  dsm_() = calcDsm();
1177  calcVelocity();
1178  }
1179 }
1180 
1181 
1182 // ************************************************************************* //
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:38
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Templated abstract base class for multiphase compressible turbulence models.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
label phasei
Definition: pEqn.H:27
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
populationBalanceModel(const phaseSystem &fluid, const word &name, HashPtrTable< volScalarField, phaseInterfaceKey, phaseInterfaceKey::hash > &pDmdt)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const phaseCompressible::momentumTransportModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
const dimensionSet dimless
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:59
fvMesh & mesh
dimensionedScalar neg(const dimensionedScalar &ds)
const dimensionedScalar eta(const label i, const dimensionedScalar &v) const
Return allocation coefficient.
const dimensionSet dimLength
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
autoPtr< populationBalanceModel > clone() const
Return clone.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
Calculate the first temporal derivative.
dimensionedScalar pos(const dimensionedScalar &ds)
stressControl lookup("compactNormalStress") >> compactNormalStress
static word groupName(Name name, const word &group)
virtual void precompute()
Precompute diameter independent expressions.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Calculate the matrix for the first temporal derivative.
const Type & value() const
Return const reference to value.
void correct()
Correct derived quantities.
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
Calculate the field for explicit evaluation of implicit and explicit sources.
Calculate the divergence of the given field.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
const dimensionSet dimVelocity
defineTypeNameAndDebug(combustionModel, 0)
HashPtrTable< volScalarField, phaseInterfaceKey, phaseInterfaceKey::hash > dmdtfTable
Definition: phaseSystem.H:93
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.
Calculate the matrix for the divergence of the given field and flux.
#define notImplemented(functionName)
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:343
phaseCompressibleMomentumTransportModel momentumTransportModel
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
int nCorr
Definition: createControls.H:3
const dimensionSet dimVolume
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const tmp< volScalarField > sigmaWithContinuousPhase(const phaseModel &dispersedPhase) const
Return the surface tension coefficient between a given dispersed.
const tmp< volScalarField::Internal > & Su
Definition: alphaSuSp.H:6
A class for managing temporary objects.
Definition: PtrList.H:53
void solve()
Solve the population balance equation.
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
Calculate the matrix for implicit and explicit sources.
Namespace for OpenFOAM.