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-2020 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"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace diameterModels
36 {
37 namespace binaryBreakupModels
38 {
39  defineTypeNameAndDebug(LuoSvendsen, 0);
41  (
42  binaryBreakupModel,
43  LuoSvendsen,
44  dictionary
45  );
46 }
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const populationBalanceModel& popBal,
56  const dictionary& dict
57 )
58 :
59  binaryBreakupModel(popBal, dict),
60  gammaUpperReg2by11_(),
61  gammaUpperReg5by11_(),
62  gammaUpperReg8by11_(),
63  C4_(dimensionedScalar::lookupOrDefault("C4", dict, dimless, 0.923)),
64  beta_(dimensionedScalar::lookupOrDefault("beta", dict, dimless, 2.05)),
65  minEddyRatio_
66  (
67  dimensionedScalar::lookupOrDefault("minEddyRatio", dict, dimless, 11.4)
68  ),
69  kolmogorovLengthScale_
70  (
71  IOobject
72  (
73  "kolmogorovLengthScale",
74  popBal_.time().timeName(),
75  popBal_.mesh()
76  ),
77  popBal_.mesh(),
79  (
80  "kolmogorovLengthScale",
81  dimLength,
82  Zero
83  )
84  )
85 {
86  List<Tuple2<scalar, scalar>> gammaUpperReg2by11Table;
87  List<Tuple2<scalar, scalar>> gammaUpperReg5by11Table;
88  List<Tuple2<scalar, scalar>> gammaUpperReg8by11Table;
89 
90  gammaUpperReg2by11Table.append(Tuple2<scalar, scalar>(0.0, 1.0));
91  gammaUpperReg5by11Table.append(Tuple2<scalar, scalar>(0.0, 1.0));
92  gammaUpperReg8by11Table.append(Tuple2<scalar, scalar>(0.0, 1.0));
93 
94  for (scalar z = 1e-2; z <= 10.0; z = z + 1e-2)
95  {
96  Tuple2<scalar, scalar> gamma2by11
97  (
98  z,
99  incGammaRatio_Q(2.0/11.0, z)
100  );
101 
102  Tuple2<scalar, scalar> gamma5by11
103  (
104  z,
105  incGammaRatio_Q(5.0/11.0, z)
106  );
107 
108  Tuple2<scalar, scalar> gamma8by11
109  (
110  z,
111  incGammaRatio_Q(8.0/11.0, z)
112  );
113 
114  gammaUpperReg2by11Table.append(gamma2by11);
115  gammaUpperReg5by11Table.append(gamma5by11);
116  gammaUpperReg8by11Table.append(gamma8by11);
117  }
118 
119  gammaUpperReg2by11_ =
120  new Function1s::Table<scalar>
121  (
122  "gamma2by11",
124  linearInterpolationWeights::typeName,
125  gammaUpperReg2by11Table
126  );
127 
128  gammaUpperReg5by11_ =
129  new Function1s::Table<scalar>
130  (
131  "gamma5by11",
133  linearInterpolationWeights::typeName,
134  gammaUpperReg5by11Table
135  );
136 
137  gammaUpperReg8by11_ =
138  new Function1s::Table<scalar>
139  (
140  "gamma8by11",
142  linearInterpolationWeights::typeName,
143  gammaUpperReg8by11Table
144  );
145 }
146 
147 
148 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
149 
151 {
152  kolmogorovLengthScale_ =
153  pow025
154  (
155  pow3
156  (
158  )
160  );
161 }
162 
163 
164 void
166 (
167  volScalarField& binaryBreakupRate,
168  const label i,
169  const label j
170 )
171 {
172  const phaseModel& continuousPhase = popBal_.continuousPhase();
173  const sizeGroup& fi = popBal_.sizeGroups()[i];
174  const sizeGroup& fj = popBal_.sizeGroups()[j];
175 
176  const dimensionedScalar cf
177  (
178  pow(fi.x()/fj.x(), 2.0/3.0) + pow((1 - fi.x()/fj.x()), 2.0/3.0) - 1
179  );
180 
181  const volScalarField b
182  (
183  12.0*cf*popBal_.sigmaWithContinuousPhase(fi.phase())
184  /(
185  beta_*continuousPhase.rho()*pow(fj.dSph(), 5.0/3.0)
186  *pow(popBal_.continuousTurbulence().epsilon(), 2.0/3.0)
187  )
188  );
189 
190  const volScalarField xiMin(minEddyRatio_*kolmogorovLengthScale_/fj.dSph());
191 
192  const volScalarField tMin(b/pow(xiMin, 11.0/3.0));
193 
194  volScalarField integral(3.0/(11.0*pow(b, 8.0/11.0)));
195 
196  forAll(integral, celli)
197  {
198  integral[celli] *=
199  2.0*pow(b[celli], 3.0/11.0)*tgamma(5.0/11.0)
200  *(
201  gammaUpperReg5by11_->value(b[celli])
202  - gammaUpperReg5by11_->value(tMin[celli])
203  );
204  }
205 
206  binaryBreakupRate +=
207  C4_*(1 - popBal_.alphas())/fj.x()
208  *cbrt
209  (
211  /sqr(fj.dSph())
212  )
213  *integral;
214 }
215 
216 
217 // ************************************************************************* //
#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:57
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
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
Definition: phaseModel.H:120
dimensionedScalar cbrt(const dimensionedScalar &ds)
const phaseCompressibleMomentumTransportModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
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 populationBalanceModel & popBal_
Reference to the populationBalanceModel.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
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.
const volScalarField & alphas() const
Return total void of phases belonging to this populationBalance.
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
Definition: fluidThermo.C:91
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:105
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
Namespace for OpenFOAM.