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-2016 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  forAll(xy_, bini)
62  {
63  xy_[bini][1] /= sumProb;
64  }
65 
66  // Calculate the mean value
67  label bini = 0;
68  forAll(xy_, i)
69  {
70  if (xy_[i][1] > 0.5)
71  {
72  bini = i;
73  break;
74  }
75  }
76 
77  meanValue_ = xy_[bini][1];
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
82 
84 (
85  const dictionary& dict,
86  Random& rndGen
87 )
88 :
89  distributionModel(typeName, dict, rndGen),
90  xy_(distributionModelDict_.lookup("distribution")),
91  meanValue_(0)
92 {
93  if (maxValue() < minValue())
94  {
96  << "Maximum value is smaller than the minimum value:"
97  << " maxValue = " << maxValue()
98  << ", minValue = " << minValue()
99  << exit(FatalError);
100  }
101 
102  initialise();
103 }
104 
105 
107 (
108  const UList<scalar>& sampleData,
109  const scalar binWidth,
110  Random& rndGen
111 )
112 :
114  xy_(),
115  meanValue_(0)
116 {
117  scalar minValue = GREAT;
118  scalar maxValue = -GREAT;
119  forAll(sampleData, i)
120  {
121  minValue = min(minValue, sampleData[i]);
122  maxValue = max(maxValue, sampleData[i]);
123  }
124 
125  const label bin0 = floor(minValue/binWidth);
126  const label bin1 = ceil(maxValue/binWidth);
127  const label nBin = bin1 - bin0;
128 
129  if (nBin == 0)
130  {
132  << "Data cannot be binned - zero bins generated" << nl
133  << " Bin width : " << binWidth << nl
134  << " Sample data : " << sampleData
135  << endl;
136 
137  return;
138  }
139 
140  // Populate bin boundaries and initialise occurrences
141  xy_.setSize(nBin);
142  forAll(xy_, bini)
143  {
144  xy_[bini][0] = (bin0 + bini)*binWidth;
145  xy_[bini][1] = 0;
146  }
147 
148  // Bin the data
149  forAll(sampleData, i)
150  {
151  // Choose the nearest bin
152  label bini = floor(sampleData[i]/binWidth) - bin0;
153  label binii = min(bini + 1, nBin - 1);
154 
155  scalar d1 = mag(sampleData[i] - xy_[bini][0]);
156  scalar d2 = mag(xy_[binii][0] - sampleData[i]);
157 
158  if (d1 < d2)
159  {
160  xy_[bini][1]++;
161  }
162  else
163  {
164  xy_[binii][1]++;
165  }
166  }
167 
168  initialise();
169 }
170 
171 
173 :
175  xy_(p.xy_),
176  meanValue_(p.meanValue_)
177 {}
178 
179 
180 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
181 
183 {}
184 
185 
186 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
187 
189 {
190  scalar y = rndGen_.sample01<scalar>();
191 
192  for (label i = 0; i < xy_.size() - 1; ++i)
193  {
194  if (xy_[i][1] > y)
195  {
196  return xy_[i][0];
197  }
198  }
199 
200  return xy_.last()[0];
201 }
202 
203 
205 {
206  return xy_.first()[0];
207 }
208 
209 
211 {
212  return xy_.last()[0];
213 }
214 
215 
217 {
218  return meanValue_;
219 }
220 
221 
223 {
224 // distributionModel::readData(is);
225  is >> xy_;
226 }
227 
228 
230 {
231 // distributionModel::writeData(os);
232  os << xy_ ;
233 }
234 
235 
237 (
238  const word& dictName
239 ) const
240 {
241 // dictionary dict = distributionModel::writeDict(dictName);
243  dict.add("distribution", xy_);
244 
245  return dict;
246 }
247 
248 
250 {
251 // distributionModel::readDict(dict);
252  dict.readEntry("distribution", xy_);
253 }
254 
255 
256 Foam::Ostream& Foam::operator<<
257 (
258  Ostream& os,
260 )
261 {
262  os.check(FUNCTION_NAME);
263 
264  b.writeData(os);
265  return os;
266 }
267 
268 
270 {
271  is.check(FUNCTION_NAME);
272 
273  b.readData(is);
274  return is;
275 }
276 
277 
278 // ************************************************************************* //
Foam::Random
Random number generator.
Definition: Random.H:59
Foam::distributionModels::binned::readDict
virtual void readDict(const dictionary &dict)
Read data from dictionary.
Definition: binned.C:249
Foam::distributionModels::binned::maxValue
virtual scalar maxValue() const
Return the maximum value.
Definition: binned.C:210
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::distributionModels::binned::writeData
virtual void writeData(Ostream &os) const
Write data to stream.
Definition: binned.C:229
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
binned.H
Foam::distributionModels::binned::readData
virtual void readData(Istream &os)
Read data from stream.
Definition: binned.C:222
Foam::distributionModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(distributionModel, binned, dictionary)
dictName
const word dictName("blockMeshDict")
Foam::distributionModels::binned
A binned distribution model where the random sample is taken from a discrete set of bins.
Definition: binned.H:68
Foam::operator>>
Istream & operator>>(Istream &, directionInfo &)
Definition: directionInfo.C:230
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::distributionModel
A library of runtime-selectable distribution models.
Definition: distributionModel.H:70
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
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::dictionary::null
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition: dictionary.H:385
Foam::distributionModels::binned::meanValue
virtual scalar meanValue() const
Return the mean value.
Definition: binned.C:216
Foam::distributionModels::binned::sample
virtual scalar sample() const
Sample the distributionModel.
Definition: binned.C:188
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::distributionModels::binned::minValue
virtual scalar minValue() const
Return the minimum value.
Definition: binned.C:204
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::distributionModels::binned::~binned
virtual ~binned()
Destructor.
Definition: binned.C:182
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:51
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:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::distributionModels::binned::header
static const char * header
Definition: binned.H:94
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::distributionModels::binned::binned
binned(const dictionary &dict, Random &rndGen)
Construct from dictionary.
Definition: binned.C:84
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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
Foam::distributionModels::binned::writeDict
virtual dictionary writeDict(const word &dictName) const
Write data in dictionary format.
Definition: binned.C:237
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:270
minValue
scalar minValue
Definition: LISASMDCalcMethod2.H:12
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
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303
y
scalar y
Definition: LISASMDCalcMethod1.H:14
maxValue
scalar maxValue
Definition: LISASMDCalcMethod1.H:5