BrownianCollisions.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) 2019-2020 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 "BrownianCollisions.H"
28 #include "fundamentalConstants.H"
29 #include "mathematicalConstants.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace diameterModels
36 {
37 namespace coalescenceModels
38 {
39  defineTypeNameAndDebug(BrownianCollisions, 0);
41  (
42  coalescenceModel,
43  BrownianCollisions,
44  dictionary
45  );
46 }
47 }
48 }
49 
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
58 (
59  const populationBalanceModel& popBal,
60  const dictionary& dict
61 )
62 :
63  coalescenceModel(popBal, dict),
64  A1_(dict.lookupOrDefault<scalar>("A1", 2.514)),
65  A2_(dict.lookupOrDefault<scalar>("A2", 0.8)),
66  A3_(dict.lookupOrDefault<scalar>("A3", 0.55)),
67  sigma_("sigma", dimLength, dict),
68  lambda_
69  (
70  IOobject
71  (
72  "lambda",
73  popBal_.time().timeName(),
74  popBal_.mesh()
75  ),
76  popBal_.mesh(),
78  (
79  "lambda",
80  dimLength,
81  Zero
82  )
83  )
84 {}
85 
86 
87 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88 
90 {
93 
94  lambda_ = k*T/(sqrt(2.0)*pi*p*sqr(sigma_));
95 }
96 
97 
98 void
101 (
102  volScalarField& coalescenceRate,
103  const label i,
104  const label j
105 )
106 {
107  const sizeGroup& fi = popBal_.sizeGroups()[i];
108  const sizeGroup& fj = popBal_.sizeGroups()[j];
109 
110  const volScalarField& T = popBal_.continuousPhase().thermo().T();
111  const volScalarField& mu = popBal_.continuousPhase().thermo().mu();
112 
113  const volScalarField Cci
114  (
115  1 + lambda_/fi.d()*(A1_ + A2_*exp(-A3_*fi.d()/lambda_))
116  );
117 
118  const volScalarField Ccj
119  (
120  1 + lambda_/fj.d()*(A1_ + A2_*exp(-A3_*fj.d()/lambda_))
121  );
122 
123  coalescenceRate += 8*k*T/(3*mu)*(fi.d() + fj.d())*(Cci/fi.d() + Ccj/fj.d());
124 }
125 
126 
127 // ************************************************************************* //
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual void precompute()
Precompute diameter independent expressions.
label k
Boltzmann constant.
virtual volScalarField & p()=0
Pressure [Pa].
Fundamental dimensioned constants.
Macros for easy insertion into run-time selection tables.
const dimensionSet dimLength
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dynamicFvMesh & mesh
dimensionedScalar exp(const dimensionedScalar &ds)
const phaseModel & continuousPhase() const
Return continuous phase.
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
Definition: phaseModel.H:121
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const dimensionedScalar mu
Atomic mass unit.
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
defineTypeNameAndDebug(combustionModel, 0)
const UPtrList< sizeGroup > & sizeGroups() const
Return the sizeGroups belonging to this populationBalance.
virtual const volScalarField & T() const =0
Temperature [K].
const dimensionedScalar k
Boltzmann constant.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
BrownianCollisions(const populationBalanceModel &popBal, const dictionary &dict)
virtual tmp< volScalarField > mu() const =0
Dynamic viscosity of mixture [kg/m/s].
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
Namespace for OpenFOAM.