general.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-2021 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 "general.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace distributionModels
34 {
35  defineTypeNameAndDebug(general, 0);
36  addToRunTimeSelectionTable(distributionModel, general, dictionary);
37 }
38 }
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
43 (
44  const dictionary& dict,
45  Random& rndGen
46 )
47 :
48  distributionModel(typeName, dict, rndGen),
49  xy_(distributionModelDict_.lookup("distribution")),
50  nEntries_(xy_.size()),
51  minValue_(xy_[0][0]),
52  maxValue_(xy_[nEntries_-1][0]),
53  meanValue_(0.0),
54  integral_(nEntries_),
55  cumulative_(distributionModelDict_.lookupOrDefault("cumulative", false))
56 {
57  check();
58 
59  // Additional sanity checks
60  if (cumulative_ && xy_[0][1] != 0)
61  {
63  << type() << "distribution: "
64  << "The cumulative distribution must start from zero."
65  << abort(FatalError);
66  }
67 
68  for (label i=0; i<nEntries_; i++)
69  {
70  if (i > 0 && xy_[i][0] <= xy_[i-1][0])
71  {
73  << type() << "distribution: "
74  << "The points must be specified in ascending order."
75  << abort(FatalError);
76  }
77 
78  if (xy_[i][1] < 0)
79  {
81  << type() << "distribution: "
82  << "The distribution can't be negative."
83  << abort(FatalError);
84  }
85 
86  if (i > 0 && cumulative_ && xy_[i][1] < xy_[i-1][1])
87  {
89  << type() << "distribution: "
90  << "Cumulative distribution must be non-decreasing."
91  << abort(FatalError);
92  }
93  }
94 
95  // Fill out the integral table (x, P(x<=0)) and calculate mean
96  // For density function: P(x<=0) = int f(x) and mean = int x*f(x)
97  // For cumulative function: mean = int 1-P(x<=0) = maxValue_ - int P(x<=0)
98  integral_[0] = 0.0;
99  for (label i=1; i<nEntries_; i++)
100  {
101  // Integrating k*x+d
102  scalar k = (xy_[i][1] - xy_[i-1][1])/(xy_[i][0] - xy_[i-1][0]);
103  scalar d = xy_[i-1][1] - k*xy_[i-1][0];
104 
105  scalar y1 = xy_[i][0]*(0.5*k*xy_[i][0] + d);
106  scalar y0 = xy_[i-1][0]*(0.5*k*xy_[i-1][0] + d);
107  scalar area = y1 - y0;
108 
109  if (cumulative_)
110  {
111  integral_[i] = xy_[i][1];
112  meanValue_ += area;
113  }
114  else
115  {
116  integral_[i] = area + integral_[i-1];
117 
118  y1 = sqr(xy_[i][0])*(1.0/3.0*k*xy_[i][0] + 0.5*d);
119  y0 = sqr(xy_[i-1][0])*(1.0/3.0*k*xy_[i-1][0] + 0.5*d);
120  meanValue_ += y1 - y0;
121  }
122  }
123 
124  // normalise the distribution
125  scalar sumArea = integral_.last();
126 
127  for (label i=0; i<nEntries_; i++)
128  {
129  xy_[i][1] /= sumArea;
130  integral_[i] /= sumArea;
131  }
132 
133  meanValue_ /= sumArea;
134  meanValue_ = cumulative_ ? (maxValue_ - meanValue_) : meanValue_;
135 
136  info();
137 }
138 
139 
141 :
143  xy_(p.xy_),
144  nEntries_(p.nEntries_),
145  minValue_(p.minValue_),
146  maxValue_(p.maxValue_),
147  meanValue_(p.meanValue_),
148  integral_(p.integral_),
149  cumulative_(p.cumulative_)
150 {}
151 
152 
153 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
154 
156 {}
157 
158 
159 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
160 
162 {
163  scalar y = rndGen_.sample01<scalar>();
164 
165  // find the interval where y is in the table
166  label n=1;
167  while (integral_[n] <= y)
168  {
169  n++;
170  }
171 
172  scalar k = (xy_[n][1] - xy_[n-1][1])/(xy_[n][0] - xy_[n-1][0]);
173  scalar d = xy_[n-1][1] - k*xy_[n-1][0];
174 
175  if (cumulative_)
176  {
177  return (y - d)/k;
178  }
179 
180  scalar alpha = y + xy_[n-1][0]*(0.5*k*xy_[n-1][0] + d) - integral_[n-1];
181  scalar x = 0.0;
182 
183  // if k is small it is a linear equation, otherwise it is of second order
184  if (mag(k) > small)
185  {
186  scalar p = 2.0*d/k;
187  scalar q = -2.0*alpha/k;
188  scalar sqrtEr = sqrt(0.25*p*p - q);
189 
190  scalar x1 = -0.5*p + sqrtEr;
191  scalar x2 = -0.5*p - sqrtEr;
192  if ((x1 >= xy_[n-1][0]) && (x1 <= xy_[n][0]))
193  {
194  x = x1;
195  }
196  else
197  {
198  x = x2;
199  }
200  }
201  else
202  {
203  x = alpha/d;
204  }
205 
206  return x;
207 }
208 
209 
211 {
212  return minValue_;
213 }
214 
215 
217 {
218  return maxValue_;
219 }
220 
221 
223 {
224  return meanValue_;
225 }
226 
227 
228 // ************************************************************************* //
virtual scalar sample() const
Sample the distributionModel.
Definition: general.C:161
virtual scalar minValue() const
Return the minimum value.
Definition: general.C:210
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
Random & rndGen_
Reference to the random number generator.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual ~general()
Destructor.
Definition: general.C:155
dimensionedScalar y0(const dimensionedScalar &ds)
label k
Boltzmann constant.
Macros for easy insertion into run-time selection tables.
scalar y
Type sample01()
Advance the state and return a sample of a given type from a.
Definition: RandomI.H:70
virtual scalar meanValue() const
Return the mean value.
Definition: general.C:222
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensionedScalar y1(const dimensionedScalar &ds)
Random number generator.
Definition: Random.H:57
defineTypeNameAndDebug(exponential, 0)
general(const dictionary &dict, Random &rndGen)
Construct from components.
Definition: general.C:43
A library of runtime-selectable distribution models.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
A general distribution model where the distribution is specified as (point, value) pairs...
Definition: general.H:58
virtual scalar maxValue() const
Return the maximum value.
Definition: general.C:216
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
addToRunTimeSelectionTable(distributionModel, exponential, dictionary)
volScalarField & p
Namespace for OpenFOAM.