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 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"
29 #include "phaseCompressibleTurbulenceModel.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace diameterModels
36 {
37 namespace coalescenceModels
38 {
39  defineTypeNameAndDebug(PrinceBlanch, 0);
41  (
42  coalescenceModel,
43  PrinceBlanch,
44  dictionary
45  );
46 }
47 }
48 }
49 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
56 (
57  const populationBalanceModel& popBal,
58  const dictionary& dict
59 )
60 :
61  coalescenceModel(popBal, dict),
62  C1_("C1", dimless, dict.lookupOrDefault<scalar>("C1", 0.356)),
63  h0_("h0", dimLength, dict.lookupOrDefault<scalar>("h0", 1e-4)),
64  hf_("hf", dimLength, dict.lookupOrDefault<scalar>("h0", 1e-8)),
65  turbulentCollisions_(dict.lookup("turbulentCollisions")),
66  buoyantCollisions_(dict.lookup("buoyantCollisions")),
67  laminarShearCollisions_(dict.lookup("laminarShearCollisions"))
68 {}
69 
70 
71 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
72 
75 (
76  volScalarField& coalescenceRate,
77  const label i,
78  const label j
79 )
80 {
81  const phaseModel& continuousPhase = popBal_.continuousPhase();
82  const sizeGroup& fi = *popBal_.sizeGroups()[i];
83  const sizeGroup& fj = *popBal_.sizeGroups()[j];
86 
87  dimensionedScalar rij(1.0/(1.0/fi.d() + 1.0/fj.d()));
88 
89  const volScalarField collisionEfficiency
90  (
91  exp
92  (
93  - sqrt
94  (
95  pow3(rij)*continuousPhase.rho()
96  /(16.0*popBal_.sigmaWithContinuousPhase(fi.phase()))
97  )
98  *log(h0_/hf_)
99  *cbrt(popBal_.continuousTurbulence().epsilon())/pow(rij, 2.0/3.0)
100  )
101  );
102 
103  if (turbulentCollisions_)
104  {
105  coalescenceRate +=
106  (
107  C1_*pi*sqr(fi.d() + fj.d())
109  *sqrt(pow(fi.d(), 2.0/3.0) + pow(fj.d(), 2.0/3.0))
110  )
111  *collisionEfficiency;
112  }
113 
114  if (buoyantCollisions_)
115  {
116  dimensionedScalar Sij(pi/4.0*sqr(fi.d() + fj.d()));
117 
118  coalescenceRate +=
119  (
120  Sij
121  *mag
122  (
123  sqrt
124  (
125  2.14*popBal_.sigmaWithContinuousPhase(fi.phase())
126  /(continuousPhase.rho()*fi.d()) + 0.505*mag(g)*fi.d()
127  )
128  - sqrt
129  (
130  2.14*popBal_.sigmaWithContinuousPhase(fi.phase())
131  /(continuousPhase.rho()*fj.d()) + 0.505*mag(g)*fj.d()
132  )
133  )
134  )
135  *collisionEfficiency;
136  }
137 
138  if (laminarShearCollisions_)
139  {
141  << "Laminar shear collision contribution not implemented for "
142  << this->type() << " coalescence model."
143  << exit(FatalError);
144  }
145 }
146 
147 
148 // ************************************************************************* //
dictionary dict
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)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
UniformDimensionedField< vector > uniformDimensionedVectorField
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
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.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensionedScalar exp(const dimensionedScalar &ds)
const phaseModel & continuousPhase() const
Return continuous phase.
PrinceBlanch(const populationBalanceModel &popBal, const dictionary &dict)
dimensionedScalar cbrt(const dimensionedScalar &ds)
const List< sizeGroup * > & sizeGroups() const
Return the sizeGroups belonging to this populationBalance.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const phaseCompressibleTurbulenceModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
defineTypeNameAndDebug(combustionModel, 0)
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
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:98
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
Namespace for OpenFOAM.