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