massRosinRammler.C
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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "massRosinRammler.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace distributionModels
36 {
38  addToRunTimeSelectionTable(distributionModel, massRosinRammler, dictionary);
39 }
40 }
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const dictionary& dict,
47  Random& rndGen
48 )
49 :
50  distributionModel(typeName, dict, rndGen),
51  minValue_(distributionModelDict_.get<scalar>("minValue")),
52  maxValue_(distributionModelDict_.get<scalar>("maxValue")),
53  d_(distributionModelDict_.get<scalar>("d")),
54  n_(distributionModelDict_.get<scalar>("n"))
55 {
56  check();
57 }
58 
59 
61 (
62  const massRosinRammler& p
63 )
64 :
66  minValue_(p.minValue_),
67  maxValue_(p.maxValue_),
68  d_(p.d_),
69  n_(p.n_)
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
74 
76 {}
77 
78 
79 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
80 
82 {
83  scalar d;
84 
85  // Re-sample if the calculated d is out of the physical range
86  do
87  {
88  const scalar a = 3/n_ + 1;
89  const scalar P = rndGen_.sample01<scalar>();
90  const scalar x = invIncGamma(a, P);
91  d = d_*pow(x, 1/n_);
92  } while (d < minValue_ || d > maxValue_);
93 
94  return d;
95 }
96 
97 
99 {
100  return minValue_;
101 }
102 
103 
105 {
106  return maxValue_;
107 }
108 
109 
111 {
112  return d_;
113 }
114 
115 
116 // ************************************************************************* //
Foam::Random
Random number generator.
Definition: Random.H:59
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::distributionModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(distributionModel, binned, dictionary)
Foam::distributionModels::massRosinRammler::meanValue
virtual scalar meanValue() const
Return the mean value.
Definition: massRosinRammler.C:110
Foam::invIncGamma
scalar invIncGamma(const scalar a, const scalar P)
Inverse normalized incomplete gamma function.
Definition: invIncGamma.C:108
Foam::distributionModels::massRosinRammler::~massRosinRammler
virtual ~massRosinRammler()
Destructor.
Definition: massRosinRammler.C:75
Foam::distributionModel
A library of runtime-selectable distribution models.
Definition: distributionModel.H:70
Foam::distributionModels::massRosinRammler
Definition: massRosinRammler.H:68
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
massRosinRammler
Mass-based Rosin-Rammler distributionModel.
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::distributionModels::massRosinRammler::sample
virtual scalar sample() const
Sample the distributionModel.
Definition: massRosinRammler.C:81
Foam::distributionModels::massRosinRammler::maxValue
virtual scalar maxValue() const
Return the maximum value.
Definition: massRosinRammler.C:104
Foam::distributionModels::massRosinRammler::massRosinRammler
massRosinRammler(const dictionary &dict, Random &rndGen)
Construct from components.
Definition: massRosinRammler.C:45
Foam::distributionModels::defineTypeNameAndDebug
defineTypeNameAndDebug(binned, 0)
massRosinRammler.H
Foam::distributionModels::massRosinRammler::minValue
virtual scalar minValue() const
Return the minimum value.
Definition: massRosinRammler.C:98
x
x
Definition: LISASMDCalcMethod2.H:52
rndGen
Random rndGen
Definition: createFields.H:23