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-2019 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  formFactor 0.5235987756;
86 
87  sizeGroups
88  (
89  f0{d 1.00e-3; value 0;}
90  f1{d 1.08e-3; value 0;}
91  f2{d 1.16e-3; value 0.25;}
92  f3{d 1.25e-3; value 0.5;}
93  f4{d 1.36e-3; value 0.25;}
94  f5{d 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 uniformBinary;
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"
155 #include "phaseCompressibleTurbulenceModelFwd.H"
156 
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 
159 namespace Foam
160 {
161 
162 class phaseSystem;
163 
164 namespace diameterModels
165 {
166 
167 class coalescenceModel;
168 class breakupModel;
169 class binaryBreakupModel;
170 class driftModel;
171 class nucleationModel;
172 
173 /*---------------------------------------------------------------------------*\
174  Class populationBalanceModel Declaration
175 \*---------------------------------------------------------------------------*/
178 :
179  public regIOobject
180 {
181  // Private Typedefs
182 
183  typedef
186 
187 
188  // Private Data
189 
190  //- Reference to the phaseSystem
191  const phaseSystem& fluid_;
192 
193  //- Interfacial Mass transfer rate between velocityGroups
195 
196  //- Reference to the mesh
197  const fvMesh& mesh_;
198 
199  //- Name of the populationBalance
200  word name_;
201 
202  //- Dictionary
203  dictionary dict_;
204 
205  //- Reference to pimpleControl
206  const pimpleControl& pimple_;
207 
208  //- Continuous phase
209  const phaseModel& continuousPhase_;
210 
211  //- velocityGroups belonging to this populationBalance
212  UPtrList<velocityGroup> velocityGroups_;
213 
214  //- sizeGroups belonging to this populationBalance
215  UPtrList<sizeGroup> sizeGroups_;
216 
217  //- List of unordered phasePairs in this populationBalance
218  phasePairTable phasePairs_;
219 
220  //- sizeGroup boundaries
222 
223  //- Section width required for binary breakup formulation
225 
226  //- Explicitly treated sources
228 
229  //- Sources treated implicitly or explicitly depending on sign
231 
232  //- Field for caching sources
233  volScalarField Sui_;
234 
235  //- Coalescence models
236  PtrList<coalescenceModel> coalescence_;
237 
238  //- Coalescence rate
239  autoPtr<volScalarField> coalescenceRate_;
240 
241  //- BreakupModels
242  PtrList<breakupModel> breakup_;
243 
244  //- Breakup rate
245  autoPtr<volScalarField> breakupRate_;
246 
247  //- Binary breakup models
248  PtrList<binaryBreakupModel> binaryBreakup_;
249 
250  //- Binary breakup rate
251  autoPtr<volScalarField> binaryBreakupRate_;
252 
253  //- Drift models
254  PtrList<driftModel> drift_;
255 
256  //- Drift rate
257  autoPtr<volScalarField> driftRate_;
258 
259  //- Ratio between successive representative volumes
261 
262  //- Ratio between successive class widths
264 
265  //- Zeroeth order models
266  PtrList<nucleationModel> nucleation_;
267 
268  //- Zeroeth order rate
269  autoPtr<volScalarField> nucleationRate_;
270 
271  //- Total void fraction
272  autoPtr<volScalarField> alphas_;
273 
274  //- Mean Sauter diameter
276 
277  //- Average velocity
279 
280  //- Counter for interval between source term updates
281  label sourceUpdateCounter_;
282 
283 
284  // Private Member Functions
285 
286  void registerVelocityGroups();
287 
288  void registerSizeGroups(sizeGroup& group);
289 
290  void createPhasePairs();
291 
292  void correct();
293 
294  void birthByCoalescence(const label j, const label k);
295 
296  void deathByCoalescence(const label i, const label j);
297 
298  void birthByBreakup(const label k, const label model);
299 
300  void deathByBreakup(const label i);
301 
302  void calcDeltas();
303 
304  void birthByBinaryBreakup(const label i, const label j);
305 
306  void deathByBinaryBreakup(const label j, const label i);
307 
308  void drift(const label i);
309 
310  void nucleation(const label i);
311 
312  void sources();
313 
314  void dmdt();
315 
316  void calcAlphas();
317 
318  tmp<volScalarField> calcDsm();
319 
320  void calcVelocity();
321 
322  //- Return whether the sources should be updated on this iteration
323  bool updateSources();
324 
325  //- Return the number of corrections
326  inline label nCorr() const;
327 
328  //- Return the interval at which the sources are updated
329  inline label sourceUpdateInterval() const;
330 
331 public:
333  friend class sizeGroup;
334  friend class velocityGroup;
335 
336  // Constructor
337 
339  (
340  const phaseSystem& fluid,
341  const word& name,
343  <
345  phasePairKey,
347  >& pDmdt
348  );
349 
350  //- Return clone
352 
353  //- Return a pointer to a new populationBalanceModel object created on
354  // freestore from Istream
355  class iNew
356  {
357  const phaseSystem& fluid_;
358 
360  pDmdt_;
361 
362  public:
363 
365  (
366  const phaseSystem& fluid,
368  pDmdt
369  )
370  :
371  fluid_(fluid),
372  pDmdt_(pDmdt)
373  {}
376  {
378  (
379  new populationBalanceModel(fluid_, word(is), pDmdt_)
380  );
381  }
382  };
383 
384 
385  //- Destructor
386  virtual ~populationBalanceModel();
387 
388  // Member Functions
389 
390  //- Dummy write for regIOobject
391  bool writeData(Ostream&) const;
392 
393  //- Return reference to the phaseSystem
394  inline const phaseSystem& fluid() const;
395 
396  //- Return reference to the mesh
397  inline const fvMesh& mesh() const;
398 
399  //- Return populationBalanceCoeffs dictionary
400  inline const dictionary& dict() const;
401 
402  //- Return continuous phase
403  inline const phaseModel& continuousPhase() const;
404 
405  //- Return the velocityGroups belonging to this populationBalance
406  inline const UPtrList<velocityGroup>& velocityGroups() const;
407 
408  //- Return the sizeGroups belonging to this populationBalance
409  inline const UPtrList<sizeGroup>& sizeGroups() const;
410 
411  //- Return list of unordered phasePairs in this populationBalance
412  inline const phasePairTable& phasePairs() const;
413 
414  //- Return the sizeGroup boundaries
415  inline const PtrList<dimensionedScalar>& v() const;
416 
417  //- Return total void of phases belonging to this populationBalance
418  inline const volScalarField& alphas() const;
419 
420  //- Return average velocity
421  inline const volVectorField& U() const;
422 
423  //- Return allocation coefficient
425  (
426  const label i,
427  const dimensionedScalar& v
428  ) const;
429 
430  //- Return the surface tension coefficient between a given dispersed
431  // and the continuous phase
433  (
434  const phaseModel& dispersedPhase
435  ) const;
436 
437  //- Return reference to turbulence model of the continuous phase
439 
440  //- Solve the population balance equation
441  void solve();
442 };
443 
444 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
445 
446 } // End namespace diameterModels
447 } // End namespace Foam
448 
449 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
450 
451 #include "populationBalanceModelI.H"
452 
453 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
454 
455 #endif
456 
457 // ************************************************************************* //
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:295
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
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:372
label k
Boltzmann constant.
bool writeData(Ostream &) const
Dummy write for regIOobject.
autoPtr< populationBalanceModel > clone() const
Return clone.
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: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
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:53
Return a pointer to a new populationBalanceModel object created on.
const phaseCompressibleTurbulenceModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
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
Namespace for OpenFOAM.