Luo.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-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::diameterModels::coalescenceModels::Luo
28
29Description
30 Model of Luo (1993). The coalescence rate is calculated by
31
32 \f[
33 \frac{\pi}{4} (d_i + d_j)^2 u_{ij}
34 \mathrm{exp}
35 \left[
36 - C_1
37 \frac
38 {[0.75(1 + \xi_{ij}^2)(1 + \xi_{ij}^3)]^{1/2}}
39 {(\rho_d/\rho_c + C_{vm})^{1/2} (1 + \xi_{ij})^3}
40 \left(\frac{\rho_c d_i u_{ij}^2}{\sigma}\right)^{1/2}
41 \right]\,,
42 \f]
43
44 where
45
46 \f[
47 u_{ij} = \sqrt{\beta} (\epsilon_c d_i)^{1/3} \sqrt{1 + \xi_{ij}^{-2/3}}
48 \f]
49
50 is the mean approach velocity of the bubbles and
51
52 \f[
53 \xi_{ij} = d_i/d_j
54 \f]
55
56 their size ratio.
57
58 \vartable
59 d_i | Diameter of bubble i [m]
60 d_j | Diameter of bubble j [m]
61 u_{ij} | Mean approach velocity [m/s]
62 \xi_{ij} | Bubble size ratio [-]
63 \rho_d | Density of dispersed phase [kg/m3]
64 \rho_c | Density of continuous phase [kg/m3]
65 \sigma | Surface tension [N/m]
66 C_{vm} | Virtual mass coefficient [-]
67 C_1 | Coefficient [-]
68 \beta | Coefficient [-]
69 \epsilon_c | Continuous phase turbulent dissipation rate [m2/s3]
70 \endvartable
71
72 Reference:
73 \verbatim
74 Luo, H. (1993).
75 Coalescence, breakup and liquid circulation in bubble column reactors.
76 Dr. Ing (Doctoral dissertation, Thesis, Department of Chemical
77 Engineering, The Norwegian Institute of Technology, Trondheim, Norway).
78 \endverbatim
79
80Usage
81 \table
82 Property | Description | Required | Default value
83 beta | Coefficient beta | no | 2.0
84 C1 | Coefficient C1 | no | 1.0
85 \endtable
86
87SourceFiles
88 Luo.C
89
90\*---------------------------------------------------------------------------*/
91
92#ifndef Luo_H
93#define Luo_H
94
95#include "coalescenceModel.H"
96
97// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
98
99namespace Foam
100{
101namespace diameterModels
102{
103namespace coalescenceModels
104{
105
106/*---------------------------------------------------------------------------*\
107 Class Luo Declaration
108\*---------------------------------------------------------------------------*/
109
110class Luo
111:
112 public coalescenceModel
113{
114 // Private data
115
116 //- Coefficient beta, defaults to 2.0
117 dimensionedScalar beta_;
118
119 //- Optional coefficient C1, defaults to 1.0
121
122
123public:
124
125 //- Runtime type information
126 TypeName("Luo");
127
128 // Constructor
129
130 Luo
131 (
132 const populationBalanceModel& popBal,
133 const dictionary& dict
134 );
135
136
137 //- Destructor
138 virtual ~Luo() = default;
139
140
141 // Member Functions
142
143 //- Add to coalescenceRate
144 virtual void addToCoalescenceRate
145 (
146 volScalarField& coalescenceRate,
147 const label i,
148 const label j
149 );
150};
151
152
153// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154
155} // End namespace coalescenceModels
156} // End namespace diameterModels
157} // End namespace Foam
158
159// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160
161#endif
162
163// ************************************************************************* //
Base class for coalescence models.
Model of Luo (1993). The coalescence rate is calculated by.
Definition: Luo.H:171
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
Definition: Luo.C:76
Luo(const populationBalanceModel &popBal, const dictionary &dict)
Definition: Luo.C:61
TypeName("Luo")
Runtime type information.
Class that solves the univariate population balance equation by means of a class method (also called ...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73