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 | Copyright (C) 2011-2013 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 "general.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33  namespace distributionModels
34  {
35  defineTypeNameAndDebug(general, 0);
36  addToRunTimeSelectionTable(distributionModel, general, dictionary);
37  }
38 }
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
43 (
44  const dictionary& dict,
45  cachedRandom& rndGen
46 )
47 :
48  distributionModel(typeName, dict, rndGen),
49  xy_(distributionModelDict_.lookup("distribution")),
50  nEntries_(xy_.size()),
51  minValue_(xy_[0][0]),
52  maxValue_(xy_[nEntries_-1][0]),
53  meanValue_(0.0),
54  integral_(nEntries_)
55 {
56  check();
57 
58  // normalize the cumulative distributionModel
59 
60  integral_[0] = 0.0;
61  for (label i=1; i<nEntries_; i++)
62  {
63 
64  scalar k = (xy_[i][1] - xy_[i-1][1])/(xy_[i][0] - xy_[i-1][0]);
65  scalar d = xy_[i-1][1] - k*xy_[i-1][0];
66  scalar y1 = xy_[i][0]*(0.5*k*xy_[i][0] + d);
67  scalar y0 = xy_[i-1][0]*(0.5*k*xy_[i-1][0] + d);
68  scalar area = y1 - y0;
69 
70  integral_[i] = area + integral_[i-1];
71  }
72 
73  scalar sumArea = integral_.last();
74 
75  meanValue_ = sumArea/(maxValue_ - minValue_);
76 
77  for (label i=0; i<nEntries_; i++)
78  {
79  xy_[i][1] /= sumArea;
80  integral_[i] /= sumArea;
81  }
82 
83 }
84 
85 
87 :
89  xy_(p.xy_),
90  nEntries_(p.nEntries_),
91  minValue_(p.minValue_),
92  maxValue_(p.maxValue_),
93  integral_(p.integral_)
94 {}
95 
96 
97 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
98 
100 {}
101 
102 
103 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104 
106 {
107  scalar y = rndGen_.sample01<scalar>();
108 
109  // find the interval where y is in the table
110  label n=1;
111  while (integral_[n] <= y)
112  {
113  n++;
114  }
115 
116  scalar k = (xy_[n][1] - xy_[n-1][1])/(xy_[n][0] - xy_[n-1][0]);
117  scalar d = xy_[n-1][1] - k*xy_[n-1][0];
118 
119  scalar alpha = y + xy_[n-1][0]*(0.5*k*xy_[n-1][0] + d) - integral_[n-1];
120  scalar x = 0.0;
121 
122  // if k is small it is a linear equation, otherwise it is of second order
123  if (mag(k) > SMALL)
124  {
125  scalar p = 2.0*d/k;
126  scalar q = -2.0*alpha/k;
127  scalar sqrtEr = sqrt(0.25*p*p - q);
128 
129  scalar x1 = -0.5*p + sqrtEr;
130  scalar x2 = -0.5*p - sqrtEr;
131  if ((x1 >= xy_[n-1][0]) && (x1 <= xy_[n][0]))
132  {
133  x = x1;
134  }
135  else
136  {
137  x = x2;
138  }
139  }
140  else
141  {
142  x = alpha/d;
143  }
144 
145  return x;
146 }
147 
148 
150 {
151  return minValue_;
152 }
153 
154 
156 {
157  return maxValue_;
158 }
159 
160 
162 {
163  return meanValue_;
164 }
165 
166 
167 // ************************************************************************* //
virtual scalar maxValue() const
Return the maximum value.
Definition: general.C:155
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
general(const dictionary &dict, cachedRandom &rndGen)
Construct from components.
Definition: general.C:43
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensionedScalar sqrt(const dimensionedScalar &ds)
defineTypeNameAndDebug(distributionModel, 0)
virtual ~general()
Destructor.
Definition: general.C:99
dimensionedScalar y0(const dimensionedScalar &ds)
Random number generator.
Definition: cachedRandom.H:63
label k
Boltzmann constant.
Macros for easy insertion into run-time selection tables.
scalar y
Type sample01()
Return a sample whose components lie in the range 0-1.
virtual scalar minValue() const
Return the minimum value.
Definition: general.C:149
dimensionedScalar y1(const dimensionedScalar &ds)
virtual scalar sample() const
Sample the distributionModel.
Definition: general.C:105
cachedRandom & rndGen_
Reference to the random number generator.
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
addToRunTimeSelectionTable(distributionModel, exponential, dictionary)
volScalarField & p
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual scalar meanValue() const
Return the mean value.
Definition: general.C:161
Namespace for OpenFOAM.