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 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 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace diameterModels
39 {
40 namespace binaryBreakupModels
41 {
42  defineTypeNameAndDebug(LuoSvendsen, 0);
44  (
45  binaryBreakupModel,
46  LuoSvendsen,
47  dictionary
48  );
49 }
50 }
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
57 (
58  const populationBalanceModel& popBal,
59  const dictionary& dict
60 )
61 :
62  binaryBreakupModel(popBal, dict),
63  gammaUpperReg2by11_(),
64  gammaUpperReg5by11_(),
65  gammaUpperReg8by11_(),
66  C4_(dimensionedScalar::getOrDefault("C4", dict, dimless, 0.923)),
67  beta_(dimensionedScalar::getOrDefault("beta", dict, dimless, 2.05)),
68  minEddyRatio_
69  (
70  dimensionedScalar::getOrDefault("minEddyRatio", dict, dimless, 11.4)
71  ),
72  kolmogorovLengthScale_
73  (
74  IOobject
75  (
76  "kolmogorovLengthScale",
77  popBal_.time().timeName(),
78  popBal_.mesh()
79  ),
80  popBal_.mesh(),
82  (
83  "kolmogorovLengthScale",
84  dimLength,
85  Zero
86  )
87  )
88 {
89  List<Tuple2<scalar, scalar>> gammaUpperReg2by11Table;
90  List<Tuple2<scalar, scalar>> gammaUpperReg5by11Table;
91  List<Tuple2<scalar, scalar>> gammaUpperReg8by11Table;
92 
93  gammaUpperReg2by11Table.append(Tuple2<scalar, scalar>(0.0, 1.0));
94  gammaUpperReg5by11Table.append(Tuple2<scalar, scalar>(0.0, 1.0));
95  gammaUpperReg8by11Table.append(Tuple2<scalar, scalar>(0.0, 1.0));
96 
97  for (scalar z = 1e-2; z <= 10.0; z = z + 1e-2)
98  {
99  Tuple2<scalar, scalar> gamma2by11
100  (
101  z,
102  incGammaRatio_Q(2.0/11.0, z)
103  );
104 
105  Tuple2<scalar, scalar> gamma5by11
106  (
107  z,
108  incGammaRatio_Q(5.0/11.0, z)
109  );
110 
111  Tuple2<scalar, scalar> gamma8by11
112  (
113  z,
114  incGammaRatio_Q(8.0/11.0, z)
115  );
116 
117  gammaUpperReg2by11Table.append(gamma2by11);
118  gammaUpperReg5by11Table.append(gamma5by11);
119  gammaUpperReg8by11Table.append(gamma8by11);
120  }
121 
122  gammaUpperReg2by11_.reset
123  (
125  (
126  gammaUpperReg2by11Table,
127  bounds::repeatableBounding::CLAMP,
128  "gamma2by11"
129  )
130  );
131 
132  gammaUpperReg5by11_.reset
133  (
135  (
136  gammaUpperReg5by11Table,
137  bounds::repeatableBounding::CLAMP,
138  "gamma5by11"
139  )
140  );
141 
142  gammaUpperReg8by11_.reset
143  (
145  (
146  gammaUpperReg8by11Table,
147  bounds::repeatableBounding::CLAMP,
148  "gamma8by11"
149  )
150  );
151 }
152 
153 
154 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
155 
157 {
158  kolmogorovLengthScale_ =
159  pow025
160  (
161  pow3
162  (
164  )
165  /popBal_.continuousTurbulence().epsilon()
166  );
167 }
168 
169 
170 void
172 (
173  volScalarField& binaryBreakupRate,
174  const label i,
175  const label j
176 )
177 {
178  const phaseModel& continuousPhase = popBal_.continuousPhase();
179  const sizeGroup& fi = popBal_.sizeGroups()[i];
180  const sizeGroup& fj = popBal_.sizeGroups()[j];
181 
182  const dimensionedScalar cf
183  (
184  pow(fi.x()/fj.x(), 2.0/3.0) + pow((1 - fi.x()/fj.x()), 2.0/3.0) - 1
185  );
186 
187  const volScalarField b
188  (
189  12.0*cf*popBal_.sigmaWithContinuousPhase(fi.phase())
190  /(
191  beta_*continuousPhase.rho()*pow(fj.d(), 5.0/3.0)
192  *pow(popBal_.continuousTurbulence().epsilon(), 2.0/3.0)
193  )
194  );
195 
196  const volScalarField xiMin(minEddyRatio_*kolmogorovLengthScale_/fj.d());
197 
198  const volScalarField tMin(b/pow(xiMin, 11.0/3.0));
199 
200  volScalarField integral(3.0/(11.0*pow(b, 8.0/11.0)));
201 
202  forAll(integral, celli)
203  {
204  integral[celli] *=
205  2.0*pow(b[celli], 3.0/11.0)*tgamma(5.0/11.0)
206  *(
207  gammaUpperReg5by11_()(b[celli])
208  - gammaUpperReg5by11_()(tMin[celli])
209  );
210  }
211 
212  binaryBreakupRate +=
213  C4_*(1 - popBal_.alphas())/fj.x()
214  *cbrt
215  (
216  popBal_.continuousTurbulence().epsilon()
217  /sqr(fj.d())
218  )
219  *integral;
220 }
221 
222 
223 // ************************************************************************* //
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:104
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
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:53
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:172
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:182
Foam::incGammaRatio_Q
scalar incGammaRatio_Q(const scalar a, const scalar x)
Normalized upper incomplete gamma function.
Definition: incGamma.C:223
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
LuoSvendsen.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:43
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:57
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:121
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::diameterModels::binaryBreakupModels::LuoSvendsen::correct
virtual void correct()
Correct diameter independent expressions.
Definition: LuoSvendsen.C:156