populationBalanceModel.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2017-2019 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::diameterModels::populationBalanceModel
28
29Description
30 Class that solves the univariate population balance equation by means of
31 a class method (also called sectional or discrete method). The internal
32 coordinate is set to the particle volume, so the equation is based on
33 a transport equation of the volume-based number density function. The
34 discretization is done using the fixed pivot technique of Kumar and
35 Ramkrishna (1996). The source terms are written in a way that particle
36 number and mass are preserved. Coalescence (aggregation), breakup, drift
37 (growth and surface loss) as well as nucleation are supported.
38 For the discrete breakup term two recipes are available, depending on the
39 model choice. For models which state a total breakup rate and a separate
40 daughter size distribution function, the formulation of Kumar and Ramkrishna
41 (1996) is applied which is applicable for binary and multiple breakup
42 events. The second formulation is given by Liao et al. (2018). It is useful
43 for binary breakup models which give the breakup rate between a sizeGroup
44 pair directly, without an explicit expression for the daughter size
45 distribution. The drift term is implemented using a finite difference upwind
46 scheme. Although it is diffusive, it ensures a stable and
47 number-conservative solution.
48
49 The implementation allows to split the population balance over multiple
50 velocity fields using the capability of reactingMultiphaseEulerFoam to solve
51 for n momentum equations. It is also possible to define multiple population
52 balances, e.g. bubbles and droplets simultaneously.
53
54 References:
55 \verbatim
56 Coalescence and breakup term formulation:
57 Kumar, S., & Ramkrishna, D. (1996).
58 On the solution of population balance equations by discretization-I. A
59 fixed pivot technique.
60 Chemical Engineering Science, 51(8), 1311-1332.
61 \endverbatim
62
63 \verbatim
64 Binary breakup term formulation:
65 Liao, Y., Oertel, R., Kriebitzsch, S., Schlegel, F., & Lucas, D. (2018).
66 A discrete population balance equation for binary breakage.
67 International Journal for Numerical Methods in Fluids, 87(4), 202-215.
68 \endverbatim
69
70Usage
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
142See also
143 Foam::diameterModels::sizeGroup
144 Foam::diameterModels::velocityGroup
145
146SourceFiles
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#include "HashPtrTable.H"
159
160// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
161
162namespace Foam
163{
164
165class phaseSystem;
166
167namespace diameterModels
168{
169
170class coalescenceModel;
171class breakupModel;
172class binaryBreakupModel;
173class driftModel;
174class nucleationModel;
175
176/*---------------------------------------------------------------------------*\
177 Class populationBalanceModel Declaration
178\*---------------------------------------------------------------------------*/
181:
182 public regIOobject
183{
184 // Private typedefs
185
186 typedef
189
190
191 // Private data
192
193 //- Reference to the phaseSystem
194 const phaseSystem& fluid_;
195
196 //- Interfacial Mass transfer rate between velocityGroups
198
199 //- Reference to the mesh
200 const fvMesh& mesh_;
201
202 //- Name of the populationBalance
203 word name_;
204
205 //- Dictionary
206 dictionary dict_;
207
208 //- Reference to pimpleControl
209 const pimpleControl& pimple_;
210
211 //- Continuous phase
212 const phaseModel& continuousPhase_;
213
214 //- velocityGroups belonging to this populationBalance
215 UPtrList<velocityGroup> velocityGroups_;
216
217 //- sizeGroups belonging to this populationBalance
218 UPtrList<sizeGroup> sizeGroups_;
219
220 //- List of unordered phasePairs in this populationBalance
221 phasePairTable phasePairs_;
222
223 //- sizeGroup boundaries
225
226 //- Section width required for binary breakup formulation
228
229 //- Explicitly treated sources
231
232 //- Sources treated implicitly or explicitly depending on sign
234
235 //- Field for caching sources
236 volScalarField Sui_;
237
238 //- Coalescence models
239 PtrList<coalescenceModel> coalescence_;
240
241 //- Coalescence rate
242 autoPtr<volScalarField> coalescenceRate_;
243
244 //- BreakupModels
245 PtrList<breakupModel> breakup_;
246
247 //- Breakup rate
248 autoPtr<volScalarField> breakupRate_;
249
250 //- Binary breakup models
251 PtrList<binaryBreakupModel> binaryBreakup_;
252
253 //- Binary breakup rate
254 autoPtr<volScalarField> binaryBreakupRate_;
255
256 //- Drift models
257 PtrList<driftModel> drift_;
258
259 //- Drift rate
260 autoPtr<volScalarField> driftRate_;
261
262 //- Ratio between successive representative volumes
264
265 //- Ratio between successive class widths
267
268 //- Zeroeth order models
269 PtrList<nucleationModel> nucleation_;
270
271 //- Zeroeth order rate
272 autoPtr<volScalarField> nucleationRate_;
273
274 //- Total void fraction
276
277 //- Mean Sauter diameter
279
280 //- Average velocity
282
283 //- Counter for interval between source term updates
284 label sourceUpdateCounter_;
285
286
287 // Private member functions
288
289 void registerVelocityGroups();
290
291 void registerSizeGroups(sizeGroup& group);
292
293 void createPhasePairs();
294
295 void correct();
296
297 void birthByCoalescence(const label j, const label k);
298
299 void deathByCoalescence(const label i, const label j);
300
301 void birthByBreakup(const label k, const label model);
302
303 void deathByBreakup(const label i);
304
305 void calcDeltas();
306
307 void birthByBinaryBreakup(const label i, const label j);
308
309 void deathByBinaryBreakup(const label j, const label i);
310
311 void drift(const label i);
312
313 void nucleation(const label i);
314
315 void sources();
316
317 void dmdt();
318
319 void calcAlphas();
320
321 tmp<volScalarField> calcDsm();
322
323 void calcVelocity();
324
325 //- Return whether the sources should be updated on this iteration
326 bool updateSources();
327
328 //- Return the number of corrections
329 inline label nCorr() const;
330
331 //- Return the interval at which the sources are updated
332 inline label sourceUpdateInterval() const;
333
334public:
336 friend class sizeGroup;
337 friend class velocityGroup;
338
339 // Constructor
340
342 (
343 const phaseSystem& fluid,
344 const word& name,
346 <
350 >& pDmdt
351 );
352
353 //- Return clone
355
356 //- Return a pointer to a new populationBalanceModel object created on
357 // freestore from Istream
358 class iNew
359 {
360 const phaseSystem& fluid_;
361
363 pDmdt_;
364
365 public:
367 iNew
368 (
369 const phaseSystem& fluid,
371 pDmdt
372 )
373 :
374 fluid_(fluid),
375 pDmdt_(pDmdt)
376 {}
379 {
381 (
382 new populationBalanceModel(fluid_, word(is), pDmdt_)
383 );
384 }
385 };
386
387
388 //- Destructor
389 virtual ~populationBalanceModel();
390
391 // Member Functions
392
393 //- Dummy write for regIOobject
394 bool writeData(Ostream&) const;
395
396 //- Return reference to the phaseSystem
397 inline const phaseSystem& fluid() const;
398
399 //- Return reference to the mesh
400 inline const fvMesh& mesh() const;
401
402 //- Return populationBalanceCoeffs dictionary
403 inline const dictionary& dict() const;
404
405 //- Return continuous phase
406 inline const phaseModel& continuousPhase() const;
407
408 //- Return the velocityGroups belonging to this populationBalance
409 inline const UPtrList<velocityGroup>& velocityGroups() const;
410
411 //- Return the sizeGroups belonging to this populationBalance
412 inline const UPtrList<sizeGroup>& sizeGroups() const;
413
414 //- Return list of unordered phasePairs in this populationBalance
415 inline const phasePairTable& phasePairs() const;
416
417 //- Return the sizeGroup boundaries
418 inline const PtrList<dimensionedScalar>& v() const;
419
420 //- Return total void of phases belonging to this populationBalance
421 inline const volScalarField& alphas() const;
422
423 //- Return average velocity
424 inline const volVectorField& U() const;
425
426 //- Return allocation coefficient
428 (
429 const label i,
430 const dimensionedScalar& v
431 ) const;
432
433 //- Return the surface tension coefficient between a given dispersed
434 // and the continuous phase
436 (
437 const phaseModel& dispersedPhase
438 ) const;
439
440 //- Return reference to turbulence model of the continuous phase
442
443 //- Solve the population balance equation
444 void solve();
445};
446
447// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
448
449} // End namespace diameterModels
450} // End namespace Foam
451
452// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
453
455
456// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
457
458#endif
459
460// ************************************************************************* //
label k
A HashTable of pointers to objects of type <T>, with deallocation management of the pointers.
Definition: HashPtrTable.H:68
A HashTable similar to std::unordered_map.
Definition: HashTable.H:123
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
word group() const
Return group (extension part of name)
Definition: IOobjectI.H:71
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:71
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Return a pointer to a new populationBalanceModel object created on.
iNew(const phaseSystem &fluid, HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > &pDmdt)
autoPtr< populationBalanceModel > operator()(Istream &is) const
Class that solves the univariate population balance equation by means of a class method (also called ...
const phaseSystem & fluid() const
Return reference to the phaseSystem.
bool writeData(Ostream &) const
Dummy write for regIOobject.
autoPtr< populationBalanceModel > clone() const
Return clone.
const volVectorField & U() const
Return average velocity.
const phasePairTable & phasePairs() const
Return list of unordered phasePairs in this populationBalance.
const UPtrList< sizeGroup > & sizeGroups() const
Return the sizeGroups belonging to this populationBalance.
const phaseModel & continuousPhase() const
Return continuous phase.
const volScalarField & alphas() const
Return total void of phases belonging to this populationBalance.
const dictionary & dict() const
Return populationBalanceCoeffs dictionary.
const phaseCompressibleTurbulenceModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
const tmp< volScalarField > sigmaWithContinuousPhase(const phaseModel &dispersedPhase) const
Return the surface tension coefficient between a given dispersed.
const PtrList< dimensionedScalar > & v() const
Return the sizeGroup boundaries.
const fvMesh & mesh() const
Return reference to the mesh.
void solve()
Solve the population balance equation.
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:99
This diameterModel is intended for use with a populationBalanceModel in order to simulate polydispers...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
An ordered or unorder pair of phase names. Typically specified as follows.
Definition: phasePairKey.H:68
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:76
PIMPLE control class to supply convergence information/checks for the PIMPLE loop.
Definition: pimpleControl.H:58
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:76
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
const scalar gamma
Definition: EEqn.H:9
Namespace for OpenFOAM.
Hashing functor for phasePairKey.
Definition: phasePairKey.H:123