massRosinRammler.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2016 OpenFOAM Foundation
9 Copyright (C) 2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::distributionModels::massRosinRammler
29
30Description
31 Particle-size distribution model wherein random samples are drawn
32 from the two-parameter Rosin-Rammler (Weibull) probability density
33 function corrected to take into account varying number of particles
34 per parcels for fixed-mass parcels.
35
36 "There is a factor of \f$x^3\f$ difference between the size distribution
37 probability density function of individual particles and modeled parcels
38 of particles, \f$ f_{parcels}(D) = x^3 f_{particles}(D) \f$, because the
39 submodel parameter, \f$ PPP \f$ (number of particles per parcel) is based
40 on a fixed mass per parcel which weights the droplet distribution by
41 a factor proportional to \f$ 1/x^3 \f$ (i.e. \f$ PPP = \dot{m}_{total}
42 \Delta_{t_{inj(per-parcel)}} / {m}_{particle} \f$)." (YHD:p. 1374-1375).
43
44 \c massRosinRammler should be preferred over \c RosinRammler
45 when \c parcelBasisType is based on \c mass:
46 \verbatim
47 parcelBasisType mass;
48 \endverbatim
49
50 The doubly-truncated mass-corrected
51 Rosin-Rammler probability density function:
52
53 \f[
54 f(x; \lambda, n, A, B) =
55 x^3
56 \frac{n}{\lambda}
57 \left( \frac{x}{\lambda} \right)^{n-1}
58 \frac{
59 \exp\{ -(\frac{x}{\lambda})^n \}
60 }
61 {
62 \exp\{ -(\frac{A}{\lambda})^n \}
63 - \exp\{ -(\frac{B}{\lambda})^n \}
64 }
65 \f]
66 where
67
68 \vartable
69 f(x; \lambda, n, A, B) | Rosin-Rammler probability density function
70 \lambda | Scale parameter
71 n | Shape parameter
72 x | Sample
73 A | Minimum of the distribution
74 B | Maximum of the distribution
75 \endvartable
76
77 Constraints:
78 - \f$ \lambda > 0 \f$
79 - \f$ n > 0 \f$
80
81 Random samples are generated by the inverse transform sampling technique
82 by using the quantile function of the doubly-truncated two-parameter
83 mass-corrected Rosin-Rammler (Weibull) probability density function:
84
85 \f[
86 d = \lambda \, Q(a, u)^{1/n}
87 \f]
88 with
89
90 \f[
91 a = \frac{3}{n} + 1
92 \f]
93 where \f$ Q(z_1, z_2) \f$ is the inverse of regularised lower incomplete
94 gamma function, and \f$ u \f$ is a sample drawn from the uniform
95 probability density function on the interval \f$ (x, y) \f$:
96
97 \f[
98 x = \gamma \left( a, \frac{A}{\lambda}^n \right)
99 \f]
100
101 \f[
102 y = \gamma \left( a, \frac{B}{\lambda}^n \right)
103 \f]
104 where \f$ \gamma(z_1, z_2) \f$ is the lower incomplete gamma function.
105
106 Reference:
107 \verbatim
108 Standard model (tag:YHD):
109 Yoon, S. S., Hewson, J. C., DesJardin, P. E.,
110 Glaze, D. J., Black, A. R., & Skaggs, R. R. (2004).
111 Numerical modeling and experimental measurements of a high speed
112 solid-cone water spray for use in fire suppression applications.
113 International Journal of Multiphase Flow, 30(11), 1369-1388.
114 DOI:10.1016/j.ijmultiphaseflow.2004.07.006
115 \endverbatim
116
117Usage
118 Minimal example by using \c constant/<CloudProperties>:
119 \verbatim
120 subModels
121 {
122 injectionModels
123 {
124 <name>
125 {
126 ...
127 parcelBasisType mass;
128 ...
129
130 sizeDistribution
131 {
132 type massRosinRammler;
133 massRosinRammlerDistribution
134 {
135 lambda <scaleParameterValue>;
136 n <shapeParameterValue>;
137 minValue <minValue>;
138 maxValue <maxValue>;
139 }
140 }
141 }
142 }
143 }
144 \endverbatim
145
146 where the entries mean:
147 \table
148 Property | Description | Type | Reqd | Deflt
149 type | Type name: massRosinRammler | word | yes | -
150 massRosinRammlerDistribution | Distribution settings | dict | yes | -
151 lambda | Scale parameter | scalar | yes | -
152 n | Shape parameter | scalar | yes | -
153 minValue | Minimum of the distribution | scalar | yes | -
154 maxValue | Maximum of the distribution | scalar | yes | -
155 \endtable
156
157Note
158 - In the current context, \c lambda (or \c d) is called "characteristic
159 droplet size", and \c n "empirical dimensionless constant to specify
160 the distribution width, sometimes referred to as the dispersion
161 coefficient." (YHD:p. 1374).
162
163SourceFiles
164 massRosinRammler.C
165
166\*---------------------------------------------------------------------------*/
167
168#ifndef massRosinRammler_H
169#define massRosinRammler_H
170
171#include "distributionModel.H"
172
173// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
174
175namespace Foam
176{
177namespace distributionModels
178{
179
180/*---------------------------------------------------------------------------*\
181 Class massRosinRammler Declaration
182\*---------------------------------------------------------------------------*/
183
184class massRosinRammler
185:
186 public distributionModel
187{
188 // Private Data
189
190 //- Scale parameter
191 scalar lambda_;
192
193 //- Shape parameter
194 scalar n_;
195
196
197public:
198
199 //- Runtime type information
200 TypeName("massRosinRammler");
201
202
203 // Constructors
204
205 //- Construct from components
206 massRosinRammler(const dictionary& dict, Random& rndGen);
207
208 //- Copy construct
209 massRosinRammler(const massRosinRammler& p);
210
211 //- Construct and return a clone
212 virtual autoPtr<distributionModel> clone() const
213 {
214 return autoPtr<distributionModel>(new massRosinRammler(*this));
215 }
216
217 //- No copy assignment
218 void operator=(const massRosinRammler&) = delete;
219
220
221 //- Destructor
222 virtual ~massRosinRammler() = default;
223
224
225 // Member Functions
226
227 //- Sample the distribution
228 virtual scalar sample() const;
229
230 //- Return the theoretical mean of the distribution
231 virtual scalar meanValue() const;
232};
233
234
235// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236
237} // End namespace distributionModels
238} // End namespace Foam
239
240// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241
242#endif
243
244// ************************************************************************* //
Random number generator.
Definition: Random.H:60
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A library of runtime-selectable doubly-truncated probability distribution models. Returns random samp...
Particle-size distribution model wherein random samples are drawn from the two-parameter Rosin-Rammle...
virtual scalar meanValue() const
Return the theoretical mean of the distribution.
massRosinRammler(const dictionary &dict, Random &rndGen)
Construct from components.
virtual autoPtr< distributionModel > clone() const
Construct and return a clone.
void operator=(const massRosinRammler &)=delete
No copy assignment.
virtual scalar sample() const
Sample the distribution.
TypeName("massRosinRammler")
Runtime type information.
virtual ~massRosinRammler()=default
Destructor.
volScalarField & p
Namespace for OpenFOAM.
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73
Random rndGen
Definition: createFields.H:23