tabulatedCumulative.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 "tabulatedCumulative.H"
27 #include "unintegrable.H"
28 #include "SubField.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace distributions
36 {
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 (
47  const dictionary& dict,
48  Random& rndGen,
49  const label sampleQ
50 )
51 :
53  (
54  typeName,
55  dict,
56  rndGen,
57  sampleQ
58  )
59 {
60  List<Pair<scalar>> values(dict.lookup("distribution"));
61 
62  // Checks
63  if (values.first().second() != 0)
64  {
66  << typeName << ": The cumulative distribution function does not "
67  << "start at zero" << abort(FatalIOError);
68  }
69 
70  for (label i = 1; i < values.size(); ++ i)
71  {
72  if (values[i].first() - values[i - 1].first() < 0)
73  {
75  << typeName << ": The cumulative distribution function is not "
76  << "in order" << abort(FatalIOError);
77  }
78 
79  if (values[i].second() - values[i - 1].second() < 0)
80  {
82  << typeName << ": The cumulative distribution function is not "
83  << "monotonic" << abort(FatalIOError);
84  }
85  }
86 
87  // Copy the coordinates
88  x_.resize(values.size());
89  forAll(values, i)
90  {
91  x_[i] = values[i].first();
92  }
93 
94  // Set the CDF. Copy if q == 0. Re-integrated if q != 0.
95  CDF_.resize(values.size());
96  CDF_[0] = 0;
97  for (label i = 1; i < values.size(); ++ i)
98  {
99  CDF_[i] =
100  CDF_[i - 1]
101  + integerPow((x_[i] + x_[i - 1])/2, q())
102  *(values[i].second() - values[i - 1].second());
103  }
104  CDF_ /= CDF_.last();
105 
106  // Compute the PDF. Equals the gradient of the CDF. Associated with the
107  // intervals, so the field is one element shorter than the coordinates and
108  // the CDF.
109  PDF_.resize(values.size() - 1);
110  forAll(PDF_, i)
111  {
112  PDF_[i] = (CDF_[i + 1] - CDF_[i])/(x_[i + 1] - x_[i]);
113  }
114 
115  // More checks
117  if (q() != 0) validatePositive(dict);
118  report();
119 }
120 
121 
123 (
124  const tabulatedCumulative& d,
125  const label sampleQ
126 )
127 :
129  x_(d.x_),
130  PDF_(d.PDF_),
131  CDF_(d.CDF_)
132 {
133  // If Q is the same then a copy is sufficient
134  if (q() == d.q()) return;
135 
136  // Re-integrate the CDF
137  CDF_[0] = 0;
138  for (label i = 1; i < x_.size(); ++ i)
139  {
140  CDF_[i] =
141  CDF_[i - 1]
142  + integerPow((d.x_[i] + d.x_[i - 1])/2, q() - d.q())
143  *(d.CDF_[i] - d.CDF_[i - 1]);
144  }
145  CDF_ /= CDF_.last();
146 
147  // Re-compute the PDF
148  forAll(PDF_, i)
149  {
150  PDF_[i] = (CDF_[i + 1] - CDF_[i])/(x_[i + 1] - x_[i]);
151  }
152 }
153 
154 
155 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
156 
158 {}
159 
160 
161 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
162 
164 {
165  return unintegrable::sample(x_, CDF_, rndGen_.sample01<scalar>());
166 }
167 
168 
170 {
171  return x_.first();
172 }
173 
174 
176 {
177  return x_.last();
178 }
179 
180 
182 {
183  const scalarField x(this->x(-1));
184  return unintegrable::integrate(x, x*PDF(x))->last();
185 }
186 
187 
189 (
190  const label
191 ) const
192 {
193  const scalar x0 = min(), x1 = max(), d = 0.1*(x1 - x0);
194 
195  tmp<scalarField> tResult(new scalarField(2*x_.size() + 2));
196  scalarField& result = tResult.ref();
197 
198  result[0] = Foam::max(x0 - d, q() < 0 ? x0/2 : rootVSmall);
199 
200  forAll(x_, i)
201  {
202  result[2*i + 1] = result[2*i + 2] = x_[i];
203  }
204 
205  result[2*x_.size() + 1] = x1 + d;
206 
207  return tResult;
208 }
209 
210 
213 {
214  tmp<scalarField> tResult(new scalarField(2*PDF_.size() + 4));
215  scalarField& result = tResult.ref();
216 
217  result[0] = result[1] = 0;
218 
219  forAll(PDF_, i)
220  {
221  result[2*i + 2] = result[2*i + 3] = PDF_[i];
222  }
223 
224  result[2*PDF_.size() + 2] = result[2*PDF_.size() + 3] = 0;
225 
226  return tResult;
227 }
228 
229 
230 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void resize(const label)
Alias for setSize(const label)
Definition: ListI.H:138
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Random number generator.
Definition: Random.H:58
T & first()
Return the first element of the list.
Definition: UListI.H:114
T & last()
Return the last element of the list.
Definition: UListI.H:128
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
void report() const
Report.
Definition: distribution.C:133
virtual void validatePositive(const dictionary &dict) const
Validate that the lower bound is positive.
Definition: distribution.C:51
label q() const
Return the effective distribution size exponent.
Definition: distribution.H:106
virtual void validateBounds(const dictionary &dict) const
Validate that the bounds are monotonic.
Definition: distribution.C:39
Distribution in which the cumulative distribution function is given as a table of values.
virtual scalar min() const
Return the minimum value.
virtual scalar sample() const
Sample the distribution.
virtual tmp< scalarField > PDF(const scalarField &x) const
Return the distribution probability density function.
virtual tmp< scalarField > x(const label n) const
Return coordinates to plot across the range of the distribution.
tabulatedCumulative(const dictionary &dict, Random &rndGen, const label sampleQ)
Construct from a dictionary.
virtual scalar max() const
Return the maximum value.
virtual scalar mean() const
Return the mean value.
scalar sample() const
Sample the distribution.
Definition: unintegrable.C:349
static tmp< scalarField > integrate(const scalarField &x, const scalarField &y)
Integrate the values y with respect to the coordinates x.
Definition: unintegrable.C:34
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
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
defineTypeNameAndDebug(exponential, 0)
addToRunTimeSelectionTable(distribution, exponential, dictionary)
Namespace for OpenFOAM.
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
errorManip< error > abort(error &err)
Definition: errorManip.H:131
scalar integerPow(const scalar x, const label e)
Compute the power of the number x to the integer e.
Definition: scalarI.H:30
labelList second(const UList< labelPair > &p)
Definition: patchToPatch.C:49
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
dictionary dict
Random rndGen(label(0))