PrinceBlanch.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) 2018-2023 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 "PrinceBlanch.H"
28 #include "fvcGrad.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace diameterModels
35 {
36 namespace coalescenceModels
37 {
40  (
44  );
45 }
46 }
47 }
48 
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
56 (
57  const populationBalanceModel& popBal,
58  const dictionary& dict
59 )
60 :
61  coalescenceModel(popBal, dict),
62  C1_(dimensionedScalar::lookupOrDefault("C1", dict, dimless, 0.356)),
63  h0_
64  (
65  dimensionedScalar::lookupOrDefault
66  (
67  "h0",
68  dict,
69  dimLength,
70  1e-4
71  )
72  ),
73  hf_
74  (
75  dimensionedScalar::lookupOrDefault
76  (
77  "hf",
78  dict,
79  dimLength,
80  1e-8
81  )
82  ),
83  turbulence_(dict.lookup("turbulence")),
84  buoyancy_(dict.lookup("buoyancy")),
85  laminarShear_(dict.lookup("laminarShear"))
86 {
87  if (laminarShear_)
88  {
89  shearStrainRate_.set
90  (
91  new volScalarField
92  (
93  IOobject
94  (
95  "shearStrainRate",
96  popBal_.time().name(),
97  popBal_.mesh()
98  ),
99  popBal_.mesh(),
101  (
102  "shearStrainRate",
104  Zero
105  )
106  )
107  );
108  }
109 }
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
115 {
116  if (laminarShear_)
117  {
118  shearStrainRate_() =
119  sqrt(2.0)*mag(symm(fvc::grad(popBal_.continuousPhase().U())));
120  }
121 }
122 
123 
126 (
127  volScalarField& coalescenceRate,
128  const label i,
129  const label j
130 )
131 {
132  const phaseModel& continuousPhase = popBal_.continuousPhase();
133  const sizeGroup& fi = popBal_.sizeGroups()[i];
134  const sizeGroup& fj = popBal_.sizeGroups()[j];
136  popBal_.mesh().lookupObject<uniformDimensionedVectorField>("g");
137 
138  const dimensionedScalar rij(1/(1/fi.dSph() + 1/fj.dSph()));
139 
140  const volScalarField collisionEfficiency
141  (
142  exp
143  (
144  - sqrt
145  (
146  pow3(rij)*continuousPhase.rho()
147  /(16*popBal_.sigmaWithContinuousPhase(fi.phase()))
148  )
149  *log(h0_/hf_)
150  *cbrt(popBal_.continuousTurbulence().epsilon())/pow(rij, 2.0/3.0)
151  )
152  );
153 
154  if (turbulence_)
155  {
156  coalescenceRate +=
157  (
158  C1_*pi*sqr(fi.dSph() + fj.dSph())
159  *cbrt(popBal_.continuousTurbulence().epsilon())
160  *sqrt(pow(fi.dSph(), 2.0/3.0) + pow(fj.dSph(), 2.0/3.0))
161  )
162  *collisionEfficiency;
163  }
164 
165  if (buoyancy_)
166  {
167  const dimensionedScalar Sij(pi/4*sqr(fi.dSph() + fj.dSph()));
168 
169  coalescenceRate +=
170  (
171  Sij
172  *mag
173  (
174  sqrt
175  (
176  2.14*popBal_.sigmaWithContinuousPhase(fi.phase())
177  /(continuousPhase.rho()*fi.dSph())
178  + 0.505*mag(g)*fi.dSph()
179  )
180  - sqrt
181  (
182  2.14*popBal_.sigmaWithContinuousPhase(fi.phase())
183  /(continuousPhase.rho()*fj.dSph())
184  + 0.505*mag(g)*fj.dSph()
185  )
186  )
187  )
188  *collisionEfficiency;
189  }
190 
191  if (laminarShear_)
192  {
193  coalescenceRate +=
194  pow3(fi.dSph() + fj.dSph())/6
195  *shearStrainRate_()*collisionEfficiency;
196  }
197 }
198 
199 
200 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const Time & time() const
Return time.
Definition: IOobject.C:318
Base class for coalescence models.
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
Model of Prince and Blanch (1990). The coalescence rate is calculated by.
Definition: PrinceBlanch.H:262
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
Definition: PrinceBlanch.C:126
virtual void precompute()
Precompute diameter independent expressions.
Definition: PrinceBlanch.C:114
PrinceBlanch(const populationBalanceModel &popBal, const dictionary &dict)
Definition: PrinceBlanch.C:56
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
const fvMesh & mesh() const
Return reference to the mesh.
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:102
const dimensionedScalar & dSph() const
Return representative spherical diameter of the sizeGroup.
Definition: sizeGroupI.H:50
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:36
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
const word & name() const
Return const reference to name.
virtual const volScalarField & rho() const =0
Return the density field.
Calculate the gradient of the given field.
const dimensionedScalar e
Elementary charge.
addToRunTimeSelectionTable(coalescenceModel, AdachiStuartFokkink, dictionary)
defineTypeNameAndDebug(AdachiStuartFokkink, 0)
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedScalar exp(const dimensionedScalar &ds)
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimless
const dimensionSet dimLength
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionSet dimVelocity
dimensionedScalar cbrt(const dimensionedScalar &ds)
dictionary dict