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-2021 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  Class that solves the univariate population balance equation by means of
29  a class method (also called sectional or discrete method). The internal
30  coordinate is set to the particle volume, so the equation is based on
31  a transport equation of the volume-based number density function. The
32  discretisation is done using the fixed pivot technique of Kumar and
33  Ramkrishna (1996). The source terms are written in a way that particle
34  number and mass are preserved. Coalescence (coagulation), breakup, drift
35  (growth and surface loss) as well as nucleation are supported.
36  For the discrete breakup term two recipes are available, depending on the
37  model choice. For models which state a total breakup rate and a separate
38  daughter size distribution function, the formulation of Kumar and Ramkrishna
39  (1996) is applied which is applicable for binary and multiple breakup
40  events. The second formulation is given by Liao et al. (2018). It is useful
41  for binary breakup models which give the breakup rate between a sizeGroup
42  pair directly, without an explicit expression for the daughter size
43  distribution. The drift term is implemented using a finite difference upwind
44  scheme. Although it is diffusive, it ensures a stable and
45  number-conservative solution.
46 
47  The implementation allows to split the population balance over multiple
48  velocity fields using the capability of multiphaseEulerFoam to solve
49  for n momentum equations. It is also possible to define multiple population
50  balances, e.g. bubbles and droplets simultaneously.
51 
52  References:
53  \verbatim
54  Coalescence and breakup term formulation:
55  Kumar, S., & Ramkrishna, D. (1996).
56  On the solution of population balance equations by discretization-I. A
57  fixed pivot technique.
58  Chemical Engineering Science, 51(8), 1311-1332.
59  \endverbatim
60 
61  \verbatim
62  Binary breakup term formulation:
63  Liao, Y., Oertel, R., Kriebitzsch, S., Schlegel, F., & Lucas, D. (2018).
64  A discrete population balance equation for binary breakage.
65  International Journal for Numerical Methods in Fluids, 87(4), 202-215.
66  \endverbatim
67 
68 Usage
69  Example excerpt from a phaseProperties dictionary.
70  \verbatim
71  type populationBalanceTwoPhaseSystem;
72 
73  phases (air water);
74 
75  populationBalances (bubbles);
76 
77  air
78  {
79  type purePhaseModel;
80  diameterModel velocityGroup;
81  velocityGroupCoeffs
82  {
83  populationBalance bubbles;
84 
85  shapeModel constant;
86 
87  sizeGroups
88  (
89  f0{dSph 1.00e-3; value 0;}
90  f1{dSph 1.08e-3; value 0;}
91  f2{dSph 1.16e-3; value 0.25;}
92  f3{dSph 1.25e-3; value 0.5;}
93  f4{dSph 1.36e-3; value 0.25;}
94  f5{dSph 1.46e-3; value 0;}
95  ...
96  );
97  }
98 
99  residualAlpha 1e-6;
100  }
101 
102  populationBalanceCoeffs
103  {
104  bubbles
105  {
106  continuousPhase water;
107 
108  coalescenceModels
109  (
110  hydrodynamic
111  {
112  C 0.25;
113  }
114  );
115 
116  binaryBreakupModels
117  ();
118 
119  breakupModels
120  (
121  exponential
122  {
123  C 0.5;
124  exponent 0.01;
125  daughterSizeDistributionModel uniform;
126  }
127  );
128 
129  driftModels
130  (
131  densityChange{}
132  );
133 
134  nucleationModels
135  ();
136  }
137  }
138  \endverbatim
139 
140 See also
141  Foam::diameterModels::sizeGroup
142  Foam::diameterModels::velocityGroup
143 
144 SourceFiles
145  populationBalanceModel.C
146 
147 \*---------------------------------------------------------------------------*/
148 
149 #ifndef populationBalanceModel_H
150 #define populationBalanceModel_H
151 
152 #include "sizeGroup.H"
153 #include "phasePair.H"
154 #include "pimpleControl.H"
156 #include "HashPtrTable.H"
157 #include "Pair.H"
158 
159 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160 
161 namespace Foam
162 {
163 
164 class phaseSystem;
165 
166 namespace diameterModels
167 {
168 
169 class coalescenceModel;
170 class breakupModel;
171 class binaryBreakupModel;
172 class driftModel;
173 class nucleationModel;
174 
175 /*---------------------------------------------------------------------------*\
176  Class populationBalanceModel Declaration
177 \*---------------------------------------------------------------------------*/
180 :
181  public regIOobject
182 {
183 private:
184 
185  // Private Typedefs
186 
187  typedef
190 
191  typedef
193  pDmdtTable;
194 
195 
196  // Private Data
197 
198  //- Reference to the phaseSystem
199  const phaseSystem& fluid_;
200 
201  //- Interfacial mass transfer rates
202  pDmdtTable& pDmdt_;
203 
204  //- Reference to the mesh
205  const fvMesh& mesh_;
206 
207  //- Name of the populationBalance
208  word name_;
209 
210  //- Dictionary
211  dictionary dict_;
212 
213  //- Reference to pimpleControl
214  const pimpleControl& pimple_;
215 
216  //- Continuous phase
217  const phaseModel& continuousPhase_;
218 
219  //- velocityGroups belonging to this populationBalance
220  UPtrList<velocityGroup> velocityGroups_;
221 
222  //- sizeGroups belonging to this populationBalance
223  UPtrList<sizeGroup> sizeGroups_;
224 
225  //- List of unordered phasePairs in this populationBalance
226  phasePairTable phasePairs_;
227 
228  //- sizeGroup boundaries
230 
231  //- Section width required for binary breakup formulation
233 
234  //- Explicitly treated sources
236 
237  //- Implicitly treated sources
239 
240  //- Field for caching sources
241  volScalarField Sui_;
242 
243  //- Coalescence models
244  PtrList<coalescenceModel> coalescenceModels_;
245 
246  //- Coalescence rate
247  autoPtr<volScalarField> coalescenceRate_;
248 
249  //- Coalescence relevant size group pairs
250  List<labelPair> coalescencePairs_;
251 
252  //- Breakup models
253  PtrList<breakupModel> breakupModels_;
254 
255  //- Breakup rate
256  autoPtr<volScalarField> breakupRate_;
257 
258  //- Binary breakup models
259  PtrList<binaryBreakupModel> binaryBreakupModels_;
260 
261  //- Binary breakup rate
262  autoPtr<volScalarField> binaryBreakupRate_;
263 
264  //- Binary breakup relevant size group pairs
265  List<labelPair> binaryBreakupPairs_;
266 
267  //- Drift models
268  PtrList<driftModel> drift_;
269 
270  //- Drift rate
271  autoPtr<volScalarField> driftRate_;
272 
273  //- Zeroeth order models
274  PtrList<nucleationModel> nucleation_;
275 
276  //- Zeroeth order rate
277  autoPtr<volScalarField> nucleationRate_;
278 
279  //- Total void fraction
280  autoPtr<volScalarField> alphas_;
281 
282  //- Mean Sauter diameter
284 
285  //- Average velocity
287 
288  //- Counter for interval between source term updates
289  label sourceUpdateCounter_;
290 
291 
292  // Private Member Functions
293 
294  void registerVelocityGroups();
295 
296  void registerSizeGroups(sizeGroup& group);
297 
298  void createPhasePairs();
299 
300  void precompute();
301 
302  void birthByCoalescence(const label j, const label k);
303 
304  void deathByCoalescence(const label i, const label j);
305 
306  void birthByBreakup(const label k, const label model);
307 
308  void deathByBreakup(const label i);
309 
310  void calcDeltas();
311 
312  void birthByBinaryBreakup(const label i, const label j);
313 
314  void deathByBinaryBreakup(const label j, const label i);
315 
316  void drift(const label i, driftModel& model);
317 
318  void nucleation(const label i, nucleationModel& model);
319 
320  void sources();
321 
322  void dmdt();
323 
324  void calcAlphas();
325 
326  tmp<volScalarField> calcDsm();
327 
328  void calcVelocity();
329 
330  //- Return whether the sources should be updated on this iteration
331  bool updateSources();
332 
333  //- Return the number of corrections
334  inline label nCorr() const;
335 
336  //- Return the interval at which the sources are updated
337  inline label sourceUpdateInterval() const;
338 
339 public:
340 
341  //- Runtime type information
342  TypeName("populationBalanceModel");
343 
344 
345  // Constructor
346 
348  (
349  const phaseSystem& fluid,
350  const word& name,
352  <
354  phasePairKey,
356  >& pDmdt
357  );
358 
359  //- Return clone
361 
362  //- Return a pointer to a new populationBalanceModel object created on
363  // freestore from Istream
364  class iNew
365  {
366  const phaseSystem& fluid_;
367 
369  pDmdt_;
370 
371  public:
372 
374  (
375  const phaseSystem& fluid,
377  pDmdt
378  )
379  :
380  fluid_(fluid),
381  pDmdt_(pDmdt)
382  {}
385  {
387  (
388  new populationBalanceModel(fluid_, word(is), pDmdt_)
389  );
390  }
391  };
392 
393 
394  //- Destructor
395  virtual ~populationBalanceModel();
396 
397  // Member Functions
398 
399  //- Dummy write for regIOobject
400  bool writeData(Ostream&) const;
401 
402  //- Return reference to the phaseSystem
403  inline const phaseSystem& fluid() const;
404 
405  //- Return reference to the mesh
406  inline const fvMesh& mesh() const;
407 
408  //- Return populationBalanceCoeffs dictionary
409  inline const dictionary& dict() const;
410 
411  //- Return continuous phase
412  inline const phaseModel& continuousPhase() const;
413 
414  //- Return the velocityGroups belonging to this populationBalance
415  inline const UPtrList<velocityGroup>& velocityGroups() const;
416 
417  //- Return the sizeGroups belonging to this populationBalance
418  inline const UPtrList<sizeGroup>& sizeGroups() const;
419 
420  //- Return list of unordered phasePairs in this populationBalance
421  inline const phasePairTable& phasePairs() const;
422 
423  //- Return implicit source terms
424  inline const volScalarField& Sp(const label i) const;
425 
426  //- Return coalescence relevant size group pairs
427  inline const List<labelPair>& coalescencePairs() const;
428 
429  //- Return binary breakup relevant size group pairs
430  inline const List<labelPair>& binaryBreakupPairs() const;
431 
432  //- Return total void of phases belonging to this populationBalance
433  inline const volScalarField& alphas() const;
434 
435  //- Return average velocity
436  inline const volVectorField& U() const;
437 
438  //- Returns true if both phases are velocity groups and
439  // belong to this populationBalance
440  inline bool isVelocityGroupPair(const phasePair& pair) const;
441 
442  //- Return allocation coefficient
443  const dimensionedScalar eta
444  (
445  const label i,
446  const dimensionedScalar& v
447  ) const;
448 
449  //- Return the surface tension coefficient between a given dispersed
450  // and the continuous phase
452  (
453  const phaseModel& dispersedPhase
454  ) const;
455 
456  //- Return reference to turbulence model of the continuous phase
458  continuousTurbulence() const;
459 
460  //- Solve the population balance equation
461  void solve();
462 
463  //- Correct derived quantities
464  void correct();
465 };
466 
467 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
468 
469 } // End namespace diameterModels
470 } // End namespace Foam
471 
472 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
473 
474 #include "populationBalanceModelI.H"
475 
476 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
477 
478 #endif
479 
480 // ************************************************************************* //
Templated abstract base class for multiphase compressible turbulence models.
Base class for drift models.
Definition: driftModel.H:50
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const word & name() const
Return name.
Definition: IOobject.H:303
const volScalarField & Sp(const label i) const
Return implicit source terms.
iNew(const phaseSystem &fluid, HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > &pDmdt)
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
Class that solves the univariate population balance equation by means of a class method (also called ...
populationBalanceModel(const phaseSystem &fluid, const word &name, HashPtrTable< volScalarField, phasePairKey, phasePairKey::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.
bool isVelocityGroupPair(const phasePair &pair) const
Returns true if both phases are velocity groups and.
const phasePairTable & phasePairs() const
Return list of unordered phasePairs in this populationBalance.
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:330
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.
const UPtrList< velocityGroup > & velocityGroups() const
Return the velocityGroups belonging to this populationBalance.
This class represents a single sizeGroup belonging to a velocityGroup. The main property of a sizeGro...
Definition: sizeGroup.H:100
const phaseModel & continuousPhase() const
Return continuous phase.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:70
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.
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
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 sizeGroups 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:78
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:52
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.
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.