binned.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) 2015-2021 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "binned.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace distributionModels
36 {
39 }
40 }
41 
42 
44  "(bin probability)";
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::distributionModels::binned::initialise()
50 {
51  const label nSample(xy_.size());
52 
53  // Convert values to integral values
54  for (label bini = 1; bini < nSample; ++bini)
55  {
56  xy_[bini][1] += xy_[bini - 1][1];
57  }
58 
59  // Normalise
60  scalar sumProb = xy_.last()[1];
61 
62  if (sumProb < VSMALL)
63  {
65  << type() << "distribution: "
66  << "The sum of elements in the second column cannot be zero." << nl
67  << "sum = " << sumProb
68  << exit(FatalError);
69  }
70 
71  forAll(xy_, bini)
72  {
73  xy_[bini][1] /= sumProb;
74  }
75 
76  // Calculate the mean value
77  label bini = 0;
78  forAll(xy_, i)
79  {
80  if (xy_[i][1] > 0.5)
81  {
82  bini = i;
83  break;
84  }
85  }
86 
87  meanValue_ = xy_[bini][1];
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
92 
94 (
95  const dictionary& dict,
96  Random& rndGen
97 )
98 :
99  distributionModel(typeName, dict, rndGen),
100  xy_(distributionModelDict_.lookup("distribution")),
101  meanValue_(0)
102 {
103  minValue_ = xy_[0][0];
104  maxValue_ = xy_[xy_.size()-1][0];
105 
106  check();
107 
108  initialise();
109 }
110 
111 
113 (
114  const UList<scalar>& sampleData,
115  const scalar binWidth,
116  Random& rndGen
117 )
118 :
119  distributionModel(typeName, dictionary::null, rndGen),
120  xy_(),
121  meanValue_(0)
122 {
123  minValue_ = GREAT;
124  maxValue_ = -GREAT;
125  forAll(sampleData, i)
126  {
127  minValue_ = min(minValue_, sampleData[i]);
128  maxValue_ = max(maxValue_, sampleData[i]);
129  }
130 
131  const label bin0 = floor(minValue_/binWidth);
132  const label bin1 = ceil(maxValue_/binWidth);
133  const label nBin = bin1 - bin0;
134 
135  if (nBin == 0)
136  {
138  << "Data cannot be binned - zero bins generated" << nl
139  << " Bin width : " << binWidth << nl
140  << " Sample data : " << sampleData
141  << endl;
142 
143  return;
144  }
145 
146  // Populate bin boundaries and initialise occurrences
147  xy_.setSize(nBin);
148  forAll(xy_, bini)
149  {
150  xy_[bini][0] = (bin0 + bini)*binWidth;
151  xy_[bini][1] = 0;
152  }
153 
154  // Bin the data
155  forAll(sampleData, i)
156  {
157  // Choose the nearest bin
158  label bini = floor(sampleData[i]/binWidth) - bin0;
159  label binii = min(bini + 1, nBin - 1);
160 
161  scalar d1 = mag(sampleData[i] - xy_[bini][0]);
162  scalar d2 = mag(xy_[binii][0] - sampleData[i]);
163 
164  if (d1 < d2)
165  {
166  xy_[bini][1]++;
167  }
168  else
169  {
170  xy_[binii][1]++;
171  }
172  }
173 
174  initialise();
175 }
176 
177 
179 :
180  distributionModel(p),
181  xy_(p.xy_),
182  meanValue_(p.meanValue_)
183 {}
184 
185 
186 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
187 
188 Foam::scalar Foam::distributionModels::binned::sample() const
189 {
190  const scalar u = rndGen_.sample01<scalar>();
191 
192  for (label i = 0; i < xy_.size() - 1; ++i)
193  {
194  if (xy_[i][1] > u)
195  {
196  return xy_[i][0];
197  }
198  }
199 
200  return maxValue_;
201 }
202 
205 {
206  return meanValue_;
207 }
208 
209 
211 {
212 // distributionModel::readData(is);
213  is >> xy_;
214 }
215 
216 
218 {
219 // distributionModel::writeData(os);
220  os << xy_ ;
221 }
222 
223 
225 (
226  const word& dictName
227 ) const
228 {
229 // dictionary dict = distributionModel::writeDict(dictName);
231  dict.add("distribution", xy_);
232 
233  return dict;
234 }
235 
236 
238 {
239 // distributionModel::readDict(dict);
240  dict.readEntry("distribution", xy_);
241 }
242 
243 
244 Foam::Ostream& Foam::operator<<
245 (
246  Ostream& os,
247  const distributionModels::binned& b
248 )
249 {
251 
252  b.writeData(os);
253  return os;
254 }
255 
256 
257 Foam::Istream& Foam::operator>>(Istream& is, distributionModels::binned& b)
258 {
259  is.check(FUNCTION_NAME);
260 
261  b.readData(is);
262  return is;
263 }
264 
265 
266 // ************************************************************************* //
dictionary dict
scalar maxValue_
Maximum of the distribution.
virtual void readData(Istream &os)
Read data from stream.
Definition: binned.C:203
virtual scalar meanValue() const
Return the arithmetic mean of the distribution data.
Definition: binned.C:197
Particle-size distribution model wherein random samples are drawn from a given discrete set of (bin...
Definition: binned.H:146
binned(const dictionary &dict, Random &rndGen)
Construct from dictionary.
Definition: binned.C:87
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:45
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Random rndGen
Definition: createFields.H:23
const word dictName("faMeshDefinition")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
addToRunTimeSelectionTable(distributionModel, binned, dictionary)
Lookup type of boundary radiation properties.
Definition: lookup.H:57
virtual void writeData(Ostream &os) const
Write data to stream.
Definition: binned.C:210
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static const char * header
Definition: binned.H:178
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
defineTypeNameAndDebug(binned, 0)
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual void check() const
Check that the distribution model is valid.
Istream & operator>>(Istream &, directionInfo &)
scalar minValue_
Minimum of the distribution.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
virtual scalar sample() const
Sample the distribution.
Definition: binned.C:181
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:105
Random number generator.
Definition: Random.H:55
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
virtual dictionary writeDict(const word &dictName) const
Write data in dictionary format.
Definition: binned.C:218
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
virtual void readDict(const dictionary &dict)
Read data from dictionary.
Definition: binned.C:230
#define WarningInFunction
Report a warning using Foam::Warning.
A library of runtime-selectable doubly-truncated probability distribution models. Returns random samp...
volScalarField & p
Namespace for OpenFOAM.