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-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 "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 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  reader_
62  (
63  TableReader<scalar>::New(word::null, {defaultUnits, unitAny}, dict)
64  )
65 {
67  reader_->read({defaultUnits, unitAny}, dict, "distribution");
68 
69  // Checks
70  forAll(values, i)
71  {
72  if (i && values[i].first() - values[i - 1].first() < 0)
73  {
75  << typeName << ": The probability density function is not "
76  << "in order" << abort(FatalIOError);
77  }
78 
79  if (values[i].second() < 0)
80  {
82  << typeName << ": The probability density function is not "
83  << "positive" << 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  // Copy the PDF. Scale if q != 0.
95  PDF_.resize(values.size());
96  forAll(values, i)
97  {
98  PDF_[i] = integerPow(x_[i], q())*values[i].second();
99  }
100 
101  // Compute the CDF
102  CDF_.resize(values.size());
103  CDF_ = unintegrable::integrate(x_, PDF_);
104 
105  // Normalise the PDF and the CDF
106  PDF_ /= CDF_.last();
107  CDF_ /= CDF_.last();
108 
109  // More checks
110  validateBounds(dict);
111  if (q() != 0) validatePositive(dict);
112 }
113 
114 
116 (
117  const tabulatedDensity& d,
118  const label sampleQ
119 )
120 :
122  reader_(d.reader_, false),
123  x_(d.x_),
124  PDF_(d.PDF_),
125  CDF_(d.CDF_)
126 {
127  if (q() == d.q()) return;
128 
129  // Scale the PDF
130  PDF_ = integerPow(x_, q() - d.q())*d.PDF_;
131 
132  // Compute the CDF
133  CDF_ = unintegrable::integrate(x_, PDF_);
134 
135  // Normalise the PDF and the CDF
136  PDF_ /= CDF_.last();
137  CDF_ /= CDF_.last();
138 }
139 
140 
141 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
142 
144 {}
145 
146 
147 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
148 
150 {
151  return unintegrable::sample(x_, PDF_, CDF_, rndGen_.sample01<scalar>());
152 }
153 
154 
156 {
157  return x_.first();
158 }
159 
160 
162 {
163  return x_.last();
164 }
165 
166 
168 {
169  return unintegrable::integrateX(x_, PDF_)->last();
170 }
171 
172 
175 (
176  const scalarField& x,
177  const label e,
178  const bool
179 ) const
180 {
181  return unintegrable::interpolateIntegrateXPow(x_, e, PDF_, x);
182 }
183 
184 
186 (
187  Ostream& os,
188  const unitConversion& units
189 ) const
190 {
192 
193  // Recover the PDF
194  scalarField PDF(integerPow(x_, -q())*PDF_);
195 
196  // Normalise the PDF
197  PDF /= unintegrable::integrate(x_, PDF)->last();
198 
199  // Construct and write the values
200  List<Tuple2<scalar, scalar>> values(PDF_.size());
201  forAll(values, i)
202  {
203  values[i].first() = x_[i];
204  values[i].second() = PDF[i];
205  }
206  reader_->write(os, {units, unitAny}, values, "distribution");
207 }
208 
209 
212 {
213  const scalar x0 = min(), x1 = max(), d = 0.1*(x1 - x0);
214 
215  tmp<scalarField> tResult(new scalarField(x_.size() + 4));
216  scalarField& result = tResult.ref();
217 
218  result[0] = Foam::max(x0 - d, q() < 0 ? x0/2 : rootVSmall);
219  result[1] = x0;
220 
221  SubField<scalar>(result, x_.size(), 2) = x_;
222 
223  result[x_.size() + 2] = x1;
224  result[x_.size() + 3] = x1 + d;
225 
226  return tResult;
227 }
228 
229 
232 {
233  tmp<scalarField> tResult(new scalarField(PDF_.size() + 4, 0));
234  scalarField& result = tResult.ref();
235 
236  SubField<scalar>(result, PDF_.size(), 2) = PDF_;
237 
238  return tResult;
239 }
240 
241 
242 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
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
Pre-declare related SubField type.
Definition: SubField.H:63
Base class to read table data for tables.
Definition: TableReader.H:51
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 keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Base class for statistical distributions.
Definition: distribution.H:76
label q() const
Return the effective distribution size exponent.
Definition: distribution.H:107
Distribution in which the probability density function is given as a table of values.
virtual tmp< scalarField > integralPDFxPow(const scalarField &x, const label e, const bool consistent=false) const
Return the integral of the PDF multiplied by an integer power of x.
virtual scalar min() const
Return the minimum value.
virtual scalar sample() const
Sample the distribution.
virtual tmp< scalarField > plotPDF(const scalarField &x) const
Return values to plot the probability density function.
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.
virtual tmp< scalarField > plotX(const label n) const
Return coordinates to plot across the range of the distribution.
tabulatedDensity(const unitConversion &units, const dictionary &dict, const label sampleQ, randomGenerator &&rndGen)
Construct from a dictionary.
virtual scalar sample() const
Sample the distribution.
Definition: unintegrable.C:416
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
static tmp< scalarField > interpolateIntegrateXPow(const scalarField &xStar, const label e, const scalarField &yStar, const scalarField &x)
Integrate the values x^e*y with respect to the coordinates x,.
Definition: unintegrable.C:79
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:197
Unit conversion structure. Contains the associated dimensions and the multiplier with which to conver...
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
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.
const doubleScalar e
Definition: doubleScalar.H:106
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
const unitConversion unitAny
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.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
dictionary dict
randomGenerator rndGen(653213)