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-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 "PrinceBlanch.H"
27 #include "fvcGrad.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace diameterModels
37 {
38 namespace coalescenceModels
39 {
42  (
46  );
47 }
48 }
49 }
50 
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
58 (
59  const populationBalanceModel& popBal,
60  const dictionary& dict
61 )
62 :
63  coalescenceModel(popBal, dict),
64  C1_("C1", dimless, dict, 0.356),
65  h0_("h0", dimLength, dict, 1e-4),
66  hf_("hf", dimLength, dict, 1e-8),
67  turbulence_(dict.lookup("turbulence")),
68  buoyancy_(dict.lookup("buoyancy")),
69  laminarShear_(dict.lookup("laminarShear"))
70 {
71  if (laminarShear_)
72  {
73  shearStrainRate_.set
74  (
76  (
77  IOobject
78  (
79  "shearStrainRate",
80  popBal_.time().name(),
81  popBal_.mesh()
82  ),
83  popBal_.mesh(),
85  (
86  "shearStrainRate",
88  Zero
89  )
90  )
91  );
92  }
93 }
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
99 {
100  if (laminarShear_)
101  {
102  shearStrainRate_() =
103  sqrt(2.0)*mag(symm(fvc::grad(popBal_.continuousPhase().U())));
104  }
105 }
106 
107 
110 (
111  volScalarField::Internal& coalescenceRate,
112  const label i,
113  const label j
114 )
115 {
116  const sizeGroup& fi = popBal_.sizeGroups()[i];
117  const sizeGroup& fj = popBal_.sizeGroups()[j];
118 
119  const volScalarField::Internal& rhoc = popBal_.continuousPhase().rho();
120 
121  tmp<volScalarField> tsigma(popBal_.sigmaWithContinuousPhase(fi.phase()));
122  const volScalarField::Internal& sigma = tsigma();
123 
124  tmp<volScalarField> tepsilonc(popBal_.continuousTurbulence().epsilon());
125  const volScalarField::Internal& epsilonc = tepsilonc();
126 
128  popBal_.mesh().lookupObject<uniformDimensionedVectorField>("g");
129 
130  const dimensionedScalar rij(1/(1/fi.dSph() + 1/fj.dSph()));
131 
132  const volScalarField::Internal collisionEfficiency
133  (
134  exp
135  (
136  - sqrt(pow3(rij)*rhoc/(16*sigma))
137  *log(h0_/hf_)
138  *cbrt(epsilonc)
139  /pow(rij, 2.0/3.0)
140  )
141  );
142 
143  if (turbulence_)
144  {
145  coalescenceRate +=
146  (
147  C1_
148  *pi
149  *sqr(fi.dSph() + fj.dSph())
150  *cbrt(epsilonc)
151  *sqrt(pow(fi.dSph(), 2.0/3.0) + pow(fj.dSph(), 2.0/3.0))
152  )
153  *collisionEfficiency;
154  }
155 
156  if (buoyancy_)
157  {
158  const dimensionedScalar Sij(pi/4*sqr(fi.dSph() + fj.dSph()));
159 
160  coalescenceRate +=
161  Sij
162  *mag
163  (
164  sqrt(2.14*sigma/(rhoc*fi.dSph()) + 0.505*mag(g)*fi.dSph())
165  - sqrt(2.14*sigma/(rhoc*fj.dSph()) + 0.505*mag(g)*fj.dSph())
166  )
167  *collisionEfficiency;
168  }
169 
170  if (laminarShear_)
171  {
172  coalescenceRate +=
173  pow3(fi.dSph() + fj.dSph())/6
174  *shearStrainRate_()*collisionEfficiency;
175  }
176 }
177 
178 
179 // ************************************************************************* //
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
const Time & time() const
Return time.
Definition: IOobject.C:315
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 precompute()
Precompute diameter independent expressions.
Definition: PrinceBlanch.C:98
virtual void addToCoalescenceRate(volScalarField::Internal &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
Definition: PrinceBlanch.C:110
PrinceBlanch(const populationBalanceModel &popBal, const dictionary &dict)
Definition: PrinceBlanch.C:58
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:96
const dimensionedScalar & dSph() const
Return representative spherical diameter of the sizeGroup.
Definition: sizeGroupI.H:51
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
const word & name() const
Return const reference to name.
A class for managing temporary objects.
Definition: tmp.H:55
Calculate the gradient of the given field.
const dimensionedScalar e
Elementary charge.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
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
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
const dimensionSet dimless
const dimensionSet dimLength
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
dimensionedScalar log(const dimensionedScalar &ds)
void symm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
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)
const dimensionSet dimVelocity
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dictionary dict