All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
BrownianCollisions.H
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-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 Class
25  Foam::diameterModels::coalescenceModels::BrownianCollisions
26 
27 Description
28  Model describing coagulation due to Brownian motion. Utilises collisional
29  diameters and the Cunningham slip correction. The slip correction
30  coefficient is implemented in the following form:
31 
32  \f[
33  C_{c_i} = 1 + \lambda [A_1 + A_2 \exp(-A_3 d_i/\lambda)]/d_i\,.
34  \f]
35 
36  The coefficients default to the values proposed by Davis (1945). The mean
37  free path is computed by
38 
39  \f[
40  \lambda = \frac{kT}{\sqrt{2} \pi p \sigma^{2}}\,.
41  \f]
42 
43  \vartable
44  A_1 | Coefficient [-]
45  A_2 | Coefficient [-]
46  A_3 | Coefficient [-]
47  \sigma | Lennard-Jones parameter [m]
48  \endvartable
49 
50  Reference:
51  \verbatim
52  Davies, C. N. (1945).
53  Definitive equations for the fluid resistance of spheres.
54  Proceedings of the Physical Society, 57(4), 259.
55  \endverbatim
56 
57 Usage
58  \table
59  Property | Description | Required | Default value
60  A1 | Coefficient A1 | no | 2.514
61  A2 | Coefficient A2 | no | 0.8
62  A3 | Coefficient A2 | no | 0.55
63  sigma | Lennard-Jones parameter | yes | none
64  \endtable
65 
66 SourceFiles
67  BrownianCollisions.C
68 
69 \*---------------------------------------------------------------------------*/
70 
71 #ifndef BrownianCollisions_H
72 #define BrownianCollisions_H
73 
74 #include "coalescenceModel.H"
75 
76 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
77 
78 namespace Foam
79 {
80 namespace diameterModels
81 {
82 namespace coalescenceModels
83 {
84 
85 /*---------------------------------------------------------------------------*\
86  Class BrownianCollisions Declaration
87 \*---------------------------------------------------------------------------*/
88 
89 class BrownianCollisions
90 :
91  public coalescenceModel
92 {
93  // Private Data
94 
95  //- Cunningham slip correction coefficient A1
96  const dimensionedScalar A1_;
97 
98  //- Cunningham slip correction coefficient A2
99  const dimensionedScalar A2_;
100 
101  //- Cunningham slip correction coefficient A3
102  const dimensionedScalar A3_;
103 
104  //- Lennard-Jones sigma parameter
105  const dimensionedScalar sigma_;
106 
107  //- Mean free path
108  volScalarField lambda_;
109 
110 
111 public:
112 
113  //- Runtime type information
114  TypeName("BrownianCollisions");
115 
116  // Constructor
117 
119  (
120  const populationBalanceModel& popBal,
121  const dictionary& dict
122  );
123 
124 
125  //- Destructor
126  virtual ~BrownianCollisions()
127  {}
128 
130  // Member Functions
131 
132  //- Precompute diameter independent expressions
133  virtual void precompute();
134 
135  //- Add to coalescenceRate
136  virtual void addToCoalescenceRate
137  (
138  volScalarField& coalescenceRate,
139  const label i,
140  const label j
141  );
142 };
143 
144 
145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146 
147 } // End namespace coalescenceModels
148 } // End namespace diameterModels
149 } // End namespace Foam
150 
151 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152 
153 #endif
154 
155 // ************************************************************************* //
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
virtual void precompute()
Precompute diameter independent expressions.
TypeName("BrownianCollisions")
Runtime type information.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
BrownianCollisions(const populationBalanceModel &popBal, const dictionary &dict)
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
Namespace for OpenFOAM.