LiaoBase.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) 2021-2024 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 "LiaoBase.H"
27 #include "fvcGrad.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
33 (
34  const populationBalanceModel& popBal,
35  const dictionary& dict
36 )
37 :
38  populationBalance_(popBal),
39  kolmogorovLengthScale_
40  (
41  IOobject
42  (
43  "kolmogorovLengthScale",
44  populationBalance_.time().name(),
45  populationBalance_.mesh()
46  ),
47  populationBalance_.mesh(),
49  (
50  "kolmogorovLengthScale",
51  dimLength,
52  Zero
53  )
54  ),
55  shearStrainRate_
56  (
57  IOobject
58  (
59  "shearStrainRate",
60  populationBalance_.time().name(),
61  populationBalance_.mesh()
62  ),
63  populationBalance_.mesh(),
65  (
66  "shearStrainRate",
68  Zero
69  )
70  ),
71  eddyStrainRate_
72  (
73  IOobject
74  (
75  "eddyStrainRate",
76  populationBalance_.time().name(),
77  populationBalance_.mesh()
78  ),
79  populationBalance_.mesh(),
81  (
82  "eddyStrainRate",
84  Zero
85  )
86  )
87 {}
88 
89 
90 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91 
93 {
94  kolmogorovLengthScale_ =
95  pow025
96  (
97  pow3(populationBalance_.continuousPhase().fluidThermo().nu())
98  /populationBalance_.continuousTurbulence().epsilon()
99  );
100 
101  shearStrainRate_ =
102  sqrt(2.0)
103  *mag(symm(fvc::grad(populationBalance_.continuousPhase().U())));
104 
105  eddyStrainRate_ =
106  sqrt
107  (
108  populationBalance_.continuousPhase().rho()
109  *populationBalance_.continuousTurbulence().epsilon()
110  /populationBalance_.continuousPhase().fluidThermo().mu()
111  );
112 
113  if (uTerminal_.empty())
114  {
115  const fvMesh& mesh = populationBalance_.mesh();
118 
119  const dimensionedScalar nuc
120  (
121  "nuc",
123  gAverage(populationBalance_.continuousPhase().fluidThermo().nu()())
124  );
125 
126  const dimensionedScalar rhoc
127  (
128  "rhoc",
129  dimDensity,
130  gAverage(populationBalance_.continuousPhase().rho())
131  );
132 
133  const dimensionedScalar rhod
134  (
135  "rhod",
136  dimDensity,
137  gAverage(populationBalance_.sizeGroups()[1].phase().rho())
138  );
139 
141  (
142  "sigma",
144  gAverage
145  (
146  populationBalance_.sigmaWithContinuousPhase
147  (
148  populationBalance_.sizeGroups()[1].phase()
149  )()
150  )
151  );
152 
153  for(int m = 0; m < populationBalance_.sizeGroups().size(); m++)
154  {
155  const sizeGroup& f = populationBalance_.sizeGroups()[m];
156 
157  dimensionedScalar uTerminal("uTerminal", dimVelocity, 0.2);
158  dimensionedScalar Cd("Cd", dimless, 0.44);
159  dimensionedScalar CdEllipse("CdEllipse", dimless, 1);
160 
161  dimensionedScalar Re(uTerminal*f.dSph()/nuc);
162  const dimensionedScalar Eo
163  (
164  mag(g)*mag(rhoc - rhod)*sqr(f.dSph())/sigma
165  );
166 
169  const dimensionedScalar uTerminalX("uTerminalX", dimVelocity, 1e-5);
170  dimensionedScalar ReX("ReX", dimless, Re.value());
171  dimensionedScalar CdX("CdX", dimless, Cd.value());
172  dimensionedScalar dCd("dCd", Cd.dimensions()/dimVelocity, Zero);
173 
174  int n = 0;
175 
176  while(mag(F.value()) >= 1.0e-05 && n++ <= 20)
177  {
178  Re = uTerminal*f.dSph()/nuc;
179 
180  Cd =
181  pos0(1000 - Re)*24/Re*(1 + 0.1*pow(Re, 0.75))
182  + neg(1000 - Re)*0.44;
183 
184  CdEllipse = 0.6666*sqrt(Eo);
185 
186  Cd =
187  pos0(CdEllipse - Cd)
188  *min(CdEllipse.value(), 8.0/3.0)
189  + neg(CdEllipse - Cd)*Cd;
190 
191  F =
192  4.0/3.0*(rhoc - rhod)*mag(g)*f.dSph()
193  - rhoc*Cd*sqr(uTerminal);
194 
195  ReX = (uTerminal + uTerminalX)*f.dSph()/nuc;
196 
197  CdX =
198  pos0(1000 - ReX)
199  *24/ReX*(1 + 0.1*pow(ReX, 0.75))
200  + neg(1000 - ReX)*0.44;
201 
202  CdX =
203  pos0(CdEllipse - CdX)
204  *min(CdEllipse.value(), 2.66667)
205  + neg(CdEllipse - CdX)*CdX;
206 
207  dCd = (CdX - Cd)/uTerminalX;
208 
209  dF = -(2*rhoc*uTerminal*Cd + rhoc*sqr(uTerminal)*dCd);
210 
211  uTerminal -= F/dF;
212  }
213 
214  uTerminal_.append(new dimensionedScalar("uTerminal", uTerminal));
215 
216  Cd_.append(new dimensionedScalar("Cd", Cd));
217  }
218  }
219 }
220 
221 
222 // ************************************************************************* //
label n
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
virtual void precompute()
Precompute diameter independent expressions.
Definition: LiaoBase.C:92
LiaoBase(const populationBalanceModel &popBal, const dictionary &dict)
Definition: LiaoBase.C:33
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:102
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const dimensionSet & dimensions() const
Return const reference to dimensions.
const Type & value() const
Return const reference to value.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:441
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
Calculate the gradient of the given field.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
const dimensionedScalar F
Faraday constant: default SI units: [C/kmol].
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
static const zero Zero
Definition: zero.H:97
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
const doubleScalar e
Definition: doubleScalar.H:106
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimKinematicViscosity
const dimensionSet dimless
const dimensionSet dimLength
const dimensionSet dimForce
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimDensity
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensionedScalar neg(const dimensionedScalar &ds)
Type gAverage(const FieldField< Field, Type > &f)
const dimensionSet dimVelocity
const dimensionSet dimArea
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pow025(const dimensionedScalar &ds)
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
labelList f(nPoints)
dictionary dict