PrinceBlanch.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) 2018-2019 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
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 
26 Class
27  Foam::diameterModels::coalescenceModels::PrinceBlanch
28 
29 Description
30  Model of Prince and Blanch (1990). The coalescence rate is calculated by
31 
32  \f[
33  \left( \theta_{ij}^{T} + \theta_{ij}^{B} + \theta_{ij}^{LS} \right)
34  \lambda_{ij}
35  \f]
36 
37  with the coalescence efficiency
38 
39  \f[
40  \lambda_{ij} =
41  \mathrm{exp}
42  \left(
43  - \sqrt{\frac{r_{ij}^3 \rho_c}{16 \sigma}}
44  \mathrm{ln} \left(\frac{h_0}{h_f}\right)
45  \epsilon_c^{1/3}/r_{ij}^{2/3}
46  \right)\;,
47  \f]
48 
49  the turbulent collision rate
50 
51  \f[
52  \theta_{ij}^{T} =
53  C_1 \pi (d_i + d_j)^{2} \epsilon_c^{1/3}
54  \sqrt{d_{i}^{2/3} + d_{j}^{2/3}}\;,
55  \f]
56 
57  and the buoyancy-driven collision rate
58 
59  \f[
60  \theta_{ij}^{B} = S_{ij} \left| u_{ri} - u_{rj} \right|\;,
61  \f]
62 
63  where the rise velocity of bubble i is calculated by
64 
65  \f[
66  u_{ri} = \sqrt{2.14 \sigma / \left(\rho_c d_i \right) + 0.505 g d_i}\;,
67  \f]
68 
69  the equivalent radius by
70 
71  \f[
72  r_{ij} = \left( \frac{1}{d_i} + \frac{1}{d_j} \right)^{-1}
73  \f]
74 
75  and the collision cross sectional area by
76 
77  \f[
78  S_{ij} = \frac{\pi}{4} \left(d_i + d_j\right)^{2}\;.
79  \f]
80 
81  Note that in equation 2, the bubble radius has been substituted by the
82  bubble diameter. Also the expression for the equivalent radius r_ij
83  (equation 19 in the paper of Prince and Blanch (1990)) was corrected.
84  The collision rate contribution due to laminar shear in the continuous phase
85  is currently neglected.
86 
87  \vartable
88  \theta_{ij}^{T} | Turbulent collision rate [m3/s]
89  \theta_{ij}^{B} | Buoyancy-driven collision rate [m3/s]
90  \theta_{ij}^{LS} | Laminar shear collision rate [m3/s]
91  \lambda_{ij} | Coalescence efficiency [-]
92  r_{ij} | Equivalent radius [m]
93  \rho_c | Density of continuous phase [kg/m3]
94  \sigma | Surface tension [N/m]
95  h_0 | Initial film thickness [m]
96  h_f | Critical film thickness [m]
97  \epsilon_c | Continuous phase turbulent dissipation rate [m2/s3]
98  d_i | Diameter of bubble i [m]
99  d_j | Diameter of bubble j [m]
100  u_{ri} | Rise velocity of bubble i [m/s]
101  S_{ij} | Collision cross sectional area [m2]
102  g | Gravitational constant [m/s2]
103  \endvartable
104 
105  Reference:
106  \verbatim
107  Prince, M. J., & Blanch, H. W. (1990).
108  Bubble coalescence and break‐up in air‐sparged bubble columns.
109  AIChE journal, 36(10), 1485-1499.
110  \endverbatim
111 
112 Usage
113  \table
114  Property | Description | Required | Default value
115  C1 | coefficient C1 | no | 0.089
116  h0 | initial film thickness | no | 1e-4m
117  hf | critical film thickness | no | 1e-8m
118  turbulence | Switch for collisions due to turbulence | yes | none
119  buoyancy | Switch for collisions due to buoyancy | yes | none
120  laminarShear | Switch for collisions due to laminar shear | yes | none
121  \endtable
122 
123 SourceFiles
124  PrinceBlanch.C
125 
126 \*---------------------------------------------------------------------------*/
127 
128 #ifndef PrinceBlanch_H
129 #define PrinceBlanch_H
130 
131 #include "coalescenceModel.H"
132 
133 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
134 
135 namespace Foam
136 {
137 namespace diameterModels
138 {
139 namespace coalescenceModels
140 {
141 
142 /*---------------------------------------------------------------------------*\
143  Class PrinceBlanch Declaration
144 \*---------------------------------------------------------------------------*/
145 
146 class PrinceBlanch
147 :
148  public coalescenceModel
149 {
150  // Private data
151 
152  //- Optional coefficient C1, defaults to 0.089
153  dimensionedScalar C1_;
154 
155  //- Initial film thickness, defaults to 1e-4m
156  dimensionedScalar h0_;
157 
158  //- Critical film thickness, defaults to 1e-8m
159  dimensionedScalar hf_;
160 
161  //- Switch for considering turbulent collisions
162  Switch turbulence_;
163 
164  //- Switch for considering buoyancy-induced collisions
165  Switch buoyancy_;
166 
167  //- Switch for considering buoyancy-induced collisions
168  Switch laminarShear_;
169 
170 
171 public:
172 
173  //- Runtime type information
174  TypeName("PrinceBlanch");
175 
176  // Constructor
177 
179  (
180  const populationBalanceModel& popBal,
181  const dictionary& dict
182  );
183 
184 
185  //- Destructor
186  virtual ~PrinceBlanch() = default;
187 
188 
189  // Member Functions
190 
191  //- Add to coalescenceRate
192  virtual void addToCoalescenceRate
193  (
194  volScalarField& coalescenceRate,
195  const label i,
196  const label j
197  );
198 };
199 
200 
201 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
202 
203 } // End namespace coalescenceModels
204 } // End namespace diameterModels
205 } // End namespace Foam
206 
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208 
209 #endif
210 
211 // ************************************************************************* //
Foam::diameterModels::coalescenceModels::PrinceBlanch
Model of Prince and Blanch (1990). The coalescence rate is calculated by.
Definition: PrinceBlanch.H:240
Foam::diameterModels::coalescenceModel
Base class for coalescence models.
Definition: coalescenceModel.H:52
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:77
Foam::diameterModels::coalescenceModels::PrinceBlanch::~PrinceBlanch
virtual ~PrinceBlanch()=default
Destructor.
Foam::diameterModels::coalescenceModels::PrinceBlanch::addToCoalescenceRate
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
Definition: PrinceBlanch.C:96
coalescenceModel.H
Foam::diameterModels::coalescenceModels::PrinceBlanch::PrinceBlanch
PrinceBlanch(const populationBalanceModel &popBal, const dictionary &dict)
Definition: PrinceBlanch.C:59
Foam::diameterModels::coalescenceModels::PrinceBlanch::TypeName
TypeName("PrinceBlanch")
Runtime type information.
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::diameterModels::populationBalanceModel
Class that solves the univariate population balance equation by means of a class method (also called ...
Definition: populationBalanceModel.H:179
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::GeometricField< scalar, fvPatchField, volMesh >