invIncGamma.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) 2016 OpenFOAM Foundation
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 Global
27  Foam::invIncGamma
28 
29 Description
30  Function to calculate the inverse of the
31  normalized incomplete gamma function.
32 
33  The algorithm is described in detain in reference:
34  \verbatim
35  DiDonato, A. R., & Morris Jr, A. H. (1986).
36  Computation of the incomplete gamma function ratios and their inverse.
37  ACM Transactions on Mathematical Software (TOMS), 12(4), 377-393.
38  \endverbatim
39 
40  All equation numbers in the following code refer to the above text and
41  the algorithm in the function 'invIncGamma' is described in section 4.
42 
43 \*---------------------------------------------------------------------------*/
44 
45 #include "mathematicalConstants.H"
46 using namespace Foam::constant::mathematical;
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 static scalar minimaxs(const scalar P)
56 {
57  // Eqn. 32
58 
59  static const scalar a[] =
60  {
61  3.31125922108741,
62  11.6616720288968,
63  4.28342155967104,
64  0.213623493715853
65  };
66 
67  static const scalar b[] =
68  {
69  6.61053765625462,
70  6.40691597760039,
71  1.27364489782223,
72  0.03611708101884203
73  };
74 
75  const scalar t = P < 0.5 ? sqrt(-2*log(P)) : sqrt(-2*log(1 - P));
76 
77  const scalar s =
78  t
79  - (a[0] + t*(a[1] + t*(a[2] + t*a[3])))
80  /(1 + t*(b[0] + t*(b[1] + t*(b[2] + t*b[3]))));
81 
82  return P < 0.5 ? -s : s;
83 }
84 
85 
86 static scalar Sn(const scalar a, const scalar x)
87 {
88  //- Eqn. 34
89 
90  scalar Sn = 1;
91  scalar Si = 1;
92 
93  for (int i=1; i<100; i++)
94  {
95  Si *= x/(a + i);
96  Sn += Si;
97 
98  if (Si < 1e-4) break;
99  }
100 
101  return Sn;
102 }
103 
104 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
105 
106 } // End namespace Foam
107 
108 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
109 
110 //- Inverse normalized incomplete gamma function
111 Foam::scalar Foam::invIncGamma(const scalar a, const scalar P)
112 {
113  const scalar Q = 1 - P;
114 
115  if (a == 1)
116  {
117  return -log(Q);
118  }
119  else if (a < 1)
120  {
121  const scalar Ga = tgamma(a);
122  const scalar B = Q*Ga;
123 
124  if (B > 0.6 || (B >= 0.45 && a >= 0.3))
125  {
126  // Eqn. 21
127  const scalar u =
128  (B*Q > 1e-8) ? pow(P*Ga*a, 1/a) : exp((-Q/a) - Eu);
129 
130  return u/(1 - (u/(a + 1)));
131  }
132  else if (a < 0.3 && B >= 0.35)
133  {
134  // Eqn. 22
135  const scalar t = exp(-Eu - B);
136  const scalar u = t*exp(t);
137 
138  return t*exp(u);
139  }
140  else if (B > 0.15 || a >= 0.3)
141  {
142  // Eqn. 23
143  const scalar y = -log(B);
144  const scalar u = y - (1 - a)*log(y);
145 
146  return y - (1 - a)*log(u) - log(1 + (1 - a)/(1 + u));
147  }
148  else if (B > 0.1)
149  {
150  // Eqn. 24
151  const scalar y = -log(B);
152  const scalar u = y - (1 - a)*log(y);
153 
154  return y
155  - (1 - a)*log(u)
156  - log
157  (
158  (sqr(u) + 2*(3 - a)*u + (2 - a)*(3 - a))
159  /(sqr(u) + (5 - a)*u + 2)
160  );
161  }
162  else
163  {
164  // Eqn. 25:
165  const scalar y = -log(B);
166  const scalar c1 = (a - 1)*log(y);
167  const scalar c12 = c1*c1;
168  const scalar c13 = c12*c1;
169  const scalar c14 = c12*c12;
170  const scalar a2 = a*a;
171  const scalar a3 = a2*a;
172  const scalar c2 = (a - 1)*(1 + c1);
173  const scalar c3 = (a - 1)*(-(c12/2) + (a - 2)*c1 + (3*a - 5)/2);
174  const scalar c4 =
175  (a - 1)
176  *(
177  (c13/3)
178  - (3*a - 5)*c12/2
179  + (a2 - 6*a + 7)*c1
180  + (11*a2 - 46*a + 47)/6
181  );
182  const scalar c5 =
183  (a - 1)*(-(c14/4)
184  + (11*a - 17)*c13/6
185  + (-3*a2 + 13*a - 13)*c12
186  + (2*a3 - 25*a2 + 72*a - 61)*c1/2
187  + (25*a3 - 195*a2 + 477*a - 379)/12);
188  const scalar y2 = y*y;
189  const scalar y3 = y2*y;
190  const scalar y4 = y2*y2;
191 
192  return y + c1 + (c2/y) + (c3/y2) + (c4/y3) + (c5/y4);
193  }
194  }
195  else
196  {
197  // Eqn. 31:
198  scalar s = minimaxs(P);
199 
200  const scalar s2 = sqr(s);
201  const scalar s3 = s*s2;
202  const scalar s4 = s2*s2;
203  const scalar s5 = s*s4;
204  const scalar sqrta = sqrt(a);
205 
206  const scalar w =
207  a + s*sqrta + (s2 - 1)/3
208  + (s3 - 7*s)/(36*sqrta)
209  - (3*s4 + 7*s2 - 16)/(810*a)
210  + (9*s5 + 256*s3 - 433*s)/(38880*a*sqrta);
211 
212  if (a >= 500 && mag(1 - w/a) < 1e-6)
213  {
214  return w;
215  }
216  else if (P > 0.5)
217  {
218  if (w < 3*a)
219  {
220  return w;
221  }
222  else
223  {
224  const scalar D = max(scalar(2), scalar(a*(a - 1)));
225  const scalar lnGa = lgamma(a);
226  const scalar lnB = log(Q) + lnGa;
227 
228  if (lnB < -2.3*D)
229  {
230  // Eqn. 25:
231  const scalar y = -lnB;
232  const scalar c1 = (a - 1)*log(y);
233  const scalar c12 = c1*c1;
234  const scalar c13 = c12*c1;
235  const scalar c14 = c12*c12;
236  const scalar a2 = a*a;
237  const scalar a3 = a2*a;
238 
239  const scalar c2 = (a - 1)*(1 + c1);
240  const scalar c3 =
241  (a - 1)
242  *(
243  - (c12/2)
244  + (a - 2)*c1
245  + (3*a - 5)/2
246  );
247  const scalar c4 =
248  (a - 1)
249  *(
250  (c13/3)
251  - (3*a - 5)*c12/2
252  + (a2 - 6*a + 7)*c1
253  + (11*a2 - 46*a + 47)/6
254  );
255  const scalar c5 =
256  (a - 1)
257  *(
258  - (c14/4)
259  + (11*a - 17)*c13/6
260  + (-3*a2 + 13*a - 13)*c12
261  + (2*a3 - 25*a2 + 72*a - 61)*c1/2
262  + (25*a3 - 195*a2 + 477*a - 379)/12
263  );
264 
265  const scalar y2 = y*y;
266  const scalar y3 = y2*y;
267  const scalar y4 = y2*y2;
268 
269  return y + c1 + (c2/y) + (c3/y2) + (c4/y3) + (c5/y4);
270  }
271  else
272  {
273  // Eqn. 33
274  const scalar u =
275  -lnB + (a - 1)*log(w) - log(1 + (1 - a)/(1 + w));
276 
277  return -lnB + (a - 1)*log(u) - log(1 + (1 - a)/(1 + u));
278  }
279  }
280  }
281  else
282  {
283  scalar z = w;
284  const scalar ap1 = a + 1;
285 
286  if (w < 0.15*ap1)
287  {
288  // Eqn. 35
289  const scalar ap2 = a + 2;
290  const scalar v = log(P) + lgamma(ap1);
291  z = exp((v + w)/a);
292  s = log1p(z/ap1*(1 + z/ap2));
293  z = exp((v + z - s)/a);
294  s = log1p(z/ap1*(1 + z/ap2));
295  z = exp((v + z - s)/a);
296  s = log1p(z/ap1*(1 + z/ap2*(1 + z/(a + 3))));
297  z = exp((v + z - s)/a);
298  }
299 
300  if (z <= 0.01*ap1 || z > 0.7*ap1)
301  {
302  return z;
303  }
304  else
305  {
306  // Eqn. 36
307  const scalar lnSn = log(Sn(a, z));
308  const scalar v = log(P) + lgamma(ap1);
309  z = exp((v + z - lnSn)/a);
310 
311  return z*(1 - (a*log(z) - z - v + lnSn)/(a - z));
312  }
313  }
314  }
315 }
316 
317 
318 // ************************************************************************* //
mathematicalConstants.H
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputSpray.H:25
B
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
Foam::constant::mathematical::e
constexpr scalar e(M_E)
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::invIncGamma
scalar invIncGamma(const scalar a, const scalar P)
Inverse normalized incomplete gamma function.
Definition: invIncGamma.C:108
Foam::lgamma
dimensionedScalar lgamma(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:278
Foam::constant::mathematical::Eu
constexpr scalar Eu(0.57721566490153286060651209)
Euler's constant.
Foam::constant::physicoChemical::c1
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
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
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::constant::physicoChemical::c2
const dimensionedScalar c2
Second radiation constant: default SI units: [m.K].
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::constant::mathematical
Mathematical constants.
Foam::Sn
static scalar Sn(const scalar a, const scalar x)
Definition: invIncGamma.C:83
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::minimaxs
static scalar minimaxs(const scalar P)
Definition: invIncGamma.C:52
D
const dimensionedScalar & D
Definition: solveBulkSurfactant.H:4
x
x
Definition: LISASMDCalcMethod2.H:52
y
scalar y
Definition: LISASMDCalcMethod1.H:14