Liao.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 "Liao.H"
28 #include "fvcGrad.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace diameterModels
35 {
36 namespace binaryBreakupModels
37 {
40  (
42  Liao,
44  );
45 }
46 }
47 }
48 
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
55 (
56  const populationBalanceModel& popBal,
57  const dictionary& dict
58 )
59 :
60  binaryBreakupModel(popBal, dict),
61  LiaoBase(popBal, dict),
62  BTurb_(dimensionedScalar::lookupOrDefault("BTurb", dict, dimless, 1)),
63  BShear_(dimensionedScalar::lookupOrDefault("BShear", dict, dimless, 1)),
64  BEddy_(dimensionedScalar::lookupOrDefault("BEddy", dict, dimless, 1)),
65  BFric_(dimensionedScalar::lookupOrDefault("BFric", dict, dimless, 0.25)),
66  turbulence_(dict.lookup("turbulence")),
67  laminarShear_(dict.lookup("laminarShear")),
68  turbulentShear_(dict.lookup("turbulentShear")),
69  interfacialFriction_(dict.lookup("interfacialFriction"))
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
76 {
78 }
79 
80 
83 (
84  volScalarField& binaryBreakupRate,
85  const label i,
86  const label j
87 )
88 {
89  const phaseModel& continuousPhase = popBal_.continuousPhase();
90  const sizeGroup& fi = popBal_.sizeGroups()[i];
91  const sizeGroup& fj = popBal_.sizeGroups()[j];
92 
93  dimensionedScalar dk(cbrt(pow3(fj.dSph()) - pow3(fi.dSph())));
94 
95  const volScalarField tauCrit1
96  (
97  6*popBal_.sigmaWithContinuousPhase(fj.phase())/fj.dSph()
98  *(sqr(fi.dSph()/fj.dSph()) + sqr(dk/fj.dSph()) - 1)
99  );
100 
101  const volScalarField tauCrit2
102  (
103  popBal_.sigmaWithContinuousPhase(fj.phase())/min(dk, fi.dSph())
104  );
105 
106  const volScalarField tauCrit(max(tauCrit1, tauCrit2));
107 
108  if (turbulence_)
109  {
110  const volScalarField tauTurb
111  (
112  pos(fj.dSph() - kolmogorovLengthScale_)*BTurb_*continuousPhase.rho()
113  *sqr(cbrt(popBal_.continuousTurbulence().epsilon()*fj.dSph()))
114  );
115 
116  binaryBreakupRate +=
117  pos(tauTurb - tauCrit)*1/fj.dSph()
118  *sqrt(mag(tauTurb - tauCrit)/continuousPhase.rho())/fj.x();
119  }
120 
121  if (laminarShear_)
122  {
123  const volScalarField tauShear
124  (
125  BShear_*continuousPhase.thermo().mu()*shearStrainRate_
126  );
127 
128  binaryBreakupRate +=
129  pos(tauShear - tauCrit)*1/fj.dSph()
130  *sqrt(mag(tauShear - tauCrit)/continuousPhase.rho())/fj.x();
131  }
132 
133  if (turbulentShear_)
134  {
135  const volScalarField tauEddy
136  (
137  pos0(kolmogorovLengthScale_ - fj.dSph())
138  *BEddy_*continuousPhase.thermo().mu()*eddyStrainRate_
139  );
140 
141  binaryBreakupRate +=
142  pos(tauEddy - tauCrit)*1/fj.dSph()
143  *sqrt(mag(tauEddy - tauCrit)/continuousPhase.rho())/fj.x();
144  }
145 
146  if (interfacialFriction_)
147  {
148  const volScalarField tauFric
149  (
150  BFric_*0.5*continuousPhase.rho()*sqr(uTerminal_[j])*Cd_[j]
151  );
152 
153  binaryBreakupRate +=
154  pos(tauFric - tauCrit)*1/fj.dSph()
155  *sqrt(mag(tauFric - tauCrit)/continuousPhase.rho())/fj.x();
156  }
157 }
158 
159 
160 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
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:91
Base class for binary breakup models that provide a breakup rate between a size class pair directly,...
Bubble breakup model of Liao et al. (2015). The terminal velocities and drag coefficients are compute...
Definition: Liao.H:131
virtual void precompute()
Precompute diameter independent expressions.
Definition: Liao.C:75
virtual void addToBinaryBreakupRate(volScalarField &binaryBreakupRate, const label i, const label j)
Add to binary breakupRate.
Definition: Liao.C:83
Liao(const populationBalanceModel &popBal, const dictionary &dict)
Definition: Liao.C:55
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:50
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition: sizeGroupI.H:57
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:36
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
virtual tmp< volScalarField > mu() const =0
Dynamic viscosity of mixture [kg/m/s].
virtual const rhoThermo & thermo() const =0
Return the thermophysical model.
virtual const volScalarField & rho() const =0
Return the density field.
Calculate the gradient of the given field.
addToRunTimeSelectionTable(binaryBreakupModel, LehrMilliesMewes, dictionary)
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedScalar pos0(const dimensionedScalar &ds)
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimless
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 cbrt(const dimensionedScalar &ds)
dictionary dict