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 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"
28 #include "phaseCompressibleTurbulenceModel.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 interpolationTable<scalar>
120  (
121  gammaUpperReg2by11Table,
123  "gamma2by11"
124  );
125 
126  gammaUpperReg5by11_ =
127  new interpolationTable<scalar>
128  (
129  gammaUpperReg5by11Table,
130  interpolationTable<scalar>::CLAMP,
131  "gamma5by11"
132  );
133 
134  gammaUpperReg8by11_ =
135  new interpolationTable<scalar>
136  (
137  gammaUpperReg8by11Table,
138  interpolationTable<scalar>::CLAMP,
139  "gamma8by11"
140  );
141 }
142 
143 
144 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
145 
147 {
148  kolmogorovLengthScale_ =
149  pow025
150  (
151  pow3
152  (
154  )
156  );
157 }
158 
159 
160 void
162 (
163  volScalarField& binaryBreakupRate,
164  const label i,
165  const label j
166 )
167 {
168  const phaseModel& continuousPhase = popBal_.continuousPhase();
169  const sizeGroup& fi = popBal_.sizeGroups()[i];
170  const sizeGroup& fj = popBal_.sizeGroups()[j];
171 
172  const dimensionedScalar cf
173  (
174  pow(fi.x()/fj.x(), 2.0/3.0) + pow((1 - fi.x()/fj.x()), 2.0/3.0) - 1
175  );
176 
177  const volScalarField b
178  (
179  12.0*cf*popBal_.sigmaWithContinuousPhase(fi.phase())
180  /(
181  beta_*continuousPhase.rho()*pow(fj.d(), 5.0/3.0)
182  *pow(popBal_.continuousTurbulence().epsilon(), 2.0/3.0)
183  )
184  );
185 
186  const volScalarField xiMin(minEddyRatio_*kolmogorovLengthScale_/fj.d());
187 
188  const volScalarField tMin(b/pow(xiMin, 11.0/3.0));
189 
190  volScalarField integral(3.0/(11.0*pow(b, 8.0/11.0)));
191 
192  forAll(integral, celli)
193  {
194  integral[celli] *=
195  2.0*pow(b[celli], 3.0/11.0)*tgamma(5.0/11.0)
196  *(
197  gammaUpperReg5by11_()(b[celli])
198  - gammaUpperReg5by11_()(tMin[celli])
199  );
200  }
201 
202  binaryBreakupRate +=
203  C4_*(1 - popBal_.alphas())/fj.x()
204  *cbrt
205  (
207  /sqr(fj.d())
208  )
209  *integral;
210 }
211 
212 
213 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
LuoSvendsen(const populationBalanceModel &popBal, const dictionary &dict)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual void correct()
Correct diameter independent expressions.
dimensionedScalar pow025(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual void addToBinaryBreakupRate(volScalarField &binaryBreakupRate, const label i, const label j)
Add to binary breakupRate.
dynamicFvMesh & mesh
const phaseModel & continuousPhase() const
Return continuous phase.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
dimensionedScalar cbrt(const dimensionedScalar &ds)
Clamp value to the start/end value.
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
scalar incGammaRatio_Q(const scalar a, const scalar x)
Normalized upper incomplete gamma function.
Definition: incGamma.C:221
const dimensionedScalar & nu() const
Return the laminar viscosity.
Definition: phaseModel.H:150
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const phaseCompressibleTurbulenceModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
defineTypeNameAndDebug(combustionModel, 0)
const UPtrList< sizeGroup > & sizeGroups() const
Return the sizeGroups belonging to this populationBalance.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
const volScalarField & alphas() const
Return total void of phases belonging to this populationBalance.
const tmp< volScalarField > sigmaWithContinuousPhase(const phaseModel &dispersedPhase) const
Return the surface tension coefficient between a given dispersed.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
Namespace for OpenFOAM.