LuoSvendsen.C
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  Copyright (C) 2020-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "LuoSvendsen.H"
31 #include "phaseCompressibleTurbulenceModel.H"
32 #include "tableBounds.H"
33 #include "MathFunctions.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace diameterModels
40 {
41 namespace binaryBreakupModels
42 {
43  defineTypeNameAndDebug(LuoSvendsen, 0);
45  (
46  binaryBreakupModel,
47  LuoSvendsen,
48  dictionary
49  );
50 }
51 }
52 }
53 
54 
55 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
56 
58 (
59  const populationBalanceModel& popBal,
60  const dictionary& dict
61 )
62 :
63  binaryBreakupModel(popBal, dict),
64  gammaUpperReg2by11_(),
65  gammaUpperReg5by11_(),
66  gammaUpperReg8by11_(),
67  C4_(dimensionedScalar::getOrDefault("C4", dict, dimless, 0.923)),
68  beta_(dimensionedScalar::getOrDefault("beta", dict, dimless, 2.05)),
69  minEddyRatio_
70  (
71  dimensionedScalar::getOrDefault("minEddyRatio", dict, dimless, 11.4)
72  ),
73  kolmogorovLengthScale_
74  (
75  IOobject
76  (
77  "kolmogorovLengthScale",
78  popBal_.time().timeName(),
79  popBal_.mesh()
80  ),
81  popBal_.mesh(),
83  (
84  "kolmogorovLengthScale",
85  dimLength,
86  Zero
87  )
88  )
89 {
90  List<Tuple2<scalar, scalar>> gammaUpperReg2by11Table;
91  List<Tuple2<scalar, scalar>> gammaUpperReg5by11Table;
92  List<Tuple2<scalar, scalar>> gammaUpperReg8by11Table;
93 
94  gammaUpperReg2by11Table.append(Tuple2<scalar, scalar>(0.0, 1.0));
95  gammaUpperReg5by11Table.append(Tuple2<scalar, scalar>(0.0, 1.0));
96  gammaUpperReg8by11Table.append(Tuple2<scalar, scalar>(0.0, 1.0));
97 
98  for (scalar z = 1e-2; z <= 10.0; z = z + 1e-2)
99  {
100  Tuple2<scalar, scalar> gamma2by11
101  (
102  z,
103  Math::incGammaRatio_Q(2.0/11.0, z)
104  );
105 
106  Tuple2<scalar, scalar> gamma5by11
107  (
108  z,
109  Math::incGammaRatio_Q(5.0/11.0, z)
110  );
111 
112  Tuple2<scalar, scalar> gamma8by11
113  (
114  z,
115  Math::incGammaRatio_Q(8.0/11.0, z)
116  );
117 
118  gammaUpperReg2by11Table.append(gamma2by11);
119  gammaUpperReg5by11Table.append(gamma5by11);
120  gammaUpperReg8by11Table.append(gamma8by11);
121  }
122 
123  gammaUpperReg2by11_.reset
124  (
126  (
127  gammaUpperReg2by11Table,
128  bounds::repeatableBounding::CLAMP,
129  "gamma2by11"
130  )
131  );
132 
133  gammaUpperReg5by11_.reset
134  (
136  (
137  gammaUpperReg5by11Table,
138  bounds::repeatableBounding::CLAMP,
139  "gamma5by11"
140  )
141  );
142 
143  gammaUpperReg8by11_.reset
144  (
146  (
147  gammaUpperReg8by11Table,
148  bounds::repeatableBounding::CLAMP,
149  "gamma8by11"
150  )
151  );
152 }
153 
154 
155 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
156 
158 {
159  kolmogorovLengthScale_ =
160  pow025
161  (
162  pow3
163  (
165  )
166  /popBal_.continuousTurbulence().epsilon()
167  );
168 }
169 
170 
171 void
173 (
174  volScalarField& binaryBreakupRate,
175  const label i,
176  const label j
177 )
178 {
179  const phaseModel& continuousPhase = popBal_.continuousPhase();
180  const sizeGroup& fi = popBal_.sizeGroups()[i];
181  const sizeGroup& fj = popBal_.sizeGroups()[j];
182 
183  const dimensionedScalar cf
184  (
185  pow(fi.x()/fj.x(), 2.0/3.0) + pow((1 - fi.x()/fj.x()), 2.0/3.0) - 1
186  );
187 
188  const volScalarField b
189  (
190  12.0*cf*popBal_.sigmaWithContinuousPhase(fi.phase())
191  /(
192  beta_*continuousPhase.rho()*pow(fj.d(), 5.0/3.0)
193  *pow(popBal_.continuousTurbulence().epsilon(), 2.0/3.0)
194  )
195  );
196 
197  const volScalarField xiMin(minEddyRatio_*kolmogorovLengthScale_/fj.d());
198 
199  const volScalarField tMin(b/pow(xiMin, 11.0/3.0));
200 
201  volScalarField integral(3.0/(11.0*pow(b, 8.0/11.0)));
202 
203  forAll(integral, celli)
204  {
205  integral[celli] *=
206  2.0*pow(b[celli], 3.0/11.0)*tgamma(5.0/11.0)
207  *(
208  gammaUpperReg5by11_()(b[celli])
209  - gammaUpperReg5by11_()(tMin[celli])
210  );
211  }
212 
213  binaryBreakupRate +=
214  C4_*(1 - popBal_.alphas())/fj.x()
215  *cbrt
216  (
217  popBal_.continuousTurbulence().epsilon()
218  /sqr(fj.d())
219  )
220  *integral;
221 }
222 
223 
224 // ************************************************************************* //
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:54
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::diameterModels::sizeGroup::phase
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:38
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::diameterModels::binaryBreakupModels::LuoSvendsen::addToBinaryBreakupRate
virtual void addToBinaryBreakupRate(volScalarField &binaryBreakupRate, const label i, const label j)
Add to binary breakupRate.
Definition: LuoSvendsen.C:173
Foam::diameterModels::binaryBreakupModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(binaryBreakupModel, LehrMilliesMewes, dictionary)
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:133
Foam::Math::incGammaRatio_Q
scalar incGammaRatio_Q(const scalar a, const scalar x)
Regularised upper incomplete gamma function.
Definition: incGamma.C:199
LuoSvendsen.H
MathFunctions.H
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::phaseModel::nu
const dimensionedScalar & nu() const
Definition: phaseModel.H:152
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::diameterModels::sizeGroup
This class represents a single sizeGroup belonging to a velocityGroup. The main property of a sizeGro...
Definition: sizeGroup.H:96
Foam::diameterModels::binaryBreakupModels::LuoSvendsen::LuoSvendsen
LuoSvendsen(const populationBalanceModel &popBal, const dictionary &dict)
Definition: LuoSvendsen.C:58
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::diameterModels::populationBalanceModel
Class that solves the univariate population balance equation by means of a class method (also called ...
Definition: populationBalanceModel.H:179
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::diameterModels::populationBalanceModel::continuousPhase
const phaseModel & continuousPhase() const
Return continuous phase.
Definition: populationBalanceModelI.H:70
Foam::diameterModels::binaryBreakupModels::defineTypeNameAndDebug
defineTypeNameAndDebug(LehrMilliesMewes, 0)
Foam::diameterModels::sizeGroup::d
const dimensionedScalar & d() const
Return representative diameter of the sizeGroup.
Definition: sizeGroupI.H:52
tableBounds.H
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::diameterModels::binaryBreakupModel::popBal_
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
Definition: binaryBreakupModel.H:61
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::diameterModels::binaryBreakupModel
Base class for binary breakup models which give the breakup rate between a sizeGroup pair directly,...
Definition: binaryBreakupModel.H:54
Foam::diameterModels::populationBalanceModel::continuousTurbulence
const phaseCompressibleTurbulenceModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
Definition: populationBalanceModel.C:1220
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::interpolationTable< scalar >
Foam::phaseModel::rho
const dimensionedScalar & rho() const
Definition: phaseModel.H:167
Foam::Tuple2< scalar, scalar >
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:155
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::diameterModels::sizeGroup::x
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition: sizeGroupI.H:59
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::diameterModels::binaryBreakupModels::LuoSvendsen::correct
virtual void correct()
Correct diameter independent expressions.
Definition: LuoSvendsen.C:157