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-2020 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. (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 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 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 
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159 
160 namespace Foam
161 {
162 
163 class phaseSystem;
164 
165 namespace diameterModels
166 {
167 
168 class coalescenceModel;
169 class breakupModel;
170 class binaryBreakupModel;
171 class driftModel;
172 class nucleationModel;
173 
174 /*---------------------------------------------------------------------------*\
175  Class populationBalanceModel Declaration
176 \*---------------------------------------------------------------------------*/
179 :
180  public regIOobject
181 {
182 private:
183 
184  // Private Typedefs
185 
186  typedef
189 
190  typedef
192  pDmdtTable;
193 
194 
195  // Private Data
196 
197  //- Reference to the phaseSystem
198  const phaseSystem& fluid_;
199 
200  //- Interfacial mass transfer rates
201  pDmdtTable& pDmdt_;
202 
203  //- Reference to the mesh
204  const fvMesh& mesh_;
205 
206  //- Name of the populationBalance
207  word name_;
208 
209  //- Dictionary
210  dictionary dict_;
211 
212  //- Reference to pimpleControl
213  const pimpleControl& pimple_;
214 
215  //- Continuous phase
216  const phaseModel& continuousPhase_;
217 
218  //- velocityGroups belonging to this populationBalance
219  UPtrList<velocityGroup> velocityGroups_;
220 
221  //- sizeGroups belonging to this populationBalance
222  UPtrList<sizeGroup> sizeGroups_;
223 
224  //- List of unordered phasePairs in this populationBalance
225  phasePairTable phasePairs_;
226 
227  //- sizeGroup boundaries
229 
230  //- Section width required for binary breakup formulation
232 
233  //- Explicitly treated sources
235 
236  //- Sources treated implicitly or explicitly depending on sign
238 
239  //- Field for caching sources
240  volScalarField Sui_;
241 
242  //- Coalescence models
243  PtrList<coalescenceModel> coalescence_;
244 
245  //- Coalescence rate
246  autoPtr<volScalarField> coalescenceRate_;
247 
248  //- BreakupModels
249  PtrList<breakupModel> breakup_;
250 
251  //- Breakup rate
252  autoPtr<volScalarField> breakupRate_;
253 
254  //- Binary breakup models
255  PtrList<binaryBreakupModel> binaryBreakup_;
256 
257  //- Binary breakup rate
258  autoPtr<volScalarField> binaryBreakupRate_;
259 
260  //- Drift models
261  PtrList<driftModel> drift_;
262 
263  //- Drift rate
264  autoPtr<volScalarField> driftRate_;
265 
266  //- Ratio between successive representative volumes
268 
269  //- Ratio between successive class widths
271 
272  //- Zeroeth order models
273  PtrList<nucleationModel> nucleation_;
274 
275  //- Zeroeth order rate
276  autoPtr<volScalarField> nucleationRate_;
277 
278  //- Total void fraction
279  autoPtr<volScalarField> alphas_;
280 
281  //- Mean Sauter diameter
283 
284  //- Average velocity
286 
287  //- Counter for interval between source term updates
288  label sourceUpdateCounter_;
289 
290 
291  // Private Member Functions
292 
293  void registerVelocityGroups();
294 
295  void registerSizeGroups(sizeGroup& group);
296 
297  void createPhasePairs();
298 
299  void correct();
300 
301  void birthByCoalescence(const label j, const label k);
302 
303  void deathByCoalescence(const label i, const label j);
304 
305  void birthByBreakup(const label k, const label model);
306 
307  void deathByBreakup(const label i);
308 
309  void calcDeltas();
310 
311  void birthByBinaryBreakup(const label i, const label j);
312 
313  void deathByBinaryBreakup(const label j, const label i);
314 
315  void drift(const label i, driftModel& model);
316 
317  void nucleation(const label i, nucleationModel& model);
318 
319  void sources();
320 
321  void dmdt();
322 
323  void calcAlphas();
324 
325  tmp<volScalarField> calcDsm();
326 
327  void calcVelocity();
328 
329  //- Return whether the sources should be updated on this iteration
330  bool updateSources();
331 
332  //- Return the number of corrections
333  inline label nCorr() const;
334 
335  //- Return the interval at which the sources are updated
336  inline label sourceUpdateInterval() const;
337 
338 public:
339 
340  //- Runtime type information
341  TypeName("populationBalanceModel");
342 
343 
344  // Constructor
345 
347  (
348  const phaseSystem& fluid,
349  const word& name,
351  <
353  phasePairKey,
355  >& pDmdt
356  );
357 
358  //- Return clone
360 
361  //- Return a pointer to a new populationBalanceModel object created on
362  // freestore from Istream
363  class iNew
364  {
365  const phaseSystem& fluid_;
366 
368  pDmdt_;
369 
370  public:
371 
373  (
374  const phaseSystem& fluid,
376  pDmdt
377  )
378  :
379  fluid_(fluid),
380  pDmdt_(pDmdt)
381  {}
384  {
386  (
387  new populationBalanceModel(fluid_, word(is), pDmdt_)
388  );
389  }
390  };
391 
392 
393  //- Destructor
394  virtual ~populationBalanceModel();
395 
396  // Member Functions
397 
398  //- Dummy write for regIOobject
399  bool writeData(Ostream&) const;
400 
401  //- Return reference to the phaseSystem
402  inline const phaseSystem& fluid() const;
403 
404  //- Return reference to the mesh
405  inline const fvMesh& mesh() const;
406 
407  //- Return populationBalanceCoeffs dictionary
408  inline const dictionary& dict() const;
409 
410  //- Return continuous phase
411  inline const phaseModel& continuousPhase() const;
412 
413  //- Return the velocityGroups belonging to this populationBalance
414  inline const UPtrList<velocityGroup>& velocityGroups() const;
415 
416  //- Return the sizeGroups belonging to this populationBalance
417  inline const UPtrList<sizeGroup>& sizeGroups() const;
418 
419  //- Return semi-implicit source terms
420  inline const volScalarField& SuSp(const label i) const;
421 
422  //- Return list of unordered phasePairs in this populationBalance
423  inline const phasePairTable& phasePairs() const;
424 
425  //- Returns true if both phases are velocity groups and
426  // belong to this populationBalance
427  inline bool isVelocityGroupPair(const phasePair& pair) const;
428 
429  //- Return the sizeGroup boundaries
430  inline const PtrList<dimensionedScalar>& v() 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  //- Return allocation coefficient
439  const dimensionedScalar eta
440  (
441  const label i,
442  const dimensionedScalar& v
443  ) const;
444 
445  //- Return the surface tension coefficient between a given dispersed
446  // and the continuous phase
448  (
449  const phaseModel& dispersedPhase
450  ) const;
451 
452  //- Return reference to turbulence model of the continuous phase
454  continuousTurbulence() const;
455 
456  //- Solve the population balance equation
457  void solve();
458 };
459 
460 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
461 
462 } // End namespace diameterModels
463 } // End namespace Foam
464 
465 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
466 
467 #include "populationBalanceModelI.H"
468 
469 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
470 
471 #endif
472 
473 // ************************************************************************* //
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
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:158
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.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
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 volScalarField & SuSp(const label i) const
Return semi-implicit source terms.
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:340
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:67
const phaseSystem & fluid() const
Return reference to the phaseSystem.
A class for handling words, derived from string.
Definition: word.H:59
const phaseCompressibleMomentumTransportModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
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:55
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
TypeName("populationBalanceModel")
Runtime type information.
Namespace for OpenFOAM.