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