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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "LuoSvendsen.H"
30 #include "phaseCompressibleTurbulenceModel.H"
31 #include "tableBounds.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace diameterModels
38 {
39 namespace binaryBreakupModels
40 {
41  defineTypeNameAndDebug(LuoSvendsen, 0);
43  (
44  binaryBreakupModel,
45  LuoSvendsen,
46  dictionary
47  );
48 }
49 }
50 }
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
56 (
57  const populationBalanceModel& popBal,
58  const dictionary& dict
59 )
60 :
61  binaryBreakupModel(popBal, dict),
62  gammaUpperReg2by11_(),
63  gammaUpperReg5by11_(),
64  gammaUpperReg8by11_(),
65  C4_(dimensionedScalar::lookupOrDefault("C4", dict, dimless, 0.923)),
66  beta_(dimensionedScalar::lookupOrDefault("beta", dict, dimless, 2.05)),
67  minEddyRatio_
68  (
69  dimensionedScalar::lookupOrDefault("minEddyRatio", dict, dimless, 11.4)
70  ),
71  kolmogorovLengthScale_
72  (
73  IOobject
74  (
75  "kolmogorovLengthScale",
76  popBal_.time().timeName(),
77  popBal_.mesh()
78  ),
79  popBal_.mesh(),
81  (
82  "kolmogorovLengthScale",
83  dimLength,
84  Zero
85  )
86  )
87 {
88  List<Tuple2<scalar, scalar>> gammaUpperReg2by11Table;
89  List<Tuple2<scalar, scalar>> gammaUpperReg5by11Table;
90  List<Tuple2<scalar, scalar>> gammaUpperReg8by11Table;
91 
92  gammaUpperReg2by11Table.append(Tuple2<scalar, scalar>(0.0, 1.0));
93  gammaUpperReg5by11Table.append(Tuple2<scalar, scalar>(0.0, 1.0));
94  gammaUpperReg8by11Table.append(Tuple2<scalar, scalar>(0.0, 1.0));
95 
96  for (scalar z = 1e-2; z <= 10.0; z = z + 1e-2)
97  {
98  Tuple2<scalar, scalar> gamma2by11
99  (
100  z,
101  incGammaRatio_Q(2.0/11.0, z)
102  );
103 
104  Tuple2<scalar, scalar> gamma5by11
105  (
106  z,
107  incGammaRatio_Q(5.0/11.0, z)
108  );
109 
110  Tuple2<scalar, scalar> gamma8by11
111  (
112  z,
113  incGammaRatio_Q(8.0/11.0, z)
114  );
115 
116  gammaUpperReg2by11Table.append(gamma2by11);
117  gammaUpperReg5by11Table.append(gamma5by11);
118  gammaUpperReg8by11Table.append(gamma8by11);
119  }
120 
121  gammaUpperReg2by11_.reset
122  (
124  (
125  gammaUpperReg2by11Table,
126  bounds::repeatableBounding::CLAMP,
127  "gamma2by11"
128  )
129  );
130 
131  gammaUpperReg5by11_.reset
132  (
134  (
135  gammaUpperReg5by11Table,
136  bounds::repeatableBounding::CLAMP,
137  "gamma5by11"
138  )
139  );
140 
141  gammaUpperReg8by11_.reset
142  (
144  (
145  gammaUpperReg8by11Table,
146  bounds::repeatableBounding::CLAMP,
147  "gamma8by11"
148  )
149  );
150 }
151 
152 
153 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
154 
156 {
157  kolmogorovLengthScale_ =
158  pow025
159  (
160  pow3
161  (
163  )
164  /popBal_.continuousTurbulence().epsilon()
165  );
166 }
167 
168 
169 void
171 (
172  volScalarField& binaryBreakupRate,
173  const label i,
174  const label j
175 )
176 {
177  const phaseModel& continuousPhase = popBal_.continuousPhase();
178  const sizeGroup& fi = popBal_.sizeGroups()[i];
179  const sizeGroup& fj = popBal_.sizeGroups()[j];
180 
181  const dimensionedScalar cf
182  (
183  pow(fi.x()/fj.x(), 2.0/3.0) + pow((1 - fi.x()/fj.x()), 2.0/3.0) - 1
184  );
185 
186  const volScalarField b
187  (
188  12.0*cf*popBal_.sigmaWithContinuousPhase(fi.phase())
189  /(
190  beta_*continuousPhase.rho()*pow(fj.d(), 5.0/3.0)
191  *pow(popBal_.continuousTurbulence().epsilon(), 2.0/3.0)
192  )
193  );
194 
195  const volScalarField xiMin(minEddyRatio_*kolmogorovLengthScale_/fj.d());
196 
197  const volScalarField tMin(b/pow(xiMin, 11.0/3.0));
198 
199  volScalarField integral(3.0/(11.0*pow(b, 8.0/11.0)));
200 
201  forAll(integral, celli)
202  {
203  integral[celli] *=
204  2.0*pow(b[celli], 3.0/11.0)*tgamma(5.0/11.0)
205  *(
206  gammaUpperReg5by11_()(b[celli])
207  - gammaUpperReg5by11_()(tMin[celli])
208  );
209  }
210 
211  binaryBreakupRate +=
212  C4_*(1 - popBal_.alphas())/fj.x()
213  *cbrt
214  (
215  popBal_.continuousTurbulence().epsilon()
216  /sqr(fj.d())
217  )
218  *integral;
219 }
220 
221 
222 // ************************************************************************* //
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:57
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.
Definition: zero.H:128
Foam::phaseModel::rho
virtual tmp< volScalarField > rho() const =0
Return the density field.
Foam::diameterModels::binaryBreakupModels::LuoSvendsen::addToBinaryBreakupRate
virtual void addToBinaryBreakupRate(volScalarField &binaryBreakupRate, const label i, const label j)
Add to binary breakupRate.
Definition: LuoSvendsen.C:171
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:290
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:133
Foam::phaseModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar kinematic viscosity.
LuoSvendsen.H
Foam::label
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:62
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::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:56
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:69
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: HashTable.H:102
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:1218
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::interpolationTable< scalar >
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:155