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-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 "PrinceBlanch.H"
28 #include "mathematicalConstants.H"
30 #include "fvcGrad.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace diameterModels
37 {
38 namespace coalescenceModels
39 {
40  defineTypeNameAndDebug(PrinceBlanch, 0);
42  (
43  coalescenceModel,
44  PrinceBlanch,
45  dictionary
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_(dimensionedScalar::lookupOrDefault("C1", dict, dimless, 0.356)),
65  h0_
66  (
67  dimensionedScalar::lookupOrDefault
68  (
69  "h0",
70  dict,
71  dimLength,
72  1e-4
73  )
74  ),
75  hf_
76  (
77  dimensionedScalar::lookupOrDefault
78  (
79  "hf",
80  dict,
81  dimLength,
82  1e-8
83  )
84  ),
85  turbulence_(dict.lookup("turbulence")),
86  buoyancy_(dict.lookup("buoyancy")),
87  laminarShear_(dict.lookup("laminarShear"))
88 {
89  if (laminarShear_)
90  {
91  shearStrainRate_.set
92  (
93  new volScalarField
94  (
95  IOobject
96  (
97  "shearStrainRate",
98  popBal_.time().timeName(),
99  popBal_.mesh()
100  ),
101  popBal_.mesh(),
103  (
104  "shearStrainRate",
106  Zero
107  )
108  )
109  );
110  }
111 }
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 
117 {
118  if (laminarShear_)
119  {
120  shearStrainRate_() =
122  }
123 }
124 
125 
128 (
129  volScalarField& coalescenceRate,
130  const label i,
131  const label j
132 )
133 {
134  const phaseModel& continuousPhase = popBal_.continuousPhase();
135  const sizeGroup& fi = popBal_.sizeGroups()[i];
136  const sizeGroup& fj = popBal_.sizeGroups()[j];
139 
140  const dimensionedScalar rij(1/(1/fi.dSph() + 1/fj.dSph()));
141 
142  const volScalarField collisionEfficiency
143  (
144  exp
145  (
146  - sqrt
147  (
148  pow3(rij)*continuousPhase.rho()
149  /(16*popBal_.sigmaWithContinuousPhase(fi.phase()))
150  )
151  *log(h0_/hf_)
152  *cbrt(popBal_.continuousTurbulence().epsilon())/pow(rij, 2.0/3.0)
153  )
154  );
155 
156  if (turbulence_)
157  {
158  coalescenceRate +=
159  (
160  C1_*pi*sqr(fi.dSph() + fj.dSph())
162  *sqrt(pow(fi.dSph(), 2.0/3.0) + pow(fj.dSph(), 2.0/3.0))
163  )
164  *collisionEfficiency;
165  }
166 
167  if (buoyancy_)
168  {
169  const dimensionedScalar Sij(pi/4*sqr(fi.dSph() + fj.dSph()));
170 
171  coalescenceRate +=
172  (
173  Sij
174  *mag
175  (
176  sqrt
177  (
178  2.14*popBal_.sigmaWithContinuousPhase(fi.phase())
179  /(continuousPhase.rho()*fi.dSph())
180  + 0.505*mag(g)*fi.dSph()
181  )
182  - sqrt
183  (
184  2.14*popBal_.sigmaWithContinuousPhase(fi.phase())
185  /(continuousPhase.rho()*fj.dSph())
186  + 0.505*mag(g)*fj.dSph()
187  )
188  )
189  )
190  *collisionEfficiency;
191  }
192 
193  if (laminarShear_)
194  {
195  coalescenceRate +=
196  pow3(fi.dSph() + fj.dSph())/6
197  *shearStrainRate_()*collisionEfficiency;
198  }
199 }
200 
201 
202 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
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
dimensionedScalar log(const dimensionedScalar &ds)
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)
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:636
const dimensionSet dimless
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Macros for easy insertion into run-time selection tables.
const dimensionSet dimLength
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensionedScalar exp(const dimensionedScalar &ds)
const phaseModel & continuousPhase() const
Return continuous phase.
Calculate the gradient of the given field.
PrinceBlanch(const populationBalanceModel &popBal, const dictionary &dict)
dimensionedScalar cbrt(const dimensionedScalar &ds)
static const zero Zero
Definition: zero.H:97
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const dimensionSet dimVelocity
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
defineTypeNameAndDebug(combustionModel, 0)
const UPtrList< sizeGroup > & sizeGroups() const
Return the sizeGroups belonging to this populationBalance.
virtual void precompute()
Precompute diameter independent expressions.
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.
const Time & time() const
Return time.
Definition: IOobject.C:318
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.
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
Namespace for OpenFOAM.