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"
27 #include "fvcGrad.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace diameterModels
36 {
37 namespace binaryBreakupModels
38 {
41  (
43  Liao,
45  );
46 }
47 }
48 }
49 
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
56 (
57  const populationBalanceModel& popBal,
58  const dictionary& dict
59 )
60 :
61  binaryBreakupModel(popBal, dict),
62  LiaoBase(popBal, dict),
63  BTurb_(dimensionedScalar::lookupOrDefault("BTurb", dict, dimless, 1)),
64  BShear_(dimensionedScalar::lookupOrDefault("BShear", dict, dimless, 1)),
65  BEddy_(dimensionedScalar::lookupOrDefault("BEddy", dict, dimless, 1)),
66  BFric_(dimensionedScalar::lookupOrDefault("BFric", dict, dimless, 0.25)),
67  turbulence_(dict.lookup("turbulence")),
68  laminarShear_(dict.lookup("laminarShear")),
69  turbulentShear_(dict.lookup("turbulentShear")),
70  interfacialFriction_(dict.lookup("interfacialFriction"))
71 {}
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
77 {
79 }
80 
81 
84 (
85  volScalarField& binaryBreakupRate,
86  const label i,
87  const label j
88 )
89 {
90  const phaseModel& continuousPhase = popBal_.continuousPhase();
91  const sizeGroup& fi = popBal_.sizeGroups()[i];
92  const sizeGroup& fj = popBal_.sizeGroups()[j];
93 
94  dimensionedScalar dk(cbrt(pow3(fj.dSph()) - pow3(fi.dSph())));
95 
96  const volScalarField tauCrit1
97  (
98  6*popBal_.sigmaWithContinuousPhase(fj.phase())/fj.dSph()
99  *(sqr(fi.dSph()/fj.dSph()) + sqr(dk/fj.dSph()) - 1)
100  );
101 
102  const volScalarField tauCrit2
103  (
104  popBal_.sigmaWithContinuousPhase(fj.phase())/min(dk, fi.dSph())
105  );
106 
107  const volScalarField tauCrit(max(tauCrit1, tauCrit2));
108 
109  if (turbulence_)
110  {
111  const volScalarField tauTurb
112  (
113  pos(fj.dSph() - kolmogorovLengthScale_)*BTurb_*continuousPhase.rho()
114  *sqr(cbrt(popBal_.continuousTurbulence().epsilon()*fj.dSph()))
115  );
116 
117  binaryBreakupRate +=
118  pos(tauTurb - tauCrit)*1/fj.dSph()
119  *sqrt(mag(tauTurb - tauCrit)/continuousPhase.rho())/fj.x();
120  }
121 
122  if (laminarShear_)
123  {
124  const volScalarField tauShear
125  (
126  BShear_*continuousPhase.fluidThermo().mu()*shearStrainRate_
127  );
128 
129  binaryBreakupRate +=
130  pos(tauShear - tauCrit)*1/fj.dSph()
131  *sqrt(mag(tauShear - tauCrit)/continuousPhase.rho())/fj.x();
132  }
133 
134  if (turbulentShear_)
135  {
136  const volScalarField tauEddy
137  (
138  pos0(kolmogorovLengthScale_ - fj.dSph())
139  *BEddy_*continuousPhase.fluidThermo().mu()*eddyStrainRate_
140  );
141 
142  binaryBreakupRate +=
143  pos(tauEddy - tauCrit)*1/fj.dSph()
144  *sqrt(mag(tauEddy - tauCrit)/continuousPhase.rho())/fj.x();
145  }
146 
147  if (interfacialFriction_)
148  {
149  const volScalarField tauFric
150  (
151  BFric_*0.5*continuousPhase.rho()*sqr(uTerminal_[j])*Cd_[j]
152  );
153 
154  binaryBreakupRate +=
155  pos(tauFric - tauCrit)*1/fj.dSph()
156  *sqrt(mag(tauFric - tauCrit)/continuousPhase.rho())/fj.x();
157  }
158 }
159 
160 
161 // ************************************************************************* //
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:92
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:76
virtual void addToBinaryBreakupRate(volScalarField &binaryBreakupRate, const label i, const label j)
Add to binary breakupRate.
Definition: Liao.C:84
Liao(const populationBalanceModel &popBal, const dictionary &dict)
Definition: Liao.C:56
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 dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition: sizeGroupI.H:58
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(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