LuoSvendsen.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) 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::binaryBreakupModels::LuoSvendsen
28 
29 Description
30  Model of Luo and Svendsen (1996). The breakup rate is calculated by
31 
32  \f[
33  C_4 \alpha_c \left(\frac{\epsilon_c}{d_j^2}\right)^{1/3}
34  \int\limits_{\xi_{min}}^{1}
35  \frac{\left(1 + \xi\right)^{2}}{\xi^{11/3}}
36  \mathrm{exp}
37  \left(
38  - \frac{12c_f\sigma}{\beta\rho_c\epsilon_c^{2/3}d_j^{5/3}\xi^{11/3}}
39  \right)
40  \mathrm{d} \xi
41  \f]
42 
43  where
44 
45  \f[
46  c_f = \left(\frac{v_i}{v_j}\right)^{2/3}
47  + \left(1 - \frac{v_i}{v_j}\right)^{2/3} - 1
48  \f]
49 
50  \f[
51  \xi_{min} = \frac{\lambda_{min}}{d_j}\,,
52  \f]
53 
54  and
55 
56  \f[
57  \lambda_{min} = C_5 \eta\,.
58  \f]
59 
60  The integral in the first expression is solved by means of incomplete Gamma
61  functions as given by Bannari et al. (2008):
62 
63  \f[
64  \frac{3}{11 b^{8/11}}
65  \left(
66  \left[\Gamma(8/11, b) - \Gamma(8/11, t_{m})\right]
67  + 2b^{3/11} \left[\Gamma(5/11, b) - \Gamma(5/11, t_{m})\right]
68  + b^{6/11} \left[\Gamma(2/11, b) - \Gamma(2/11, t_{m})\right]
69  \right)
70  \f]
71 
72  where
73 
74  \f[
75  b = \frac{12c_f\sigma}{\beta\rho_c\epsilon_c^{2/3}d_j^{5/3}}
76  \f]
77 
78  and
79 
80  \f[
81  t_{min} = b \xi_{min}^{-11/3}\,.
82  \f]
83 
84  Note that in the code, the upper incomplete gamma function is expressed as
85 
86  \f[
87  \Gamma(a,z) = Q(a,z) \Gamma(a)
88  \f]
89 
90  \vartable
91  \alpha_c | Void fraction of continuous phase [-]
92  \epsilon_c | Turbulent dissipation rate of continuous phase [m2/s3]
93  d_j | Diameter of mother bubble j [m3]
94  v_i | Volume of daughter bubble i [m3]
95  v_j | Volume of mother bubble j [m3]
96  \xi | Integration variable [-]
97  \xi_{min} | Lower bound of integral [-]
98  c_f | Increase coefficient of surface area [-]
99  \sigma | Surface tension [N/m]
100  \rho_c | Density of continuous phase [kg/m3]
101  \eta | Kolmogorov length scale [m]
102  \Gamma(a,z) | Upper incomplete gamma function
103  Q(a,z) | Regularized upper incomplete gamma function
104  \Gamma(a) | Gamma function
105  \endvartable
106 
107  References:
108  \verbatim
109  Luo, H., & Svendsen, H. F. (1996).
110  Theoretical model for drop and bubble breakup in turbulent dispersions.
111  AIChE Journal, 42(5), 1225-1233.
112  Eq. 27, p. 1229.
113  \endverbatim
114 
115  \verbatim
116  Bannari, R., Kerdouss, F., Selma, B., Bannari, A., & Proulx, P. (2008).
117  Three-dimensional mathematical modeling of dispersed two-phase flow
118  using class method of population balance in bubble columns.
119  Computers & chemical engineering, 32(12), 3224-3237.
120  Eq. 49, p. 3230.
121  \endverbatim
122 
123 Usage
124  \table
125  Property | Description | Required | Default value
126  C4 | Coefficient C4 | no | 0.923
127  beta | Coefficient beta | no | 2.05
128  C5 | Minimum eddy ratio | no | 11.4
129  \endtable
130 
131 SourceFiles
132  LuoSvendsen.C
133 
134 \*---------------------------------------------------------------------------*/
135 
136 #ifndef LuoSvendsen_H
137 #define LuoSvendsen_H
138 
139 #include "binaryBreakupModel.H"
140 #include "interpolationTable.H"
141 
142 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
143 
144 namespace Foam
145 {
146 namespace diameterModels
147 {
148 namespace binaryBreakupModels
149 {
150 
151 /*---------------------------------------------------------------------------*\
152  Class LuoSvendsen Declaration
153 \*---------------------------------------------------------------------------*/
154 
155 class LuoSvendsen
156 :
157  public binaryBreakupModel
158 {
159 private:
160 
161  // Private data
162 
163  //- Interpolation table of Q(a,z) for a=2/11
164  autoPtr<interpolationTable<scalar>> gammaUpperReg2by11_;
165 
166  //- Interpolation table of Q(a,z) for a=5/11
167  autoPtr<interpolationTable<scalar>> gammaUpperReg5by11_;
168 
169  //- Interpolation table of Q(a,z) for a=8/11
170  autoPtr<interpolationTable<scalar>> gammaUpperReg8by11_;
171 
172  //- Empirical constant, defaults to 0.923
173  dimensionedScalar C4_;
174 
175  //- Empirical constant, defaults to 2.05
176  dimensionedScalar beta_;
177 
178  //- Ratio between minimum size of eddies in the inertial subrange
179  // and Kolmogorov length scale, defaults to 11.4
180  dimensionedScalar minEddyRatio_;
181 
182  //- Kolmogorov length scale
183  volScalarField kolmogorovLengthScale_;
184 
185 
186 public:
187 
188  //- Runtime type information
189  TypeName("LuoSvendsen");
190 
191  // Constructor
192 
194  (
195  const populationBalanceModel& popBal,
196  const dictionary& dict
197  );
198 
199 
200  //- Destructor
201  virtual ~LuoSvendsen() = default;
202 
203 
204  // Member Functions
205 
206  //- Correct diameter independent expressions
207  virtual void correct();
208 
209  //- Add to binary breakupRate
210  virtual void addToBinaryBreakupRate
211  (
212  volScalarField& binaryBreakupRate,
213  const label i,
214  const label j
215  );
216 };
217 
218 
219 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220 
221 } // End namespace binaryBreakupModels
222 } // End namespace diameterModels
223 } // End namespace Foam
224 
225 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
227 #endif
228 
229 // ************************************************************************* //
dictionary dict
LuoSvendsen(const populationBalanceModel &popBal, const dictionary &dict)
Definition: LuoSvendsen.C:51
virtual void addToBinaryBreakupRate(volScalarField &binaryBreakupRate, const label i, const label j)
Add to binary breakupRate.
Definition: LuoSvendsen.C:166
TypeName("LuoSvendsen")
Runtime type information.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
virtual void correct()
Correct diameter independent expressions.
Definition: LuoSvendsen.C:150
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Namespace for OpenFOAM.