All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 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"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
32 (
33  const populationBalanceModel& popBal,
34  const dictionary& dict
35 )
36 :
37  populationBalance_(popBal),
38  kolmogorovLengthScale_
39  (
40  IOobject
41  (
42  "kolmogorovLengthScale",
43  populationBalance_.time().timeName(),
44  populationBalance_.mesh()
45  ),
46  populationBalance_.mesh(),
48  (
49  "kolmogorovLengthScale",
50  dimLength,
51  Zero
52  )
53  ),
54  shearStrainRate_
55  (
56  IOobject
57  (
58  "shearStrainRate",
59  populationBalance_.time().timeName(),
60  populationBalance_.mesh()
61  ),
62  populationBalance_.mesh(),
64  (
65  "shearStrainRate",
67  Zero
68  )
69  ),
70  eddyStrainRate_
71  (
72  IOobject
73  (
74  "eddyStrainRate",
75  populationBalance_.time().timeName(),
76  populationBalance_.mesh()
77  ),
78  populationBalance_.mesh(),
80  (
81  "eddyStrainRate",
83  Zero
84  )
85  )
86 {}
87 
88 
89 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
90 
92 {
94  pow025
95  (
98  );
99 
100  shearStrainRate_ =
101  sqrt(2.0)
103 
104  eddyStrainRate_ =
105  sqrt
106  (
110  );
111 
112  if (uTerminal_.empty())
113  {
114  const fvMesh& mesh = populationBalance_.mesh();
116  mesh.lookupObject<uniformDimensionedVectorField>("g");
117 
118  const dimensionedScalar nuc
119  (
120  "nuc",
121  dimViscosity,
123  );
124 
125  const dimensionedScalar rhoc
126  (
127  "rhoc",
128  dimDensity,
130  );
131 
132  const dimensionedScalar rhod
133  (
134  "rhod",
135  dimDensity,
136  gAverage(populationBalance_.sizeGroups()[1].phase().rho()())
137  );
138 
140  (
141  "sigma",
143  gAverage
144  (
146  (
147  populationBalance_.sizeGroups()[1].phase()
148  )()
149  )
150  );
151 
152  for(int m = 0; m < populationBalance_.sizeGroups().size(); m++)
153  {
154  const sizeGroup& f = populationBalance_.sizeGroups()[m];
155 
156  dimensionedScalar uTerminal("uTerminal", dimVelocity, 0.2);
157  dimensionedScalar Cd("Cd", dimless, 0.44);
158  dimensionedScalar CdEllipse("CdEllipse", dimless, 1);
159 
160  dimensionedScalar Re(uTerminal*f.dSph()/nuc);
161  const dimensionedScalar Eo
162  (
163  mag(g)*mag(rhoc - rhod)*sqr(f.dSph())/sigma
164  );
165 
168  const dimensionedScalar uTerminalX("uTerminalX", dimVelocity, 1e-5);
169  dimensionedScalar ReX("ReX", dimless, Re.value());
170  dimensionedScalar CdX("CdX", dimless, Cd.value());
171  dimensionedScalar dCd("dCd", Cd.dimensions()/dimVelocity, Zero);
172 
173  int n = 0;
174 
175  while(mag(F.value()) >= 1.0e-05 && n++ <= 20)
176  {
177  Re = uTerminal*f.dSph()/nuc;
178 
179  Cd =
180  pos0(1000 - Re)*24/Re*(1 + 0.1*pow(Re, 0.75))
181  + neg(1000 - Re)*0.44;
182 
183  CdEllipse = 0.6666*sqrt(Eo);
184 
185  Cd =
186  pos0(CdEllipse - Cd)
187  *min(CdEllipse.value(), 8.0/3.0)
188  + neg(CdEllipse - Cd)*Cd;
189 
190  F =
191  4.0/3.0*(rhoc - rhod)*mag(g)*f.dSph()
192  - rhoc*Cd*sqr(uTerminal);
193 
194  ReX = (uTerminal + uTerminalX)*f.dSph()/nuc;
195 
196  CdX =
197  pos0(1000 - ReX)
198  *24/ReX*(1 + 0.1*pow(ReX, 0.75))
199  + neg(1000 - ReX)*0.44;
200 
201  CdX =
202  pos0(CdEllipse - CdX)
203  *min(CdEllipse.value(), 2.66667)
204  + neg(CdEllipse - CdX)*CdX;
205 
206  dCd = (CdX - Cd)/uTerminalX;
207 
208  dF = -(2*rhoc*uTerminal*Cd + rhoc*sqr(uTerminal)*dCd);
209 
210  uTerminal -= F/dF;
211  }
212 
213  uTerminal_.append(new dimensionedScalar("uTerminal", uTerminal));
214 
215  Cd_.append(new dimensionedScalar("Cd", Cd));
216  }
217  }
218 }
219 
220 
221 // ************************************************************************* //
const populationBalanceModel & populationBalance_
Reference to the populationBalanceModel.
Definition: LiaoBase.H:67
const dimensionSet dimViscosity
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
const dimensionSet dimArea
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
UniformDimensionedField< vector > uniformDimensionedVectorField
virtual tmp< volVectorField > U() const =0
Return the velocity.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const phaseCompressible::momentumTransportModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual tmp< volScalarField > rho() const =0
Return the density field.
dimensionedScalar pow025(const dimensionedScalar &ds)
const dimensionSet dimless
dimensionedScalar neg(const dimensionedScalar &ds)
const dimensionSet dimLength
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
const phaseModel & continuousPhase() const
Return continuous phase.
Calculate the gradient of the given field.
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
Definition: phaseModel.H:121
LiaoBase(const populationBalanceModel &popBal, const dictionary &dict)
const dimensionSet dimDensity
PtrList< dimensionedScalar > Cd_
Drag coefficients.
Definition: LiaoBase.H:82
const Type & value() const
Return const reference to value.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
const dimensionSet dimForce
volScalarField rhoc(IOobject(rhocValue.name(), runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, rhocValue)
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionSet dimVelocity
const UPtrList< sizeGroup > & sizeGroups() const
Return the size groups belonging to this populationBalance.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Type gAverage(const FieldField< Field, Type > &f)
virtual tmp< volScalarField > mu() const =0
Dynamic viscosity of mixture [kg/m/s].
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
Definition: fluidThermo.C:73
dimensioned< scalar > mag(const dimensioned< Type > &)
const fvMesh & mesh() const
Return reference to the mesh.
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.
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
PtrList< dimensionedScalar > uTerminal_
Terminal velocities.
Definition: LiaoBase.H:79
virtual void precompute()
Precompute diameter independent expressions.
volScalarField kolmogorovLengthScale_
Kolmogorov length scale.
Definition: LiaoBase.H:70