massRosinRammler.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) 2016-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 "massRosinRammler.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace distributionModels
34 {
35  defineTypeNameAndDebug(massRosinRammler, 0);
36  addToRunTimeSelectionTable(distributionModel, massRosinRammler, dictionary);
37 }
38 }
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
43 (
44  const dictionary& dict,
45  Random& rndGen
46 )
47 :
48  distributionModel(typeName, dict, rndGen),
49  minValue_(readScalar(distributionModelDict_.lookup("minValue"))),
50  maxValue_(readScalar(distributionModelDict_.lookup("maxValue"))),
51  d_(readScalar(distributionModelDict_.lookup("d"))),
52  n_(readScalar(distributionModelDict_.lookup("n")))
53 {
54  check();
55 }
56 
57 
59 (
60  const massRosinRammler& p
61 )
62 :
64  minValue_(p.minValue_),
65  maxValue_(p.maxValue_),
66  d_(p.d_),
67  n_(p.n_)
68 {}
69 
70 
71 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
72 
74 {}
75 
76 
77 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78 
80 {
81  scalar d;
82 
83  // Re-sample if the calculated d is out of the physical range
84  do
85  {
86  const scalar a = 3/n_ + 1;
87  const scalar P = rndGen_.sample01<scalar>();
88  const scalar x = invIncGamma(a, P);
89  d = d_*pow(x, 1/n_);
90  } while (d < minValue_ || d > maxValue_);
91 
92  return d;
93 }
94 
95 
97 {
98  return minValue_;
99 }
100 
101 
103 {
104  return maxValue_;
105 }
106 
107 
109 {
110  return d_;
111 }
112 
113 
114 // ************************************************************************* //
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual scalar meanValue() const
Return the mean value.
Macros for easy insertion into run-time selection tables.
massRosinRammler(const dictionary &dict, Random &rndGen)
Construct from components.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
Random number generator.
Definition: Random.H:57
scalar invIncGamma(const scalar a, const scalar P)
Inverse normalized incomplete gamma function.
Definition: invIncGamma.C:106
defineTypeNameAndDebug(exponential, 0)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
A library of runtime-selectable distribution models.
addToRunTimeSelectionTable(distributionModel, exponential, dictionary)
virtual scalar maxValue() const
Return the maximum value.
virtual scalar minValue() const
Return the minimum value.
virtual scalar sample() const
Sample the distributionModel.
Mass-based Rosin-Rammler distributionModel.
Namespace for OpenFOAM.