DahnekeInterpolation.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-2023 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::DahnekeInterpolation
26 
27 Description
28  Interpolation formula of Dahneke (1983) as presented by Otto et al. (1999).
29  Utilises collisional diameters.
30 
31  References:
32  \verbatim
33  Dahneke, B. (1983).
34  Simple kinetic theory of Brownian diffusion in vapors and aerosols.
35  In Theory of dispersed multiphase flow (pp. 97-133). Academic Press.
36  \endverbatim
37 
38  \verbatim
39  Otto, E., Fissan, H., Park, S. H., & Lee, K. W. (1999).
40  The log-normal size distribution theory of Brownian aerosol coagulation
41  for the entire particle size range: part II—analytical solution using
42  Dahneke’s coagulation kernel.
43  Journal of aerosol science, 30(1), 17-34.
44  \endverbatim
45 
46 SourceFiles
47  DahnekeInterpolation.C
48 
49 \*---------------------------------------------------------------------------*/
50 
51 #ifndef DahnekeInterpolation_H
52 #define DahnekeInterpolation_H
53 
54 #include "coalescenceModel.H"
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 namespace Foam
59 {
60 namespace diameterModels
61 {
62 namespace coalescenceModels
63 {
64 
65 class BrownianCollisions;
66 class ballisticCollisions;
67 
68 /*---------------------------------------------------------------------------*\
69  Class DahnekeInterpolation Declaration
70 \*---------------------------------------------------------------------------*/
71 
73 :
74  public coalescenceModel
75 {
76  // Private Data
77 
78  //- Model for coagulation due to Brownian collisions
80 
81  //- Rate for coagulation due to Brownian collisions
82  volScalarField BrownianRate_;
83 
84  //- Model for coagulation due to ballistic collisions
86 
87  //- Rate for coagulation due to ballistic collisions
88  volScalarField ballisticRate_;
89 
90 
91 public:
92 
93  //- Runtime type information
94  TypeName("DahnekeInterpolation");
95 
96  // Constructor
97 
99  (
100  const populationBalanceModel& popBal,
101  const dictionary& dict
102  );
103 
104 
105  //- Destructor
106  virtual ~DahnekeInterpolation()
107  {}
108 
109 
110  // Member Functions
111 
112  //- Precompute diameter independent expressions
113  virtual void precompute();
114 
115  //- Add to coalescenceRate
116  virtual void addToCoalescenceRate
117  (
118  volScalarField& coalescenceRate,
119  const label i,
120  const label j
121  );
122 };
123 
124 
125 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
126 
127 } // End namespace coalescenceModels
128 } // End namespace diameterModels
129 } // End namespace Foam
130 
131 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
132 
133 #endif
134 
135 // ************************************************************************* //
Generic GeometricField class.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Base class for coalescence models.
Interpolation formula of Dahneke (1983) as presented by Otto et al. (1999). Utilises collisional diam...
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
virtual void precompute()
Precompute diameter independent expressions.
DahnekeInterpolation(const populationBalanceModel &popBal, const dictionary &dict)
TypeName("DahnekeInterpolation")
Runtime type information.
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Namespace for OpenFOAM.
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
dictionary dict