multiFixedValue.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 "multiFixedValue.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace distributions
34 {
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const unitConversion& defaultUnits,
46  const dictionary& dict,
47  const label sampleQ,
49 )
50 :
52  (
53  typeName,
54  defaultUnits,
55  dict,
56  sampleQ,
57  std::move(rndGen)
58  ),
59  reader_
60  (
61  TableReader<scalar>::New(word::null, {defaultUnits, unitAny}, dict)
62  )
63 {
65  reader_->read({defaultUnits, unitAny}, dict, "values");
66 
67  // Sort
69  (
70  values,
72  {
73  return a.first() < b.first();
74  }
75  );
76 
77  // Checks
78  for (label i = 1; i < values.size(); ++ i)
79  {
80  if (values[i].second() < 0)
81  {
83  << typeName << ": The probabilities are not all positive "
84  << abort(FatalIOError);
85  }
86  }
87 
88  // Copy the coordinates
89  x_.resize(values.size());
90  forAll(values, i)
91  {
92  x_[i] = values[i].first();
93  }
94 
95  // Copy the probabilities. Scale if q != 0.
96  P_.resize(values.size());
97  forAll(values, i)
98  {
99  P_[i] = integerPow(x_[i], q())*values[i].second();
100  }
101 
102  // Compute the cumulative sum of the probabilities
103  sumP_.resize(values.size() + 1);
104  sumP_[0] = 0;
105  forAll(values, i)
106  {
107  sumP_[i + 1] = sumP_[i] + P_[i];
108  }
109 
110  // Normalise
111  P_ /= sumP_.last();
112  sumP_ /= sumP_.last();
113 }
114 
115 
117 (
118  const multiFixedValue& d,
119  const label sampleQ
120 )
121 :
123  reader_(d.reader_, false),
124  x_(d.x_),
125  P_(d.P_),
126  sumP_(d.sumP_)
127 {
128  if (q() == d.q()) return;
129 
130  // Scale the probabilities
131  P_ = integerPow(x_, q() - d.q())*d.P_;
132 
133  // Compute the cumulative sum of the probabilities
134  sumP_.resize(x_.size() + 1);
135  sumP_[0] = 0;
136  forAll(x_, i)
137  {
138  sumP_[i + 1] = sumP_[i] + P_[i];
139  }
140 
141  // Normalise
142  P_ /= sumP_.last();
143  sumP_ /= sumP_.last();
144 }
145 
146 
147 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
148 
150 {}
151 
152 
153 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
154 
156 {
157  const scalar S = rndGen_.sample01<scalar>();
158 
159  label i = 0;
160  for (; i < sumP_.size() - 1 && sumP_[i + 1] < S; ++ i);
161 
162  return x_[i];
163 }
164 
165 
167 {
168  return x_.first();
169 }
170 
171 
173 {
174  return x_.last();
175 }
176 
177 
179 {
180  return sum(x_*P_);
181 }
182 
183 
186 (
187  const scalarField& x,
188  const label e,
189  const bool
190 ) const
191 {
192  tmp<scalarField> tResult(new scalarField(x.size()));
193  scalarField& result = tResult.ref();
194 
195  label i = 0;
196 
197  while (i < x.size() && x[i] < x_[0])
198  {
199  result[i] = 0;
200  i ++;
201  }
202 
203  scalar integral_PDFxPowE_0_j = P_[0]*integerPow(x_[0], e);
204 
205  for (label j = 0; j < x_.size() - 1; ++ j)
206  {
207  while (i < x.size() && x[i] < x_[j + 1])
208  {
209  result[i] = integral_PDFxPowE_0_j;
210  i ++;
211  }
212 
213  integral_PDFxPowE_0_j += P_[j + 1]*integerPow(x_[j + 1], e);
214  }
215 
216  while (i < x.size())
217  {
218  result[i] = integral_PDFxPowE_0_j;
219  i ++;
220  }
221 
222  return tResult;
223 }
224 
225 
227 (
228  Ostream& os,
229  const unitConversion& units
230 ) const
231 {
233 
234  // Recover the probabilities
235  scalarField P(integerPow(x_, -q())*P_);
236 
237  // Normalise the probabilities
238  P /= sum(P);
239 
240  // Construct and write the values
241  List<Tuple2<scalar, scalar>> values(P_.size());
242  forAll(values, i)
243  {
244  values[i].first() = x_[i];
245  values[i].second() = P[i];
246  }
247  reader_->write(os, {units, unitAny}, values, "values");
248 }
249 
250 
253 {
254  tmp<scalarField> tResult(new scalarField(3*x_.size()));
255  scalarField& result = tResult.ref();
256 
257  forAll(x_, i)
258  {
259  result[3*i] = result[3*i + 1] = result[3*i + 2] = x_[i];
260  }
261 
262  return tResult;
263 }
264 
265 
268 {
269  tmp<scalarField> tResult(new scalarField(3*x_.size()));
270  scalarField& result = tResult.ref();
271 
272  forAll(x_, i)
273  {
274  result[3*i] = result[3*i + 2] = 0;
275  result[3*i + 1] = P_[i]/rootVSmall;
276  }
277 
278  return tResult;
279 }
280 
281 
282 // ************************************************************************* //
#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
Base class to read table data for tables.
Definition: TableReader.H:51
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
const Type1 & first() const
Return first.
Definition: Tuple2.H:119
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 which comprises a list of fixed values with given probabilities. The probabilities are n...
multiFixedValue(const unitConversion &units, const dictionary &dict, const label sampleQ, randomGenerator &&rndGen)
Construct from a dictionary.
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.
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
volScalarField & b
Definition: createFields.H:27
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
defineTypeNameAndDebug(exponential, 0)
addToRunTimeSelectionTable(distribution, exponential, dictionary)
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
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.
const HashTable< unitConversion > & units()
Get the table of unit conversions.
void sort(UList< T > &)
Definition: UList.C:115
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh, PrimitiveField > &df)
IOerror FatalIOError
dictionary dict
randomGenerator rndGen(653213)