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