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-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 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 exemplary 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"
210 #include "HashPtrTable.H"
211 #include "Pair.H"
212 
213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214 
215 namespace Foam
216 {
217 namespace diameterModels
218 {
219 
220 class coalescenceModel;
221 class breakupModel;
222 class binaryBreakupModel;
223 class driftModel;
224 class nucleationModel;
225 
226 /*---------------------------------------------------------------------------*\
227  Class populationBalanceModel Declaration
228 \*---------------------------------------------------------------------------*/
229 
231 :
232  public regIOobject
233 {
234  // Private Data
235 
236  //- Reference to the phaseSystem
237  const phaseSystem& fluid_;
238 
239  //- Interfacial mass transfer rates
240  phaseSystem::dmdtfTable dmdtfs_;
241 
242  //- Reference to the mesh
243  const fvMesh& mesh_;
244 
245  //- Name of the populationBalance
246  word name_;
247 
248  //- Dictionary
249  dictionary dict_;
250 
251  //- Continuous phase
252  const phaseModel& continuousPhase_;
253 
254  //- Velocity groups belonging to this populationBalance
255  HashTable<const velocityGroup*> velocityGroupPtrs_;
256 
257  //- Size groups belonging to this populationBalance
258  UPtrList<sizeGroup> sizeGroups_;
259 
260  //- sizeGroup boundaries
262 
263  //- Section width required for binary breakup formulation
265 
266  //- Explicitly treated sources
268 
269  //- Implicitly treated sources
271 
272  //- Dilatation errors
273  HashTable<volScalarField> dilatationErrors_;
274 
275  //- Field for caching sources
276  volScalarField Sui_;
277 
278  //- Coalescence models
279  PtrList<coalescenceModel> coalescenceModels_;
280 
281  //- Coalescence rate
282  autoPtr<volScalarField> coalescenceRate_;
283 
284  //- Coalescence relevant size group pairs
285  List<labelPair> coalescencePairs_;
286 
287  //- Breakup models
288  PtrList<breakupModel> breakupModels_;
289 
290  //- Breakup rate
291  autoPtr<volScalarField> breakupRate_;
292 
293  //- Binary breakup models
294  PtrList<binaryBreakupModel> binaryBreakupModels_;
295 
296  //- Binary breakup rate
297  autoPtr<volScalarField> binaryBreakupRate_;
298 
299  //- Binary breakup relevant size group pairs
300  List<labelPair> binaryBreakupPairs_;
301 
302  //- Drift models
303  PtrList<driftModel> drift_;
304 
305  //- Drift rate
306  autoPtr<volScalarField> driftRate_;
307 
308  //- Zeroeth order models
309  PtrList<nucleationModel> nucleation_;
310 
311  //- Zeroeth order rate
312  autoPtr<volScalarField> nucleationRate_;
313 
314  //- Total void fraction
315  autoPtr<volScalarField> alphas_;
316 
317  //- Mean Sauter diameter
319 
320  //- Average velocity
322 
323  //- Counter for interval between source term updates
324  label sourceUpdateCounter_;
325 
326 
327  // Private Member Functions
328 
329  void registerVelocityGroups();
330 
331  void registerSizeGroups(sizeGroup& group);
332 
333  void initialiseDmdtfs();
334 
335  void createPhasePairs();
336 
337  void precompute();
338 
339  void birthByCoalescence(const label j, const label k);
340 
341  void deathByCoalescence(const label i, const label j);
342 
343  void birthByBreakup(const label k, const label model);
344 
345  void deathByBreakup(const label i);
346 
347  void calcDeltas();
348 
349  void birthByBinaryBreakup(const label i, const label j);
350 
351  void deathByBinaryBreakup(const label j, const label i);
352 
353  void drift(const label i, driftModel& model);
354 
355  void nucleation(const label i, nucleationModel& model);
356 
357  void sources();
358 
359  void correctDilatationError();
360 
361  void calcAlphas();
362 
363  tmp<volScalarField> calcDsm();
364 
365  void calcVelocity();
366 
367  //- Return whether the sources should be updated on this iteration
368  bool updateSources();
369 
370  //- Return the interval at which the sources are updated
371  inline label sourceUpdateInterval() const;
372 
373 public:
374 
375  //- Runtime type information
376  TypeName("populationBalanceModel");
377 
378 
379  // Constructors
380 
381  //- Construct for a fluid
383 
384  //- Return clone
386 
387  //- Return a pointer to a new populationBalanceModel object created on
388  // freestore from Istream
389  class iNew
390  {
391  const phaseSystem& fluid_;
392 
393  public:
394 
395  iNew(const phaseSystem& fluid)
396  :
397  fluid_(fluid)
398  {}
399 
401  {
402  const word name(is);
403 
404  Info << "Setting up population balance: " << name << endl;
405 
407  (
408  new populationBalanceModel(fluid_, name)
409  );
410  }
411  };
412 
413 
414  //- Destructor
415  virtual ~populationBalanceModel();
416 
417 
418  // Member Functions
419 
420  //- Dummy write for regIOobject
421  bool writeData(Ostream&) const;
422 
423  //- Return reference to the phaseSystem
424  inline const phaseSystem& fluid() const;
425 
426  //- Return reference to the interfacial mass transfer rates
427  inline const phaseSystem::dmdtfTable& dmdtfs() const;
428 
429  //- Return reference to the mesh
430  inline const fvMesh& mesh() const;
431 
432  //- Return populationBalanceCoeffs dictionary
433  inline const dictionary& dict() const;
434 
435  //- Return the number of corrections
436  inline label nCorr() const;
437 
438  //- Solve on final pimple iteration only
439  inline Switch solveOnFinalIterOnly() const;
440 
441  //- Return continuous phase
442  inline const phaseModel& continuousPhase() const;
443 
444  //- Return the velocity groups belonging to this populationBalance
446 
447  //- Return the size groups belonging to this populationBalance
448  inline const UPtrList<sizeGroup>& sizeGroups() const;
449 
450  //- Return implicit source terms
451  inline const volScalarField& Sp(const label i) const;
452 
453  //- Return dilatation obtained from source terms
454  inline const HashTable<volScalarField>& sourceDilatation() const;
455 
456  //- Return coalescence relevant size group pairs
457  inline const List<labelPair>& coalescencePairs() const;
458 
459  //- Return binary breakup relevant size group pairs
460  inline const List<labelPair>& binaryBreakupPairs() const;
461 
462  //- Return total void of phases belonging to this populationBalance
463  inline const volScalarField& alphas() const;
464 
465  //- Return average velocity
466  inline const volVectorField& U() const;
467 
468  //- Return allocation coefficient
469  const dimensionedScalar eta
470  (
471  const label i,
472  const dimensionedScalar& v
473  ) const;
474 
475  //- Return the surface tension coefficient between a given dispersed
476  // and the continuous phase
478  (
479  const phaseModel& dispersedPhase
480  ) const;
481 
482  //- Return reference to momentumTransport model of the continuous phase
484  continuousTurbulence() const;
485 
486  //- Solve the population balance equation
487  void solve();
488 
489  //- Correct derived quantities
490  void correct();
491 };
492 
493 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
494 
495 } // End namespace diameterModels
496 } // End namespace Foam
497 
498 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
499 
500 #include "populationBalanceModelI.H"
501 
502 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
503 
504 #endif
505 
506 // ************************************************************************* //
label k
Generic GeometricField class.
An STL-conforming hash table.
Definition: HashTable.H:127
word group() const
Return group (extension part of name)
Definition: IOobject.C:324
const word & name() const
Return name.
Definition: IOobject.H:310
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
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
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.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:66
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Base class for drift models.
Definition: driftModel.H:54
Base class for nucleation models.
Return a pointer to a new populationBalanceModel object created on.
autoPtr< populationBalanceModel > operator()(Istream &is) const
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
const HashTable< volScalarField > & sourceDilatation() const
Return dilatation obtained from source terms.
const phaseSystem & fluid() const
Return reference to the phaseSystem.
const dimensionedScalar eta(const label i, const dimensionedScalar &v) const
Return allocation coefficient.
bool writeData(Ostream &) const
Dummy write for regIOobject.
const List< labelPair > & coalescencePairs() const
Return coalescence relevant size group pairs.
autoPtr< populationBalanceModel > clone() const
Return clone.
const List< labelPair > & binaryBreakupPairs() const
Return binary breakup relevant size group pairs.
const volVectorField & U() const
Return average velocity.
populationBalanceModel(const phaseSystem &fluid, const word &name)
Construct for a fluid.
const UPtrList< sizeGroup > & sizeGroups() const
Return the size groups belonging to this populationBalance.
Switch solveOnFinalIterOnly() const
Solve on final pimple iteration only.
const phaseModel & continuousPhase() const
Return continuous phase.
TypeName("populationBalanceModel")
Runtime type information.
const HashTable< const velocityGroup * > & velocityGroupPtrs() const
Return the velocity groups belonging to this populationBalance.
const volScalarField & alphas() const
Return total void of phases belonging to this populationBalance.
const dictionary & dict() const
Return populationBalanceCoeffs dictionary.
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.
label nCorr() const
Return the number of corrections.
const phaseCompressible::momentumTransportModel & continuousTurbulence() const
Return reference to momentumTransport model of the continuous phase.
const volScalarField & Sp(const label i) const
Return implicit source terms.
const phaseSystem::dmdtfTable & dmdtfs() const
Return reference to the interfacial mass transfer rates.
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
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Templated abstract base class for multiphase compressible turbulence models.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
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
Namespace for OpenFOAM.
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:251
messageStream Info