general.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2019 OpenFOAM Foundation
9  Copyright (C) 2016-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "general.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace distributionModels
37 {
39  addToRunTimeSelectionTable(distributionModel, general, dictionary);
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 void Foam::distributionModels::general::initialise()
47 {
48  static scalar eps = ROOTVSMALL;
49 
50  integral_.setSize(nEntries_);
51 
52  // Fill out the integral table (x, P(x<=0)) and calculate mean
53  // For density function: P(x<=0) = int f(x) and mean = int x*f(x)
54  // For cumulative function: mean = int 1-P(x<=0) = maxValue_ - int P(x<=0)
55  integral_[0] = 0;
56 
57  for (label i = 1; i < nEntries_; ++i)
58  {
59  scalar k = (xy_[i][1] - xy_[i-1][1])/(xy_[i][0] - xy_[i-1][0] + eps);
60  scalar d = xy_[i-1][1] - k*xy_[i-1][0];
61  scalar y1 = xy_[i][0]*(0.5*k*xy_[i][0] + d);
62  scalar y0 = xy_[i-1][0]*(0.5*k*xy_[i-1][0] + d);
63  scalar area = y1 - y0;
64 
65  if (cumulative_)
66  {
67  integral_[i] = xy_[i][1];
68  meanValue_ += area;
69  }
70  else
71  {
72  integral_[i] = area + integral_[i-1];
73 
74  y1 = sqr(xy_[i][0])*(1.0/3.0*k*xy_[i][0] + 0.5*d);
75  y0 = sqr(xy_[i-1][0])*(1.0/3.0*k*xy_[i-1][0] + 0.5*d);
76  meanValue_ += y1 - y0;
77  }
78  }
79 
80  // normalize the distribution
81  const scalar sumArea = integral_.last();
82 
83  for (label i = 0; i < nEntries_; ++i)
84  {
85  xy_[i][1] /= sumArea + eps;
86  integral_[i] /= sumArea + eps;
87  }
88 
89  meanValue_ /= sumArea + eps;
90  meanValue_ = cumulative_ ? (maxValue_ - meanValue_) : meanValue_;
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
95 
97 (
98  const dictionary& dict,
99  Random& rndGen
100 )
101 :
102  distributionModel(typeName, dict, rndGen),
103  xy_(distributionModelDict_.lookup("distribution")),
104  nEntries_(xy_.size()),
105  meanValue_(0),
106  integral_(nEntries_),
107  cumulative_(distributionModelDict_.getOrDefault("cumulative", false))
108 {
109  minValue_ = xy_[0][0];
110  maxValue_ = xy_[nEntries_-1][0];
111 
112  check();
113 
114  // Additional sanity checks
115  if (cumulative_ && xy_[0][1] != 0)
116  {
118  << type() << "distribution: "
119  << "Elements in the second column for cumulative "
120  << "distribution functions must start from zero." << nl
121  << "First element = " << xy_[0][1]
122  << exit(FatalError);
123  }
124 
125  for (label i = 0; i < nEntries_; ++i)
126  {
127  if (i > 0 && xy_[i][0] <= xy_[i-1][0])
128  {
130  << type() << "distribution: "
131  << "Elements in the first column must "
132  << "be specified in an ascending order." << nl
133  << "Please see the row i = " << i << nl
134  << "xy[i-1] = " << xy_[i-1] << nl
135  << "xy[i] = " << xy_[i]
136  << exit(FatalError);
137  }
138 
139  if (xy_[i][0] < 0 || xy_[i][1] < 0)
140  {
142  << type() << "distribution: "
143  << "Input pairs cannot contain any negative element." << nl
144  << "Please see the row i = " << i << nl
145  << "xy[i] = " << xy_[i]
146  << exit(FatalError);
147  }
148 
149  if (cumulative_ && i > 0 && xy_[i][1] < xy_[i-1][1])
150  {
152  << type() << "distribution: "
153  << "Elements in the second column for cumulative "
154  << "distribution functions must be non-decreasing." << nl
155  << "Please see the row i = " << i << nl
156  << "xy[i-1] = " << xy_[i-1] << nl
157  << "xy[i] = " << xy_[i]
158  << exit(FatalError);
159  }
160  }
161 
162  initialise();
163 }
164 
165 
167 (
168  const UList<scalar>& sampleData,
169  const scalar binWidth,
170  Random& rndGen
171 )
172 :
174  xy_(),
175  nEntries_(0),
176  meanValue_(0),
177  integral_(),
178  cumulative_(false)
179 {
180  minValue_ = GREAT;
181  maxValue_ = -GREAT;
182  forAll(sampleData, i)
183  {
184  minValue_ = min(minValue_, sampleData[i]);
185  maxValue_ = max(maxValue_, sampleData[i]);
186  }
187 
188  label bin0 = floor(minValue_/binWidth);
189  label bin1 = ceil(maxValue_/binWidth);
190  nEntries_ = bin1 - bin0;
191 
192  if (nEntries_ == 0)
193  {
195  << "Data cannot be binned - zero bins generated" << nl
196  << " Bin width : " << binWidth << nl
197  << " Sample data : " << sampleData
198  << endl;
199 
200  return;
201  }
202 
203  xy_.setSize(nEntries_);
204 
205  // Populate bin boundaries and initialise occurrences
206  for (label bini = 0; bini < nEntries_; ++bini)
207  {
208  xy_[bini][0] = (bin0 + bini)*binWidth;
209  xy_[bini][1] = 0;
210  }
211 
212  // Populate occurrences
213  forAll(sampleData, i)
214  {
215  label bini = floor(sampleData[i]/binWidth) - bin0;
216  xy_[bini][1]++;
217  }
218 
219  initialise();
220 }
221 
222 
224 :
226  xy_(p.xy_),
227  nEntries_(p.nEntries_),
228  meanValue_(p.meanValue_),
229  integral_(p.integral_),
230  cumulative_(p.cumulative_)
231 {}
232 
233 
234 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
235 
237 {
238  const scalar u = rndGen_.sample01<scalar>();
239 
240  // Find the interval where u is in the table
241  label n = 1;
242  while (integral_[n] <= u)
243  {
244  n++;
245  }
246 
247  const scalar k = (xy_[n][1] - xy_[n-1][1])/(xy_[n][0] - xy_[n-1][0]);
248  const scalar d = xy_[n-1][1] - k*xy_[n-1][0];
249 
250  if (cumulative_)
251  {
252  return (u - d)/k;
253  }
254 
255  const scalar alpha =
256  u + xy_[n-1][0]*(0.5*k*xy_[n-1][0] + d) - integral_[n-1];
257 
258  // If k is small it is a linear equation, otherwise it is of second order
259  if (mag(k) > SMALL)
260  {
261  const scalar p = 2.0*d/k;
262  const scalar q = -2.0*alpha/k;
263  const scalar sqrtEr = sqrt(0.25*p*p - q);
264 
265  const scalar x1 = -0.5*p + sqrtEr;
266  const scalar x2 = -0.5*p - sqrtEr;
267  if ((x1 >= xy_[n-1][0]) && (x1 <= xy_[n][0]))
268  {
269  return x1;
270  }
271  else
272  {
273  return x2;
274  }
275  }
276 
277  return alpha/d;
278 }
279 
280 
282 {
283  return meanValue_;
284 }
285 
286 
288 {
289 // distributionModel::readData(is);
290  is >> xy_;
291  initialise();
292 }
293 
294 
296 {
297 // distributionModel::writeData(os);
298  os << xy_;
299 }
300 
301 
303 (
304  const word& dictName
305 ) const
306 {
307  // dictionary dict = distributionModel::writeDict(dictName);
309  dict.add("x", x());
310  dict.add("y", y());
311 
312  return dict;
313 }
314 
315 
317 {
318  // distributionModel::readDict(dict);
319  List<scalar> x(dict.lookup("x"));
320  List<scalar> y(dict.lookup("y"));
321 
322  xy_.setSize(x.size());
323  forAll(xy_, i)
324  {
325  xy_[i][0] = x[i];
326  xy_[i][1] = y[i];
327  }
328 
329  initialise();
330 }
331 
332 
335 {
336  auto tx = tmp<scalarField>::New(xy_.size());
337  auto& x = tx.ref();
338  forAll(xy_, i)
339  {
340  x[i] = xy_[i][0];
341  }
342 
343  return tx;
344 }
345 
346 
349 {
350  auto ty = tmp<scalarField>::New(xy_.size());
351  auto& y = ty.ref();
352  forAll(xy_, i)
353  {
354  y[i] = xy_[i][1];
355  }
356 
357  return ty;
358 }
359 
360 
361 Foam::Ostream& Foam::operator<<
362 (
363  Ostream& os,
365 )
366 {
367  os.check(FUNCTION_NAME);
368 
369  b.writeData(os);
370  return os;
371 }
372 
373 
375 {
376  is.check(FUNCTION_NAME);
377 
378  b.readData(is);
379  return is;
380 }
381 
382 
383 // ************************************************************************* //
Foam::Random
Random number generator.
Definition: Random.H:59
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::distributionModels::general::x
virtual tmp< scalarField > x() const
Bin boundaries.
Definition: general.C:334
Foam::y1
dimensionedScalar y1(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:282
Foam::distributionModels::general::meanValue
virtual scalar meanValue() const
Return the arithmetic mean of the distribution data.
Definition: general.C:281
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
dictName
const word dictName("faMeshDefinition")
Foam::distributionModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(distributionModel, binned, dictionary)
general.H
Foam::operator>>
Istream & operator>>(Istream &, directionInfo &)
Definition: directionInfo.C:230
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::distributionModel
A library of runtime-selectable doubly-truncated probability distribution models. Returns random samp...
Definition: distributionModel.H:72
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::distributionModel::maxValue_
scalar maxValue_
Maximum of the distribution.
Definition: distributionModel.H:88
Foam::distributionModels::general
Particle-size distribution model wherein random samples are drawn from a given arbitrary probability ...
Definition: general.H:173
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::distributionModels::general::y
virtual tmp< scalarField > y() const
Probabilities.
Definition: general.C:348
Foam::dictionary::null
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition: dictionary.H:392
Foam::distributionModels::general::sample
virtual scalar sample() const
Sample the distribution.
Definition: general.C:236
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::check
static void check(const int retVal, const char *what)
Definition: ptscotchDecomp.C:80
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
general
General relative velocity model.
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::distributionModels::general::readDict
virtual void readDict(const dictionary &dict)
Read data from dictionary.
Definition: general.C:316
Foam::y0
dimensionedScalar y0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:281
Foam::distributionModels::general::writeDict
virtual dictionary writeDict(const word &dictName) const
Write data in dictionary format.
Definition: general.C:303
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::distributionModels::general::general
general(const dictionary &dict, Random &rndGen)
Construct from components.
Definition: general.C:97
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::List< scalar >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::distributionModels::defineTypeNameAndDebug
defineTypeNameAndDebug(binned, 0)
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:103
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::distributionModels::general::readData
virtual void readData(Istream &os)
Read data from stream.
Definition: general.C:287
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
x
x
Definition: LISASMDCalcMethod2.H:52
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:295
Foam::fieldTypes::area
const wordList area
Standard area field types (scalar, vector, tensor, etc)
rndGen
Random rndGen
Definition: createFields.H:23
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::distributionModels::general::writeData
virtual void writeData(Ostream &os) const
Write data to stream.
Definition: general.C:295
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
y
scalar y
Definition: LISASMDCalcMethod1.H:14