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-2018 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  discretization 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 (aggregation), 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. (2017). 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 reactingMultiphaseEulerFoam 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
57  discretization - I. A fixed pivot technique.
58  Chemical Engineering Science,
59  51(8), 1311-1332.
60  \endverbatim
61 
62  \verbatim
63  Binary breakup term formulation:
64  Liao, Y.; Oertel, R.; Kriebitzsch, S.; Schlegel, F.; Lucas, D.
65  (2018). A Discrete Population Balance Equation for Binary Breakage.
66  International Journal for Numerical Methods in Fluids,
67  DOI 10.1002/fld.4491.
68  \endverbatim
69 
70 Usage
71  Example excerpt from a phaseProperties dictionary.
72  \verbatim
73  type populationBalanceTwoPhaseSystem;
74 
75  phases (air water);
76 
77  populationBalances (bubbles);
78 
79  air
80  {
81  type purePhaseModel;
82  diameterModel velocityGroup;
83  velocityGroupCoeffs
84  {
85  populationBalance bubbles;
86 
87  formFactor 0.5235987756;
88 
89  sizeGroups
90  (
91  f0{d 1.00e-3; value 0;}
92  f1{d 1.08e-3; value 0;}
93  f2{d 1.16e-3; value 0.25;}
94  f3{d 1.25e-3; value 0.5;}
95  f4{d 1.36e-3; value 0.25;}
96  f5{d 1.46e-3; value 0;}
97  ...
98  );
99  }
100 
101  residualAlpha 1e-6;
102  }
103 
104  populationBalanceCoeffs
105  {
106  bubbles
107  {
108  continuousPhase water;
109 
110  coalescenceModels
111  (
112  hydrodynamic
113  {
114  C 0.25;
115  }
116  );
117 
118  binaryBreakupModels
119  ();
120 
121  breakupModels
122  (
123  exponential
124  {
125  C 0.5;
126  exponent 0.01;
127  daughterSizeDistributionModel uniformBinary;
128  }
129  );
130 
131  driftModels
132  (
133  densityChange{}
134  );
135 
136  nucleationModels
137  ();
138  }
139  }
140  \endverbatim
141 
142 See also
143  Foam::diameterModels::sizeGroup
144  Foam::diameterModels::velocityGroup
145 
146 SourceFiles
147  populationBalanceModel.C
148 
149 \*---------------------------------------------------------------------------*/
150 
151 #ifndef populationBalanceModel_H
152 #define populationBalanceModel_H
153 
154 #include "sizeGroup.H"
155 #include "phasePair.H"
156 #include "pimpleControl.H"
157 #include "phaseCompressibleTurbulenceModelFwd.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 typedefs
184 
185  typedef
188 
189 
190  // Private data
191 
192  //- Reference to the phaseSystem
193  const phaseSystem& fluid_;
194 
195  //- Interfacial Mass transfer rate between velocityGroups
197 
198  //- Reference to the mesh
199  const fvMesh& mesh_;
200 
201  //- Name of the populationBalance
202  word name_;
203 
204  //- Dictionary
205  dictionary dict_;
206 
207  //- Reference to pimpleControl
208  const pimpleControl& pimple_;
209 
210  //- Continuous phase
211  const phaseModel& continuousPhase_;
212 
213  //- velocityGroups belonging to this populationBalance
214  List<velocityGroup*> velocityGroups_;
215 
216  //- sizeGroups belonging to this populationBalance
217  List<sizeGroup*> sizeGroups_;
218 
219  //- List of unordered phasePairs in this populationBalance
220  phasePairTable phasePairs_;
221 
222  //- sizeGroup boundaries
224 
225  //- Section width required for binary breakup formulation
227 
228  //- Explicitly treated sources
230 
231  //- Sources treated implicitly or explicitly depending on sign
233 
234  //- Field for caching sources
235  volScalarField Sui_;
236 
237  //- Coalescence models
238  PtrList<coalescenceModel> coalescence_;
239 
240  //- Coalescence rate
241  autoPtr<volScalarField> coalescenceRate_;
242 
243  //- BreakupModels
244  PtrList<breakupModel> breakup_;
245 
246  //- Breakup rate
247  autoPtr<volScalarField> breakupRate_;
248 
249  //- Binary breakup models
250  PtrList<binaryBreakupModel> binaryBreakup_;
251 
252  //- Binary breakup rate
253  autoPtr<volScalarField> binaryBreakupRate_;
254 
255  //- Drift models
256  PtrList<driftModel> drift_;
257 
258  //- Drift rate
259  autoPtr<volScalarField> driftRate_;
260 
261  //- Ratio between successive representative volumes
263 
264  //- Ratio between successive class widths
266 
267  //- Zeroeth order models
268  PtrList<nucleationModel> nucleation_;
269 
270  //- Zeroeth order rate
271  autoPtr<volScalarField> nucleationRate_;
272 
273  //- Total void fraction of phases belonging to this populationBalance
274  volScalarField alphas_;
275 
276  //- Average velocity
277  volVectorField U_;
278 
279 
280  // Private member functions
281 
282  void registerVelocityAndSizeGroups();
283 
284  void add(sizeGroup* group);
285 
286  void createPhasePairs();
287 
288  void correct();
289 
290  void birthByCoalescence(const label j, const label k);
291 
292  void deathByCoalescence(const label i, const label j);
293 
294  void birthByBreakup(const label k, const label model);
295 
296  void deathByBreakup(const label i);
297 
298  void calcDeltas();
299 
300  void birthByBinaryBreakup(const label i, const label j);
301 
302  void deathByBinaryBreakup(const label j, const label i);
303 
304  void drift(const label i);
305 
306  void nucleation(const label i);
307 
308  void sources();
309 
310  void dmdt();
311 
312  void calcAlphas();
313 
314  void calcVelocity();
315 
316 
317 public:
319  friend class sizeGroup;
320  friend class velocityGroup;
321 
322  // Constructor
323 
325  (
326  const phaseSystem& fluid,
327  const word& name,
329  <
331  phasePairKey,
333  >& pDmdt
334  );
335 
336  //- Return clone
338 
339  //- Return a pointer to a new populationBalanceModel object created on
340  // freestore from Istream
341  class iNew
342  {
343  const phaseSystem& fluid_;
344 
346  pDmdt_;
347 
348  public:
349 
351  (
352  const phaseSystem& fluid,
354  pDmdt
355  )
356  :
357  fluid_(fluid),
358  pDmdt_(pDmdt)
359  {}
362  {
364  (
365  new populationBalanceModel(fluid_, word(is), pDmdt_)
366  );
367  }
368  };
369 
370 
371  //- Destructor
372  virtual ~populationBalanceModel();
373 
374  // Member Functions
375 
376  //- Dummy write for regIOobject
377  bool writeData(Ostream&) const;
378 
379  //- Return reference to the phaseSystem
380  inline const phaseSystem& fluid() const;
381 
382  //- Return reference to the mesh
383  inline const fvMesh& mesh() const;
384 
385  //- Return populationBalanceCoeffs dictionary
386  inline const dictionary& dict() const;
387 
388  //- Return continuous phase
389  inline const phaseModel& continuousPhase() const;
390 
391  //- Return the velocityGroups belonging to this populationBalance
392  inline const List<velocityGroup*>& velocityGroups() const;
393 
394  //- Return the sizeGroups belonging to this populationBalance
395  inline const List<sizeGroup*>& sizeGroups() const;
396 
397  //- Return list of unordered phasePairs in this populationBalance
398  inline const phasePairTable& phasePairs() const;
399 
400  //- Return the sizeGroup boundaries
401  inline const PtrList<dimensionedScalar>& v() const;
402 
403  //- Return total void of phases belonging to this populationBalance
404  inline const volScalarField& alphas() const;
405 
406  //- Return average velocity
407  inline const volVectorField& U() const;
408 
409  //- Return allocation coefficient
411  (
412  const label i,
413  const dimensionedScalar& v
414  ) const;
415 
416  //- Return the surface tension coefficient between a given dispersed
417  // and the continuous phase
419  (
420  const phaseModel& dispersedPhase
421  ) const;
422 
423  //- Return reference to turbulence model of the continuous phase
425 
426  //- Solve the population balance equation
427  void solve();
428 };
429 
430 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
431 
432 } // End namespace diameterModels
433 } // End namespace Foam
434 
435 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
436 
437 #include "populationBalanceModelI.H"
438 
439 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
440 
441 #endif
442 
443 // ************************************************************************* //
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:297
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:137
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)
const PtrList< dimensionedScalar > & v() const
Return the sizeGroup boundaries.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
const phasePairTable & phasePairs() const
Return list of unordered phasePairs in this populationBalance.
const dimensionedScalar gamma(const label i, const dimensionedScalar &v) const
Return allocation coefficient.
const dictionary & dict() const
Return populationBalanceCoeffs dictionary.
autoPtr< populationBalanceModel > operator()(Istream &is) const
This diameterModel is intended for use with a populationBalanceModel in order to simulate polydispers...
word group() const
Return group (extension part of name)
Definition: IOobject.C:379
label k
Boltzmann constant.
bool writeData(Ostream &) const
Dummy write for regIOobject.
const List< velocityGroup * > & velocityGroups() const
Return the velocityGroups belonging to this populationBalance.
autoPtr< populationBalanceModel > clone() const
Return clone.
This class represents a single sizeGroup belonging to a velocityGroup. The main property of a sizeGro...
Definition: sizeGroup.H:94
const phaseModel & continuousPhase() const
Return continuous phase.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:67
const phaseSystem & fluid() const
Return reference to the phaseSystem.
A class for handling words, derived from string.
Definition: word.H:59
const List< sizeGroup * > & sizeGroups() const
Return the sizeGroups belonging to this populationBalance.
An STL-conforming hash table.
Definition: HashTable.H:62
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
Return a pointer to a new populationBalanceModel object created on.
const phaseCompressibleTurbulenceModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
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:65
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.
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
Namespace for OpenFOAM.