LiaoCoalescence.C
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) 2021-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 \*---------------------------------------------------------------------------*/
25 
26 #include "LiaoCoalescence.H"
27 #include "phaseSystem.H"
29 #include "fvcGrad.H"
30 #include "dragModel.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace diameterModels
38 {
39 namespace coalescenceModels
40 {
43  (
47  );
48 }
49 }
50 }
51 
53 
54 
55 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
56 
59 (
60  const populationBalanceModel& popBal,
61  const dictionary& dict
62 )
63 :
64  coalescenceModel(popBal, dict),
65  LiaoBase(popBal, dict),
66  PMax_(dimensionedScalar::lookupOrDefault("PMax", dict, dimless, 0.8)),
67  AH_(dimensionedScalar::lookupOrDefault("AH", dict, dimEnergy, 3.7e-20)),
68  CEff_(dimensionedScalar::lookupOrDefault("CEff", dict, dimless, 2.5)),
69  CTurb_(dimensionedScalar::lookupOrDefault("CTurb", dict, dimless, 1)),
70  CBuoy_(dimensionedScalar::lookupOrDefault("CBuoy", dict, dimless, 1)),
71  CShear_(dimensionedScalar::lookupOrDefault("CShear", dict, dimless, 1)),
72  CEddy_(dimensionedScalar::lookupOrDefault("CEddy", dict, dimless, 1)),
73  CWake_(dimensionedScalar::lookupOrDefault("CWake", dict, dimless, 1)),
74  turbulence_(dict.lookup("turbulence")),
75  buoyancy_(dict.lookup("buoyancy")),
76  laminarShear_(dict.lookup("laminarShear")),
77  eddyCapture_(dict.lookup("eddyCapture")),
78  wakeEntrainment_(dict.lookup("wakeEntrainment")),
79  CPack_
80  (
81  IOobject
82  (
83  "CPack",
84  popBal_.time().name(),
85  popBal_.mesh()
86  ),
87  popBal_.mesh(),
89  (
90  "CPack",
91  dimless,
92  Zero
93  )
94  ),
95  CPackMax_
96  (
97  dimensionedScalar::lookupOrDefault("CPackMax", dict, dimless, 1e5)
98  ),
99  dCrit_
100  (
101  IOobject
102  (
103  "dCrit",
104  popBal_.time().name(),
105  popBal_.mesh()
106  ),
107  popBal_.mesh(),
109  (
110  "dCrit",
111  dimLength,
112  Zero
113  )
114  ),
115  uRelTurb_
116  (
117  IOobject
118  (
119  "uRelTurb",
120  popBal_.time().name(),
121  popBal_.mesh()
122  ),
123  popBal_.mesh(),
125  (
126  "uRelTurb",
127  dimVelocity,
128  Zero
129  )
130  ),
131  uRelBuoy_
132  (
133  IOobject
134  (
135  "uRelBuoy",
136  popBal_.time().name(),
137  popBal_.mesh()
138  ),
139  popBal_.mesh(),
141  (
142  "uRelBuoy",
143  dimVelocity,
144  Zero
145  )
146  ),
147  uRelShear_
148  (
149  IOobject
150  (
151  "uRelShear",
152  popBal_.time().name(),
153  popBal_.mesh()
154  ),
155  popBal_.mesh(),
157  (
158  "uRelShear",
159  dimVelocity,
160  Zero
161  )
162  )
163 {}
164 
165 
166 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
167 
169 {
171 
172  CPack_ = min(PMax_/max(PMax_ - popBal_.alphas(), SMALL), CPackMax_);
173 
175  popBal_.mesh().lookupObject<uniformDimensionedVectorField>("g");
176 
177  dCrit_ =
178  4*sqrt
179  (
180  popBal_.sigmaWithContinuousPhase(popBal_.sizeGroups()[1].phase())()
181  /(mag(g)*(popBal_.continuousPhase().rho()
182  - popBal_.sizeGroups()[1].phase().rho()))
183  );
184 }
185 
186 
189 (
190  volScalarField& coalescenceRate,
191  const label i,
192  const label j
193 )
194 {
195  const phaseModel& continuousPhase = popBal_.continuousPhase();
196  const sizeGroup& fi = popBal_.sizeGroups()[i];
197  const sizeGroup& fj = popBal_.sizeGroups()[j];
198 
199  dimensionedScalar dEq(2*fi.dSph()*fj.dSph()/(fi.dSph() + fj.dSph()));
200  dimensionedScalar Aij(pi*0.25*sqr(fi.dSph() + fj.dSph()));
201 
202  if (turbulence_)
203  {
204  uRelTurb_ =
205  CTurb_*sqrt(2.0)
206  *sqrt(sqr(cbrt(fi.dSph())) + sqr(cbrt(fj.dSph())))
207  *cbrt(popBal_.continuousTurbulence().epsilon());
208  }
209 
210  if (buoyancy_)
211  {
212  uRelBuoy_ = CBuoy_*mag(uTerminal_[i] - uTerminal_[j]);
213  }
214 
215  if (laminarShear_)
216  {
217  uRelShear_ = CShear_*0.5/pi*(fi.dSph() + fj.dSph())*shearStrainRate_;
218  }
219 
220  const volScalarField collisionEfficiency
221  (
222  neg(kolmogorovLengthScale_ - (fi.dSph() + fj.dSph()))
223  *exp
224  (
225  - CEff_*sqrt
226  (
227  continuousPhase.rho()*dEq
228  /popBal_.sigmaWithContinuousPhase(fi.phase())
229  *sqr(max(uRelTurb_, max(uRelBuoy_, uRelShear_)))
230  )
231  )
232  + pos0(kolmogorovLengthScale_ - (fi.dSph() + fj.dSph()))
233  *exp
234  (
235  - 3*continuousPhase.fluidThermo().mu()*dEq*eddyStrainRate_
236  /(4*popBal_.sigmaWithContinuousPhase(fi.phase()))
237  *log
238  (
239  cbrt
240  (
241  pi*popBal_.sigmaWithContinuousPhase(fi.phase())*sqr(dEq)
242  /(32*AH_)
243  )
244  )
245  )
246  );
247 
248  if (turbulence_)
249  {
250  coalescenceRate +=
251  neg(kolmogorovLengthScale_ - (fi.dSph() + fj.dSph()))
252  *CPack_*Aij*uRelTurb_*collisionEfficiency;
253  }
254 
255  if (buoyancy_)
256  {
257  coalescenceRate += CPack_*0.5*Aij*uRelBuoy_*collisionEfficiency;
258  }
259 
260  if (laminarShear_)
261  {
262  coalescenceRate += CPack_*0.5*Aij*uRelShear_*collisionEfficiency;
263  }
264 
265  if (eddyCapture_)
266  {
267  const volScalarField uRelEddy
268  (
269  CEddy_*0.5/pi*(fi.dSph() + fj.dSph())*eddyStrainRate_
270  );
271 
272  coalescenceRate +=
273  pos0(kolmogorovLengthScale_ - (fi.dSph() + fj.dSph()))
274  *CPack_*0.5*Aij*uRelEddy*collisionEfficiency;
275  }
276 
277  if (wakeEntrainment_)
278  {
279  const dimensionedScalar uRelWakeI
280  (
281  CWake_*uTerminal_[i]*cbrt(Cd_[i])
282  );
283 
284  const dimensionedScalar uRelWakeJ
285  (
286  CWake_*uTerminal_[j]*cbrt(Cd_[j])
287  );
288 
289  coalescenceRate +=
290  CPack_*0.125*pi
291  *(
292  sqr(fi.dSph())*uRelWakeI
293  *pos0(fi.dSph() - 0.5*dCrit_)
294  *(
295  pow6(fi.dSph() - 0.5*dCrit_)
296  /(pow6(fi.dSph() - 0.5*dCrit_) + pow6(0.5*dCrit_))
297  )
298  + sqr(fj.dSph())*uRelWakeJ
299  *pos0(fj.dSph() - 0.5*dCrit_)
300  *(
301  pow6(fj.dSph() - 0.5*dCrit_)
302  /(pow6(fj.dSph() - 0.5*dCrit_) + pow6(0.5*dCrit_))
303  )
304  );
305  }
306 }
307 
308 
309 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Base class for coalescence and breakup models of Liao et al. (2015).
Definition: LiaoBase.H:61
virtual void precompute()
Precompute diameter independent expressions.
Definition: LiaoBase.C:92
Base class for coalescence models.
Bubble coalescence model of Liao et al. (2015). The terminal velocities and drag coefficients are com...
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
virtual void precompute()
Precompute diameter independent expressions.
LiaoCoalescence(const populationBalanceModel &popBal, const dictionary &dict)
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:102
const dimensionedScalar & dSph() const
Return representative spherical diameter of the sizeGroup.
Definition: sizeGroupI.H:51
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:37
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
virtual const volScalarField & mu() const =0
Dynamic viscosity of mixture [kg/m/s].
virtual const rhoFluidThermo & fluidThermo() const =0
Return the thermophysical model.
virtual const volScalarField & rho() const =0
Return the density field.
Calculate the gradient of the given field.
addToRunTimeSelectionTable(coalescenceModel, AdachiStuartFokkink, dictionary)
defineTypeNameAndDebug(AdachiStuartFokkink, 0)
Namespace for OpenFOAM.
dimensionedScalar pow6(const dimensionedScalar &ds)
static const zero Zero
Definition: zero.H:97
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar pos0(const dimensionedScalar &ds)
static const scalar SMALL
Definition: scalar.H:115
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
const dimensionSet dimEnergy
dimensionedSymmTensor sqr(const dimensionedVector &dv)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
const dimensionSet dimLength
dimensionedScalar log(const dimensionedScalar &ds)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar neg(const dimensionedScalar &ds)
const dimensionSet dimVelocity
dimensionedScalar cbrt(const dimensionedScalar &ds)
dictionary dict