populationBalanceModel.H
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 Class
25  Foam::diameterModels::populationBalanceModel
26 
27 Description
28  Model for tracking the evolution of a dispersed phase size distribution due
29  to coalescence (synonymous with coagulation, aggregation, agglomeration)
30  and breakup events as well as density or phase changes. Provides an
31  approximate solution of the population balance equation by means of a class
32  method. The underlying theory is described in the article of Lehnigk et al.
33  (2021).
34 
35  The size distribution, expressed through a volume-based number density
36  function, is discretised using the fixot pivot technique of Kumar and
37  Ramkrishna (1996). Thereby, the population balance equation is transformed
38  into a series of transport equations for the particle (bubble, droplet)
39  number concentrations in separate size classes that are coupled through
40  their source terms. The discretisation is based on representative particle
41  volumes, which are provided by the user through the corresponding sphere
42  equivalent diameters.
43 
44  Since the representative volumes are fixed a priori and the total dispersed
45  phase volume already available from solving the phase continuity equation,
46  the model only determines the evolution of the individual size class
47  fractions
48 
49  \f[
50  f_{i,\varphi} = \frac{\alpha_{i,\varphi}}{\alpha_{\varphi}}\,,
51  \f]
52 
53  where \f$\alpha_{i,\varphi}\f$ is the volume fraction of the size class and
54  \f$\alpha_{\varphi}\f$ the total phase fraction of phase \f$\varphi\f$.
55 
56  The source terms are formulated such that the first and second moment of
57  the distribution, i.e. the total particle number and volume, are conserved
58  irrespective of the discretisation of the size domain. The treatment of
59  particle breakup depends on the selected breakup submodels. For models
60  which provide a total breakup frequency and a separate daughter size
61  distribution function, the formulation provided Kumar and Ramkrishna (1996)
62  is utilised, which is applicable both for binary and multiple breakup
63  events. Currently, only field-independent daughter size distribution models
64  are allowed. In case of binary breakup models that provide the breakup
65  frequency between a size class pair directly, the formulation of Liao et
66  al. (2018) is adopted, which is computationally more efficient compared to
67  first extracting the field-dependent daughter size distribution and then
68  consuming it in the formulation of Kumar and Ramkrishna. The source terms
69  describing a drift of the size distribution through particle growth or
70  shrinkage are derived using upwind differencing, thus ensuring conservation
71  of the total particle number and volume. Note that due to the volume-based
72  implementation, both density as well as phase change lead to a drift of the
73  size distribution function. Further, users can specify multiple submodels
74  for each mechanism, whose contributions are added up.
75 
76  The model also allows to distribute the size classes over multiple
77  representative phases with identical physical properties that collectively
78  define the dispersed phase. Thereby, size class fields can be transported
79  with different velocity fields in order to account for the size dependency
80  of the particle motion. A possible mass transfer between representative
81  phases by means of coalescence, breakup and drift is taken into account.
82  Similarly, the spatial evolution of secondary particle properties such as
83  the particle surface area can be tracked.
84 
85  The key variable during a simulation is the Sauter diameter, which is
86  computed from the size class fractions of the corresponding phase. The
87  underlying size distribution can be extracted from the simulation using the
88  functionObject 'sizeDistribution'. Integral and mean properties of a size
89  distribution can be computed with the functionObject 'moments'.
90 
91  Verification cases for the population balance modeling functionality are
92  provided in test/multiphase/multiphaseEulerFoam/populationBalance.
93 
94  References:
95  \verbatim
96  Lehnigk, R., Bainbridge, W., Liao, Y., Lucas, D., Niemi, T.,
97  Peltola, J., & Schlegel, F. (2021).
98  An openā€source population balance modeling framework for the simulation
99  of polydisperse multiphase flows.
100  AIChE Journal, 68(3), e17539.
101  \endverbatim
102 
103  \verbatim
104  Coalescence and breakup term formulation:
105  Kumar, S., & Ramkrishna, D. (1996).
106  On the solution of population balance equations by discretization-I. A
107  fixed pivot technique.
108  Chemical Engineering Science, 51(8), 1311-1332.
109  \endverbatim
110 
111  \verbatim
112  Binary breakup term formulation:
113  Liao, Y., Oertel, R., Kriebitzsch, S., Schlegel, F., & Lucas, D. (2018).
114  A discrete population balance equation for binary breakage.
115  International Journal for Numerical Methods in Fluids, 87(4), 202-215.
116  \endverbatim
117 
118 Usage
119  Excerpt from an examplary phaseProperties dictionary:
120  \verbatim
121  type populationBalanceMultiphaseSystem;
122 
123  ...
124 
125  populationBalances (bubbles);
126 
127  air
128  {
129  type pureIsothermalPhaseModel;
130 
131  diameterModel velocityGroup;
132 
133  velocityGroupCoeffs
134  {
135  populationBalance bubbles;
136 
137  shapeModel spherical;
138 
139  sizeGroups
140  (
141  f1 {dSph 1e-3; value 1.0;}
142  f2 {dSph 2e-3; value 0.0;}
143  f3 {dSph 3e-3; value 0.0;}
144  f4 {dSph 4e-3; value 0.0;}
145  f5 {dSph 5e-3; value 0.0;}
146  ...
147  );
148  }
149 
150  residualAlpha 1e-6;
151  }
152 
153  ...
154 
155  populationBalanceCoeffs
156  {
157  bubbles
158  {
159  continuousPhase water;
160 
161  coalescenceModels
162  (
163  LehrMilliesMewes{}
164  );
165 
166  binaryBreakupModels
167  (
168  LehrMilliesMewes{}
169  );
170 
171  breakupModels
172  ();
173 
174  driftModels
175  (
176  densityChange{}
177  );
178 
179  nucleationModels
180  ();
181  }
182  }
183  \endverbatim
184 
185 See also
186  Foam::PopulationBalancePhaseSystem
187  Foam::diameterModels::sizeGroup
188  Foam::diameterModels::velocityGroup
189  Foam::diameterModels::SecondaryPropertyModel
190  Foam::diameterModels::coalescenceModel
191  Foam::diameterModels::breakupModel
192  Foam::diameterModels::daughterSizeDistributionModel
193  Foam::diameterModels::binaryBreakupModel
194  Foam::diameterModels::driftModel
195  Foam::diameterModels::nucleationModel
196  Foam::functionObjects::sizeDistribution
197  Foam::functionObjects::moments
198 
199 SourceFiles
200  populationBalanceModel.C
201 
202 \*---------------------------------------------------------------------------*/
203 
204 #ifndef populationBalanceModel_H
205 #define populationBalanceModel_H
206 
207 #include "sizeGroup.H"
208 #include "phaseSystem.H"
209 #include "pimpleControl.H"
211 #include "HashPtrTable.H"
212 #include "Pair.H"
213 
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 
216 namespace Foam
217 {
218 namespace diameterModels
219 {
220 
221 class coalescenceModel;
222 class breakupModel;
223 class binaryBreakupModel;
224 class driftModel;
225 class nucleationModel;
226 
227 /*---------------------------------------------------------------------------*\
228  Class populationBalanceModel Declaration
229 \*---------------------------------------------------------------------------*/
232 :
233  public regIOobject
234 {
235 private:
236 
237  // Private Data
238 
239  //- Reference to the phaseSystem
240  const phaseSystem& fluid_;
241 
242  //- Interfacial mass transfer rates
243  phaseSystem::dmdtfTable& pDmdt_;
244 
245  //- Reference to the mesh
246  const fvMesh& mesh_;
247 
248  //- Name of the populationBalance
249  word name_;
250 
251  //- Dictionary
252  dictionary dict_;
253 
254  //- Reference to pimpleControl
255  const pimpleControl& pimple_;
256 
257  //- Continuous phase
258  const phaseModel& continuousPhase_;
259 
260  //- Velocity groups belonging to this populationBalance
261  HashTable<const velocityGroup*> velocityGroupPtrs_;
262 
263  //- Size groups belonging to this populationBalance
264  UPtrList<sizeGroup> sizeGroups_;
265 
266  //- sizeGroup boundaries
268 
269  //- Section width required for binary breakup formulation
271 
272  //- Explicitly treated sources
274 
275  //- Implicitly treated sources
277 
278  //- Dilatation errors
279  HashTable<volScalarField> dilatationErrors_;
280 
281  //- Field for caching sources
282  volScalarField Sui_;
283 
284  //- Coalescence models
285  PtrList<coalescenceModel> coalescenceModels_;
286 
287  //- Coalescence rate
288  autoPtr<volScalarField> coalescenceRate_;
289 
290  //- Coalescence relevant size group pairs
291  List<labelPair> coalescencePairs_;
292 
293  //- Breakup models
294  PtrList<breakupModel> breakupModels_;
295 
296  //- Breakup rate
297  autoPtr<volScalarField> breakupRate_;
298 
299  //- Binary breakup models
300  PtrList<binaryBreakupModel> binaryBreakupModels_;
301 
302  //- Binary breakup rate
303  autoPtr<volScalarField> binaryBreakupRate_;
304 
305  //- Binary breakup relevant size group pairs
306  List<labelPair> binaryBreakupPairs_;
307 
308  //- Drift models
309  PtrList<driftModel> drift_;
310 
311  //- Drift rate
312  autoPtr<volScalarField> driftRate_;
313 
314  //- Zeroeth order models
315  PtrList<nucleationModel> nucleation_;
316 
317  //- Zeroeth order rate
318  autoPtr<volScalarField> nucleationRate_;
319 
320  //- Total void fraction
321  autoPtr<volScalarField> alphas_;
322 
323  //- Mean Sauter diameter
325 
326  //- Average velocity
328 
329  //- Counter for interval between source term updates
330  label sourceUpdateCounter_;
331 
332 
333  // Private Member Functions
334 
335  void registerVelocityGroups();
336 
337  void registerSizeGroups(sizeGroup& group);
338 
339  void createPhasePairs();
340 
341  void precompute();
342 
343  void birthByCoalescence(const label j, const label k);
344 
345  void deathByCoalescence(const label i, const label j);
346 
347  void birthByBreakup(const label k, const label model);
348 
349  void deathByBreakup(const label i);
350 
351  void calcDeltas();
352 
353  void birthByBinaryBreakup(const label i, const label j);
354 
355  void deathByBinaryBreakup(const label j, const label i);
356 
357  void drift(const label i, driftModel& model);
358 
359  void nucleation(const label i, nucleationModel& model);
360 
361  void sources();
362 
363  void correctDilatationError();
364 
365  void calcAlphas();
366 
367  tmp<volScalarField> calcDsm();
368 
369  void calcVelocity();
370 
371  //- Return whether the sources should be updated on this iteration
372  bool updateSources();
373 
374  //- Return the interval at which the sources are updated
375  inline label sourceUpdateInterval() const;
376 
377 public:
378 
379  //- Runtime type information
380  TypeName("populationBalanceModel");
381 
382 
383  // Constructor
384 
386  (
387  const phaseSystem& fluid,
388  const word& name,
390  <
394  >& pDmdt
395  );
396 
397  //- Return clone
399 
400  //- Return a pointer to a new populationBalanceModel object created on
401  // freestore from Istream
402  class iNew
403  {
404  const phaseSystem& fluid_;
405 
406  phaseSystem::dmdtfTable& pDmdt_;
407 
408  public:
409 
411  (
412  const phaseSystem& fluid,
414  )
415  :
416  fluid_(fluid),
417  pDmdt_(pDmdt)
418  {}
421  {
422  word name(is);
423 
424  Info << "Setting up population balance: " << name << endl;
425 
427  (
428  new populationBalanceModel(fluid_, name, pDmdt_)
429  );
430  }
431  };
432 
433 
434  //- Destructor
435  virtual ~populationBalanceModel();
436 
437  // Member Functions
438 
439  //- Dummy write for regIOobject
440  bool writeData(Ostream&) const;
441 
442  //- Return reference to the phaseSystem
443  inline const phaseSystem& fluid() const;
444 
445  //- Return reference to the mesh
446  inline const fvMesh& mesh() const;
447 
448  //- Return populationBalanceCoeffs dictionary
449  inline const dictionary& dict() const;
450 
451  //- Return the number of corrections
452  inline label nCorr() const;
453 
454  //- Solve on final pimple iteration only
455  inline Switch solveOnFinalIterOnly() const;
456 
457  //- Return continuous phase
458  inline const phaseModel& continuousPhase() const;
459 
460  //- Return the velocity groups belonging to this populationBalance
462 
463  //- Return the size groups belonging to this populationBalance
464  inline const UPtrList<sizeGroup>& sizeGroups() const;
465 
466  //- Return implicit source terms
467  inline const volScalarField& Sp(const label i) const;
468 
469  //- Return dilatation obtained from source terms
470  inline const HashTable<volScalarField>& sourceDilatation() const;
471 
472  //- Return coalescence relevant size group pairs
473  inline const List<labelPair>& coalescencePairs() const;
474 
475  //- Return binary breakup relevant size group pairs
476  inline const List<labelPair>& binaryBreakupPairs() const;
477 
478  //- Return total void of phases belonging to this populationBalance
479  inline const volScalarField& alphas() const;
480 
481  //- Return average velocity
482  inline const volVectorField& U() const;
483 
484  //- Return allocation coefficient
485  const dimensionedScalar eta
486  (
487  const label i,
488  const dimensionedScalar& v
489  ) const;
490 
491  //- Return the surface tension coefficient between a given dispersed
492  // and the continuous phase
494  (
495  const phaseModel& dispersedPhase
496  ) const;
497 
498  //- Return reference to turbulence model of the continuous phase
500  continuousTurbulence() const;
501 
502  //- Solve the population balance equation
503  void solve();
504 
505  //- Correct derived quantities
506  void correct();
507 };
508 
509 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
510 
511 } // End namespace diameterModels
512 } // End namespace Foam
513 
514 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
515 
516 #include "populationBalanceModelI.H"
517 
518 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
519 
520 #endif
521 
522 // ************************************************************************* //
Templated abstract base class for multiphase compressible turbulence models.
Base class for drift models.
Definition: driftModel.H:53
const word & name() const
Return name.
Definition: IOobject.H:315
const volScalarField & Sp(const label i) const
Return implicit source terms.
const volVectorField & U() const
Return average velocity.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
populationBalanceModel(const phaseSystem &fluid, const word &name, HashPtrTable< volScalarField, phaseInterfaceKey, phaseInterfaceKey::hash > &pDmdt)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
const List< labelPair > & coalescencePairs() const
Return coalescence relevant size group pairs.
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
Word-pair based class used for keying interface models in hash tables.
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
const dictionary & dict() const
Return populationBalanceCoeffs dictionary.
autoPtr< populationBalanceModel > operator()(Istream &is) const
word group() const
Return group (extension part of name)
Definition: IOobject.C:324
label k
Boltzmann constant.
bool writeData(Ostream &) const
Dummy write for regIOobject.
const dimensionedScalar eta(const label i, const dimensionedScalar &v) const
Return allocation coefficient.
autoPtr< populationBalanceModel > clone() const
Return clone.
Base class for nucleation models.
label nCorr() const
Return the number of corrections.
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:99
const phaseModel & continuousPhase() const
Return continuous phase.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:68
const phaseSystem & fluid() const
Return reference to the phaseSystem.
A class for handling words, derived from string.
Definition: word.H:59
void correct()
Correct derived quantities.
const HashTable< volScalarField > & sourceDilatation() const
Return dilatation obtained from source terms.
An STL-conforming hash table.
Definition: HashTable.H:61
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:54
iNew(const phaseSystem &fluid, phaseSystem::dmdtfTable &pDmdt)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
Return a pointer to a new populationBalanceModel object created on.
const UPtrList< sizeGroup > & sizeGroups() const
Return the size groups belonging to this populationBalance.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Pimple control class. Provides time-loop control methods which exit the simulation once convergence c...
Definition: pimpleControl.H:74
const volScalarField & alphas() const
Return total void of phases belonging to this populationBalance.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:52
messageStream Info
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.
const List< labelPair > & binaryBreakupPairs() const
Return binary breakup relevant size group pairs.
Switch solveOnFinalIterOnly() const
Solve on final pimple iteration only.
A class for managing temporary objects.
Definition: PtrList.H:53
void solve()
Solve the population balance equation.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
TypeName("populationBalanceModel")
Runtime type information.
Namespace for OpenFOAM.
const HashTable< const velocityGroup * > & velocityGroupPtrs() const
Return the velocity groups belonging to this populationBalance.