All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
LuoSvendsen.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) 2019-2022 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 "LuoSvendsen.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace diameterModels
35 {
36 namespace binaryBreakupModels
37 {
38  defineTypeNameAndDebug(LuoSvendsen, 0);
40  (
41  binaryBreakupModel,
42  LuoSvendsen,
43  dictionary
44  );
45 }
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const populationBalanceModel& popBal,
55  const dictionary& dict
56 )
57 :
58  binaryBreakupModel(popBal, dict),
59  gammaUpperReg2by11_(),
60  gammaUpperReg5by11_(),
61  gammaUpperReg8by11_(),
62  C4_(dimensionedScalar::lookupOrDefault("C4", dict, dimless, 0.923)),
63  beta_(dimensionedScalar::lookupOrDefault("beta", dict, dimless, 2.05)),
64  minEddyRatio_
65  (
66  dimensionedScalar::lookupOrDefault("minEddyRatio", dict, dimless, 11.4)
67  ),
68  kolmogorovLengthScale_
69  (
70  IOobject
71  (
72  "kolmogorovLengthScale",
73  popBal_.time().timeName(),
74  popBal_.mesh()
75  ),
76  popBal_.mesh(),
78  (
79  "kolmogorovLengthScale",
80  dimLength,
81  Zero
82  )
83  )
84 {
85  List<Tuple2<scalar, scalar>> gammaUpperReg2by11Table;
86  List<Tuple2<scalar, scalar>> gammaUpperReg5by11Table;
87  List<Tuple2<scalar, scalar>> gammaUpperReg8by11Table;
88 
89  gammaUpperReg2by11Table.append(Tuple2<scalar, scalar>(0.0, 1.0));
90  gammaUpperReg5by11Table.append(Tuple2<scalar, scalar>(0.0, 1.0));
91  gammaUpperReg8by11Table.append(Tuple2<scalar, scalar>(0.0, 1.0));
92 
93  for (scalar z = 1e-2; z <= 10.0; z = z + 1e-2)
94  {
95  Tuple2<scalar, scalar> gamma2by11
96  (
97  z,
98  incGammaRatio_Q(2.0/11.0, z)
99  );
100 
101  Tuple2<scalar, scalar> gamma5by11
102  (
103  z,
104  incGammaRatio_Q(5.0/11.0, z)
105  );
106 
107  Tuple2<scalar, scalar> gamma8by11
108  (
109  z,
110  incGammaRatio_Q(8.0/11.0, z)
111  );
112 
113  gammaUpperReg2by11Table.append(gamma2by11);
114  gammaUpperReg5by11Table.append(gamma5by11);
115  gammaUpperReg8by11Table.append(gamma8by11);
116  }
117 
118  gammaUpperReg2by11_ =
119  new Function1s::Table<scalar>
120  (
121  "gamma2by11",
123  linearInterpolationWeights::typeName,
124  autoPtr<TableReader<scalar>>(nullptr),
125  gammaUpperReg2by11Table
126  );
127 
128  gammaUpperReg5by11_ =
129  new Function1s::Table<scalar>
130  (
131  "gamma5by11",
133  linearInterpolationWeights::typeName,
134  autoPtr<TableReader<scalar>>(nullptr),
135  gammaUpperReg5by11Table
136  );
137 
138  gammaUpperReg8by11_ =
139  new Function1s::Table<scalar>
140  (
141  "gamma8by11",
143  linearInterpolationWeights::typeName,
144  autoPtr<TableReader<scalar>>(nullptr),
145  gammaUpperReg8by11Table
146  );
147 }
148 
149 
150 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151 
153 {
155  pow025
156  (
157  pow3
158  (
159  popBal_.continuousPhase().thermo().nu()
160  )
161  /popBal_.continuousTurbulence().epsilon()
162  );
163 }
164 
165 
166 void
168 (
169  volScalarField& binaryBreakupRate,
170  const label i,
171  const label j
172 )
173 {
174  const phaseModel& continuousPhase = popBal_.continuousPhase();
175  const sizeGroup& fi = popBal_.sizeGroups()[i];
176  const sizeGroup& fj = popBal_.sizeGroups()[j];
177 
178  const dimensionedScalar cf
179  (
180  pow(fi.x()/fj.x(), 2.0/3.0) + pow((1 - fi.x()/fj.x()), 2.0/3.0) - 1
181  );
182 
183  const volScalarField b
184  (
185  12*cf*popBal_.sigmaWithContinuousPhase(fi.phase())
186  /(
187  beta_*continuousPhase.rho()*pow(fj.dSph(), 5.0/3.0)
188  *pow(popBal_.continuousTurbulence().epsilon(), 2.0/3.0)
189  )
190  );
191 
192  const volScalarField xiMin(minEddyRatio_*kolmogorovLengthScale_/fj.dSph());
193 
194  const volScalarField tMin(b/pow(xiMin, 11.0/3.0));
195 
196  volScalarField integral(3/(11*pow(b, 8.0/11.0)));
197 
198  forAll(integral, celli)
199  {
200  integral[celli] *=
201  2*pow(b[celli], 3.0/11.0)*tgamma(5.0/11.0)
202  *(
203  gammaUpperReg5by11_->value(b[celli])
204  - gammaUpperReg5by11_->value(tMin[celli])
205  );
206  }
207 
208  binaryBreakupRate +=
209  C4_*(1 - popBal_.alphas())/fj.x()
210  *cbrt
211  (
212  popBal_.continuousTurbulence().epsilon()
213  /sqr(fj.dSph())
214  )
215  *integral;
216 }
217 
218 
219 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
LuoSvendsen(const populationBalanceModel &popBal, const dictionary &dict)
virtual void precompute()
Precompute diameter independent expressions.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
dimensionedScalar pow025(const dimensionedScalar &ds)
const dimensionSet dimless
fvMesh & mesh
Macros for easy insertion into run-time selection tables.
const dimensionSet dimLength
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
virtual void addToBinaryBreakupRate(volScalarField &binaryBreakupRate, const label i, const label j)
Add to binary breakupRate.
dimensionedScalar cbrt(const dimensionedScalar &ds)
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
scalar incGammaRatio_Q(const scalar a, const scalar x)
Normalised upper incomplete gamma function.
Definition: incGamma.C:221
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
Namespace for OpenFOAM.
volScalarField kolmogorovLengthScale_
Kolmogorov length scale.
Definition: LiaoBase.H:70