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-2025 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 {
41  (
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_("C4", dimless, dict, 0.923),
64  beta_("beta", dimless, dict, 2.05),
65  minEddyRatio_("minEddyRatio", dimless, dict, 11.4),
66  kolmogorovLengthScale_
67  (
68  IOobject
69  (
70  "kolmogorovLengthScale",
71  popBal_.time().name(),
72  popBal_.mesh()
73  ),
74  popBal_.mesh(),
76  (
77  "kolmogorovLengthScale",
78  dimLength,
79  Zero
80  )
81  )
82 {
83  List<Tuple2<scalar, scalar>> gammaUpperReg2by11Table;
84  List<Tuple2<scalar, scalar>> gammaUpperReg5by11Table;
85  List<Tuple2<scalar, scalar>> gammaUpperReg8by11Table;
86 
87  gammaUpperReg2by11Table.append(Tuple2<scalar, scalar>(0.0, 1.0));
88  gammaUpperReg5by11Table.append(Tuple2<scalar, scalar>(0.0, 1.0));
89  gammaUpperReg8by11Table.append(Tuple2<scalar, scalar>(0.0, 1.0));
90 
91  for (scalar z = 1e-2; z <= 10.0; z = z + 1e-2)
92  {
93  Tuple2<scalar, scalar> gamma2by11
94  (
95  z,
96  incGammaRatio_Q(2.0/11.0, z)
97  );
98 
99  Tuple2<scalar, scalar> gamma5by11
100  (
101  z,
102  incGammaRatio_Q(5.0/11.0, z)
103  );
104 
105  Tuple2<scalar, scalar> gamma8by11
106  (
107  z,
108  incGammaRatio_Q(8.0/11.0, z)
109  );
110 
111  gammaUpperReg2by11Table.append(gamma2by11);
112  gammaUpperReg5by11Table.append(gamma5by11);
113  gammaUpperReg8by11Table.append(gamma8by11);
114  }
115 
116  gammaUpperReg2by11_ =
118  (
119  "gamma2by11",
121  linearInterpolationWeights::typeName,
122  autoPtr<TableReader<scalar>>(nullptr),
123  gammaUpperReg2by11Table
124  );
125 
126  gammaUpperReg5by11_ =
128  (
129  "gamma5by11",
131  linearInterpolationWeights::typeName,
132  autoPtr<TableReader<scalar>>(nullptr),
133  gammaUpperReg5by11Table
134  );
135 
136  gammaUpperReg8by11_ =
138  (
139  "gamma8by11",
141  linearInterpolationWeights::typeName,
142  autoPtr<TableReader<scalar>>(nullptr),
143  gammaUpperReg8by11Table
144  );
145 }
146 
147 
148 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
149 
151 {
152  kolmogorovLengthScale_ =
153  pow025
154  (
155  pow3
156  (
157  popBal_.continuousPhase().fluidThermo().nu()
158  )
159  /popBal_.continuousTurbulence().epsilon()
160  );
161 }
162 
163 
164 void
166 (
167  volScalarField::Internal& binaryBreakupRate,
168  const label i,
169  const label j
170 )
171 {
172  const sizeGroup& fi = popBal_.sizeGroups()[i];
173  const sizeGroup& fj = popBal_.sizeGroups()[j];
174 
175  const volScalarField::Internal& rhoc = popBal_.continuousPhase().rho();
176 
177  tmp<volScalarField> tsigma(popBal_.sigmaWithContinuousPhase(fi.phase()));
178  const volScalarField::Internal& sigma = tsigma();
179 
180  tmp<volScalarField> tepsilonc(popBal_.continuousTurbulence().epsilon());
181  const volScalarField::Internal& epsilonc = tepsilonc();
182 
183  const dimensionedScalar cf
184  (
185  pow(fi.x()/fj.x(), 2.0/3.0) + pow((1 - fi.x()/fj.x()), 2.0/3.0) - 1
186  );
187 
189  (
190  12*cf*sigma
191  /(beta_*rhoc*pow(fj.dSph(), 5.0/3.0)*pow(epsilonc, 2.0/3.0))
192  );
193 
194  const volScalarField::Internal xiMin
195  (
196  minEddyRatio_*kolmogorovLengthScale_/fj.dSph()
197  );
198 
199  const volScalarField::Internal tMin(b/pow(xiMin, 11.0/3.0));
200 
201  volScalarField::Internal integral(3/(11*pow(b, 8.0/11.0)));
202  forAll(integral, celli)
203  {
204  integral[celli] *=
205  2
206  *pow(b[celli], 3.0/11.0)
207  *tgamma(5.0/11.0)
208  *(
209  gammaUpperReg5by11_->value(b[celli])
210  - gammaUpperReg5by11_->value(tMin[celli])
211  );
212  }
213 
214  binaryBreakupRate +=
215  C4_
216  *(1 - popBal_.alphas()())
217  /fj.x()
218  *cbrt(epsilonc/sqr(fj.dSph()))
219  *integral;
220 }
221 
222 
223 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
Base class to read table data for tables.
Definition: TableReader.H:51
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Base class for binary breakup models that provide a breakup rate between a size class pair directly,...
Model of Luo and Svendsen (1996). The breakup rate is calculated by.
Definition: LuoSvendsen.H:231
virtual void precompute()
Precompute diameter independent expressions.
Definition: LuoSvendsen.C:150
virtual void addToBinaryBreakupRate(volScalarField::Internal &binaryBreakupRate, const label i, const label j)
Add to binary breakupRate.
Definition: LuoSvendsen.C:166
LuoSvendsen(const populationBalanceModel &popBal, const dictionary &dict)
Definition: LuoSvendsen.C:54
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:96
const dimensionedScalar & dSph() const
Return representative spherical diameter of the sizeGroup.
Definition: sizeGroupI.H:51
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition: sizeGroupI.H:58
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:37
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
A class for managing temporary objects.
Definition: tmp.H:55
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
volScalarField & b
Definition: createFields.H:27
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
addToRunTimeSelectionTable(binaryBreakupModel, LehrMilliesMewes, dictionary)
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
const doubleScalar e
Definition: doubleScalar.H:106
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
void pow025(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
const dimensionSet dimless
const dimensionSet dimLength
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void cbrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
scalar incGammaRatio_Q(const scalar a, const scalar x)
Normalised upper incomplete gamma function.
Definition: incGamma.C:221
dictionary dict