tabulatedDensity.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 "tabulatedDensity.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  forAll(values, i)
64  {
65  if (i && values[i].first() - values[i - 1].first() < 0)
66  {
68  << typeName << ": The probability density function is not "
69  << "in order" << abort(FatalIOError);
70  }
71 
72  if (values[i].second() < 0)
73  {
75  << typeName << ": The probability density function is not "
76  << "positive" << abort(FatalIOError);
77  }
78  }
79 
80  // Copy the coordinates
81  x_.resize(values.size());
82  forAll(values, i)
83  {
84  x_[i] = values[i].first();
85  }
86 
87  // Copy the PDF. Scale if q != 0.
88  PDF_.resize(values.size());
89  forAll(values, i)
90  {
91  PDF_[i] = integerPow(x_[i], q())*values[i].second();
92  }
93 
94  // Compute the CDF
95  CDF_.resize(values.size());
96  CDF_ = unintegrable::integrate(x_, PDF_);
97 
98  // Normalise the PDF and the CDF
99  PDF_ /= CDF_.last();
100  CDF_ /= CDF_.last();
101 
102  // More checks
104  if (q() != 0) validatePositive(dict);
105  report();
106 }
107 
108 
110 (
111  const tabulatedDensity& d,
112  const label sampleQ
113 )
114 :
116  x_(d.x_),
117  PDF_(d.PDF_),
118  CDF_(d.CDF_)
119 {
120  if (q() == d.q()) return;
121 
122  // Scale the PDF
123  PDF_ = integerPow(x_, q() - d.q())*d.PDF_;
124 
125  // Compute the CDF
126  CDF_ = unintegrable::integrate(x_, PDF_);
127 
128  // Normalise the PDF and the CDF
129  PDF_ /= CDF_.last();
130  CDF_ /= CDF_.last();
131 }
132 
133 
134 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
135 
137 {}
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
143 {
144  return unintegrable::sample(x_, PDF_, CDF_, rndGen_.sample01<scalar>());
145 }
146 
147 
149 {
150  return x_.first();
151 }
152 
153 
155 {
156  return x_.last();
157 }
158 
159 
161 {
162  return unintegrable::integrateX(x_, PDF_)->last();
163 }
164 
165 
167 (
168  const label
169 ) const
170 {
171  const scalar x0 = min(), x1 = max(), d = 0.1*(x1 - x0);
172 
173  tmp<scalarField> tResult(new scalarField(x_.size() + 4));
174  scalarField& result = tResult.ref();
175 
176  result[0] = Foam::max(x0 - d, q() < 0 ? x0/2 : rootVSmall);
177  result[1] = x0;
178 
179  SubField<scalar>(result, x_.size(), 2) = x_;
180 
181  result[x_.size() + 2] = x1;
182  result[x_.size() + 3] = x1 + d;
183 
184  return tResult;
185 }
186 
187 
190 {
191  tmp<scalarField> tResult(new scalarField(PDF_.size() + 4, 0));
192  scalarField& result = tResult.ref();
193 
194  SubField<scalar>(result, PDF_.size(), 2) = PDF_;
195 
196  return tResult;
197 }
198 
199 
200 // ************************************************************************* //
#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
Pre-declare related SubField type.
Definition: SubField.H:63
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 probability density 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.
tabulatedDensity(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
static tmp< scalarField > integrateX(const scalarField &x, const scalarField &y)
Integrate the values x*y with respect to the coordinates x.
Definition: unintegrable.C:54
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))