RosinRammler.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) 2011-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 \*---------------------------------------------------------------------------*/
25 
26 #include "RosinRammler.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace distributions
34 {
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42 
44 (
45  const label q,
46  const scalarField& x
47 ) const
48 {
49  const scalarField xByDPowNm1(pow(x/d_, n_ - 1));
50 
51  return integerPow(x, q)*(n_/d_)*xByDPowNm1*exp(- xByDPowNm1*x/d_);
52 }
53 
54 
56 (
57  const label q,
58  const scalarField& x
59 ) const
60 {
61  if (q == 0)
62  {
63  return - exp(- pow(x/d_, n_));
64  }
65  else
66  {
68  }
69 }
70 
71 
72 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
73 
75 (
76  const dictionary& dict,
77  Random& rndGen,
78  const label sampleQ
79 )
80 :
82  (
83  typeName,
84  dict,
85  rndGen,
86  sampleQ
87  ),
88  min_(dict.lookupBackwardsCompatible<scalar>({"min", "minValue"})),
89  max_(dict.lookupBackwardsCompatible<scalar>({"max", "maxValue"})),
90  d_(dict.lookup<scalar>("d")),
91  n_(dict.lookup<scalar>("n"))
92 {
93  validateBounds(dict);
94  validatePositive(dict);
95  mean();
96  report();
97 }
98 
99 
101 (
102  const RosinRammler& d,
103  const label sampleQ
104 )
105 :
107  min_(d.min_),
108  max_(d.max_),
109  d_(d.d_),
110  n_(d.n_)
111 {}
112 
113 
114 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
115 
117 {}
118 
119 
120 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
121 
123 {
124  if (q() == 0)
125  {
126  const scalar s = rndGen_.sample01<scalar>();
127  const Pair<scalar>& Phi01 = this->Phi01();
128  const scalar PhiS = (1 - s)*Phi01[0] + s*Phi01[1];
129  return d_*pow(- log(- PhiS), 1/n_);
130  }
131  else
132  {
134  }
135 }
136 
137 
139 {
140  return min_;
141 }
142 
143 
145 {
146  return max_;
147 }
148 
149 
152 {
154  tx.ref()[0] = Foam::max(tx.ref()[0], q() < 0 ? min_/2 : rootVSmall);
155  return tx;
156 }
157 
158 
159 // ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
Random number generator.
Definition: Random.H:58
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Accumulating histogram of values. Specified bin resolution automatic generation of bins.
Definition: distribution.H:61
virtual tmp< scalarField > x(const label n) const
Return coordinates to plot across the range of the distribution.
Definition: distribution.C:140
Rosin-Rammler distribution, scaled so that it spans between a specified minimum and maximum value,...
Definition: RosinRammler.H:77
virtual tmp< scalarField > Phi(const label q, const scalarField &x) const
Return values of the un-normalised CDF for the given size exponent.
Definition: RosinRammler.C:56
virtual scalar min() const
Return the minimum value.
Definition: RosinRammler.C:138
virtual scalar sample() const
Sample the distribution.
Definition: RosinRammler.C:122
virtual tmp< scalarField > x(const label n) const
Return coordinates to plot across the range of the distribution.
Definition: RosinRammler.C:151
virtual ~RosinRammler()
Destructor.
Definition: RosinRammler.C:116
RosinRammler(const dictionary &dict, Random &rndGen, const label sampleQ)
Construct from a dictionary.
Definition: RosinRammler.C:75
virtual scalar max() const
Return the maximum value.
Definition: RosinRammler.C:144
virtual tmp< scalarField > phi(const label q, const scalarField &x) const
Return values of the un-normalised PDF for the given size exponent.
Definition: RosinRammler.C:44
virtual tmp< scalarField > Phi(const label q, const scalarField &x) const
Return values of the un-normalised CDF for the given size exponent.
Definition: unintegrable.C:253
scalar sample() const
Sample the distribution.
Definition: unintegrable.C:349
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
defineTypeNameAndDebug(exponential, 0)
addToRunTimeSelectionTable(distribution, exponential, dictionary)
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
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
scalar integerPow(const scalar x, const label e)
Compute the power of the number x to the integer e.
Definition: scalarI.H:30
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dictionary dict
Random rndGen(label(0))