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-2024 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 unitConversion& defaultUnits,
48  const dictionary& dict,
49  const label sampleQ,
51 )
52 :
54  (
55  typeName,
56  defaultUnits,
57  dict,
58  sampleQ,
59  std::move(rndGen)
60  )
61 {
62  List<Tuple2<scalar, scalar>> values(dict.lookup("distribution"));
63 
64  // Checks
65  if (values.first().second() != 0)
66  {
68  << typeName << ": The cumulative distribution function does not "
69  << "start at zero" << abort(FatalIOError);
70  }
71 
72  for (label i = 1; i < values.size(); ++ i)
73  {
74  if (values[i].first() - values[i - 1].first() < 0)
75  {
77  << typeName << ": The cumulative distribution function is not "
78  << "in order" << abort(FatalIOError);
79  }
80 
81  if (values[i].second() - values[i - 1].second() < 0)
82  {
84  << typeName << ": The cumulative distribution function is not "
85  << "monotonic" << abort(FatalIOError);
86  }
87  }
88 
89  // Optionally read units
90  unitConversion units(defaultUnits);
91  units.readIfPresent("units", dict);
92 
93  // Copy the coordinates
94  x_.resize(values.size());
95  forAll(values, i)
96  {
97  x_[i] = units.toStandard(values[i].first());
98  }
99 
100  // Set the CDF. Copy if q == 0. Re-integrated if q != 0.
101  CDF_.resize(values.size());
102  CDF_[0] = 0;
103  for (label i = 1; i < values.size(); ++ i)
104  {
105  CDF_[i] =
106  CDF_[i - 1]
107  + integerPow((x_[i] + x_[i - 1])/2, q())
108  *(values[i].second() - values[i - 1].second());
109  }
110  CDF_ /= CDF_.last();
111 
112  // Compute the PDF. Equals the gradient of the CDF. Associated with the
113  // intervals, so the field is one element shorter than the coordinates and
114  // the CDF.
115  PDF_.resize(values.size() - 1);
116  forAll(PDF_, i)
117  {
118  PDF_[i] = (CDF_[i + 1] - CDF_[i])/(x_[i + 1] - x_[i]);
119  }
120 
121  // More checks
123  if (q() != 0) validatePositive(dict);
124 }
125 
126 
128 (
129  const tabulatedCumulative& d,
130  const label sampleQ
131 )
132 :
134  x_(d.x_),
135  PDF_(d.PDF_),
136  CDF_(d.CDF_)
137 {
138  // If Q is the same then a copy is sufficient
139  if (q() == d.q()) return;
140 
141  // Re-integrate the CDF
142  CDF_[0] = 0;
143  for (label i = 1; i < x_.size(); ++ i)
144  {
145  CDF_[i] =
146  CDF_[i - 1]
147  + integerPow((d.x_[i] + d.x_[i - 1])/2, q() - d.q())
148  *(d.CDF_[i] - d.CDF_[i - 1]);
149  }
150  CDF_ /= CDF_.last();
151 
152  // Re-compute the PDF
153  forAll(PDF_, i)
154  {
155  PDF_[i] = (CDF_[i + 1] - CDF_[i])/(x_[i + 1] - x_[i]);
156  }
157 }
158 
159 
160 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
161 
163 {}
164 
165 
166 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
167 
169 {
170  return unintegrable::sample(x_, CDF_, rndGen_.sample01<scalar>());
171 }
172 
173 
175 {
176  return x_.first();
177 }
178 
179 
181 {
182  return x_.last();
183 }
184 
185 
187 {
188  const scalarField x(this->x(-1));
189  return unintegrable::integrate(x, x*PDF(x))->last();
190 }
191 
192 
194 (
195  Ostream& os,
196  const unitConversion& units
197 ) const
198 {
200 
201  // Recover the CDF
202  scalarField CDF(CDF_.size());
203  CDF[0] = 0;
204  for (label i = 1; i < CDF_.size(); ++ i)
205  {
206  CDF[i] =
207  CDF[i - 1]
208  + integerPow((x_[i] + x_[i - 1])/2, -q())
209  *(CDF_[i] - CDF_[i - 1]);
210  }
211 
212  // Normalise the CDF
213  CDF /= CDF.last();
214 
215  // Construct and write the values
216  List<Tuple2<scalar, scalar>> values(CDF_.size());
217  forAll(values, i)
218  {
219  values[i].first() = units.toUser(x_[i]);
220  values[i].second() = CDF[i];
221  }
222  writeEntry(os, "distribution", values);
223 }
224 
225 
227 (
228  const label
229 ) const
230 {
231  const scalar x0 = min(), x1 = max(), d = 0.1*(x1 - x0);
232 
233  tmp<scalarField> tResult(new scalarField(2*x_.size() + 2));
234  scalarField& result = tResult.ref();
235 
236  result[0] = Foam::max(x0 - d, q() < 0 ? x0/2 : rootVSmall);
237 
238  forAll(x_, i)
239  {
240  result[2*i + 1] = result[2*i + 2] = x_[i];
241  }
242 
243  result[2*x_.size() + 1] = x1 + d;
244 
245  return tResult;
246 }
247 
248 
251 {
252  tmp<scalarField> tResult(new scalarField(2*PDF_.size() + 4));
253  scalarField& result = tResult.ref();
254 
255  result[0] = result[1] = 0;
256 
257  forAll(PDF_, i)
258  {
259  result[2*i + 2] = result[2*i + 3] = PDF_[i];
260  }
261 
262  result[2*PDF_.size() + 2] = result[2*PDF_.size() + 3] = 0;
263 
264  return tResult;
265 }
266 
267 
268 // ************************************************************************* //
#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
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
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:162
Base class for statistical distributions.
Definition: distribution.H:76
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:107
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.
tabulatedCumulative(const unitConversion &units, const dictionary &dict, const label sampleQ, randomGenerator &&rndGen)
Construct from a dictionary.
virtual tmp< scalarField > x(const label n) const
Return coordinates to plot across the range of the distribution.
virtual void write(Ostream &os, const unitConversion &units) const
Write to a stream.
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:351
static tmp< scalarField > integrate(const scalarField &x, const scalarField &y)
Integrate the values y with respect to the coordinates x.
Definition: unintegrable.C:34
Random number generator.
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
Unit conversion structure. Contains the associated dimensions and the multiplier with which to conver...
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
defineTypeNameAndDebug(exponential, 0)
addToRunTimeSelectionTable(distribution, exponential, dictionary)
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
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)
const HashTable< unitConversion > & units()
Get the table of unit conversions.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
dictionary dict
randomGenerator rndGen(653213)