incGamma.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) 2019 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::incGamma
28 
29 Description
30  Calculates the upper and lower incomplete gamma functions as well as their
31  normalized versions.
32 
33  The algorithm is described in detail in DiDonato et al. (1986).
34 
35  \verbatim
36  DiDonato, A. R., & Morris Jr, A. H. (1986).
37  Computation of the incomplete gamma function ratios and their inverse.
38  ACM Transactions on Mathematical Software (TOMS), 12(4), 377-393.
39  \endverbatim
40 
41  All equation numbers in the following code refer to the above paper.
42  The algorithm in function 'incGammaRatio_Q' is described in section 3.
43  The accuracy parameter IND is set to a value of 1.
44 
45 \*---------------------------------------------------------------------------*/
46 
47 #include "mathematicalConstants.H"
48 
49 using namespace Foam::constant::mathematical;
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 // Eqn. (13)
59 static scalar calcQE11(const scalar a, const scalar x, const int e = 30)
60 {
61  scalar a_2n = 0;
62  scalar b_2n = 1;
63 
64  scalar a_2np1 = 1;
65  scalar b_2np1 = x;
66 
67  int n = 1;
68  for (n = 1; (2*n) <= e; n++)
69  {
70  const scalar a_2nm1 = a_2np1;
71  const scalar b_2nm1 = b_2np1;
72 
73  a_2n = a_2nm1 + (n - a)*a_2n;
74  b_2n = b_2nm1 + (n - a)*b_2n;
75 
76  a_2np1 = x*a_2n + n*a_2nm1;
77  b_2np1 = x*b_2n + n*b_2nm1;
78  }
79 
80  if (2*(n - 1) < e)
81  {
82  return a_2np1/b_2np1;
83  }
84  else
85  {
86  return a_2n/b_2n;
87  }
88 }
89 
90 
91 // Eqn. (15)
92 static scalar calcPE15(const scalar a, const scalar x, const int nmax = 20)
93 {
94  scalar prod = 1;
95  scalar sum = 0;
96 
97  for (int n = 1; n <= nmax; n++)
98  {
99  prod *= (a + n);
100  sum += pow(x, n)/prod;
101  }
102 
103  const scalar R = (exp(-x)*pow(x, a))/tgamma(a);
104 
105  return R/a*(1 + sum);
106 }
107 
108 
109 // Eq. (16)
110 static scalar calcQE16(const scalar a, const scalar x, const int N = 20)
111 {
112  scalar an = 1;
113  scalar sum = 0;
114 
115  for (int n = 1; n <= (N - 1); n++)
116  {
117  an *= (a - n);
118  sum += an/pow(x, n);
119  }
120 
121  const scalar R = (exp(-x)*pow(x, a))/tgamma(a);
122 
123  return R/x*(1 + sum);
124 }
125 
126 
127 // Eq. (18)
128 static scalar calcTE18
129 (
130  const scalar a,
131  const scalar e0,
132  const scalar x,
133  const scalar lambda,
134  const scalar sigma,
135  const scalar phi
136 )
137 {
138  static const scalar D0[] =
139  {
140  -0.333333333333333E-00,
141  0.833333333333333E-01,
142  -0.148148148148148E-01,
143  0.115740740740741E-02,
144  0.352733686067019E-03,
145  -0.178755144032922E-03,
146  0.391926317852244E-04,
147  -0.218544851067999E-05,
148  -0.185406221071516E-05,
149  0.829671134095309E-06,
150  -0.176659527368261E-06,
151  0.670785354340150E-08,
152  0.102618097842403E-07,
153  -0.438203601845335E-08
154  };
155 
156  static const scalar D1[] =
157  {
158  -0.185185185185185E-02,
159  -0.347222222222222E-02,
160  0.264550264550265E-02,
161  -0.990226337448560E-03,
162  0.205761316872428E-03,
163  -0.401877572016461E-06,
164  -0.180985503344900E-04,
165  0.764916091608111E-05,
166  -0.161209008945634E-05,
167  0.464712780280743E-08,
168  0.137863344691572E-06,
169  -0.575254560351770E-07,
170  0.119516285997781E-07
171  };
172 
173  static const scalar D2[] =
174  {
175  0.413359788359788E-02,
176  -0.268132716049383E-02,
177  0.771604938271605E-03,
178  0.200938786008230E-05,
179  -0.107366532263652E-03,
180  0.529234488291201E-04,
181  -0.127606351886187E-04,
182  0.342357873409614E-07,
183  0.137219573090629E-05,
184  -0.629899213838006E-06,
185  0.142806142060642E-06
186  };
187 
188  const scalar u = 1/a;
189  scalar z = sqrt(2*phi);
190 
191  if (lambda < 1)
192  {
193  z = -z;
194  }
195 
196  if (sigma > (e0/sqrt(a)))
197  {
198  const scalar C0 =
199  D0[6]*pow6(z) + D0[5]*pow5(z) + D0[4]*pow4(z)
200  + D0[3]*pow3(z) + D0[2]*sqr(z) + D0[1]*z + D0[0];
201 
202  const scalar C1 =
203  D1[4]*pow4(z) + D1[3]*pow3(z) + D1[2]*sqr(z) + D1[1]*z + D1[0];
204 
205  const scalar C2 = D2[1]*z + D2[0];
206 
207  return C2*sqr(u) + C1*u + C0;
208  }
209  else
210  {
211  const scalar C0 = D0[2]*sqr(z) + D0[1]*z + D0[0];
212  const scalar C1 = D1[1]*z + D1[0];
213  const scalar C2 = D2[1]*z + D2[0];
214 
215  return C2*sqr(u) + C1*u + C0;
216  }
217 }
218 
219 
220 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
221 
222 } // End namespace Foam
223 
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
225 
226 Foam::scalar Foam::incGammaRatio_Q(const scalar a, const scalar x)
227 {
228  const scalar BIG = 14;
229  const scalar x0 = 17;
230  const scalar e0 = 0.025;
231 
232  if (a < 1)
233  {
234  if (a == 0.5)
235  {
236  // Eqn. (8)
237  if (x < 0.25)
238  {
239  return 1 - erf(sqrt(x));
240  }
241  else
242  {
243  return erfc(sqrt(x));
244  }
245  }
246  else if ( x < 1.1)
247  {
248  // Eqn. (12)
249  scalar alpha = x/2.59;
250 
251  if (x < 0.5)
252  {
253  alpha = log(sqrt(0.765))/log(x);
254  }
255 
256  scalar sum = 0;
257 
258  for (int n = 1; n <= 10; n++)
259  {
260  sum += pow((-x), n)/((a + n)*factorial(n));
261  }
262 
263  const scalar J = -a*sum;
264 
265  if (a > alpha || a == alpha)
266  {
267  // Eqn. (9)
268  return 1 - (pow(x, a)*(1 - J))/tgamma(a + 1);
269  }
270  else
271  {
272  // Eqn. (10)
273  const scalar L = exp(a*log(x)) - 1;
274  const scalar H = 1/(tgamma(a + 1)) - 1;
275 
276  return (pow(x, a)*J - L)/tgamma(a + 1) - H;
277  }
278  }
279  else
280  {
281  // Eqn. (11)
282  const scalar R = (exp(-x)*pow(x, a))/tgamma(a);
283 
284  return R*calcQE11(a, x);
285  }
286  }
287  else if (a >= BIG)
288  {
289  const scalar sigma = fabs(1 - x/a);
290 
291  if (sigma <= e0/sqrt(a))
292  {
293  // Eqn. (19)
294  const scalar lambda = x/a;
295  const scalar phi = lambda - 1 - log(lambda);
296  const scalar y = a*phi;
297 
298  const scalar E = 0.5 - (1 - y/3)*sqrt(y/pi);
299 
300  if (lambda <= 1)
301  {
302  return
303  1
304  - (
305  E
306  - (1 - y)/sqrt(2*pi*a)
307  *calcTE18(a, e0, x, lambda, sigma, phi)
308  );
309  }
310  else
311  {
312  return
313  E
314  + (1 - y)/sqrt(2*pi*a)
315  *calcTE18(a, e0, x, lambda, sigma, phi);
316  }
317  }
318  else
319  {
320  if (sigma <= 0.4)
321  {
322  // Eqn. (17)
323  const scalar lambda = x/a;
324  const scalar phi = lambda - 1 - log(lambda);
325  const scalar y = a*phi;
326 
327  if (lambda <= 1)
328  {
329  return
330  1
331  - (0.5*erfc(sqrt(y))
332  - exp(-y)/sqrt(2*pi*a)
333  *calcTE18(a, e0, x, lambda, sigma, phi));
334  }
335  else
336  {
337  return
338  0.5*erfc(sqrt(y))
339  + exp(-y)/sqrt(2*pi*a)
340  *calcTE18(a, e0, x, lambda, sigma, phi);
341  }
342  }
343  else
344  {
345  if (x <= max(a, log(10.0)))
346  {
347  // Eqn. (15)
348  return 1 - calcPE15(a, x);
349  }
350  else if (x < x0)
351  {
352  // Eqn. (11)
353  const scalar R = (exp(-x)*pow(x, a))/tgamma(a);
354 
355  return R*calcQE11(a, x);
356  }
357  else
358  {
359  // Eqn. (16)
360  return calcQE16(a, x);
361  }
362  }
363  }
364  }
365  else
366  {
367  if (a > x || x >= x0)
368  {
369  if (x <= max(a, log(10.0)))
370  {
371  // Eqn. (15)
372  return 1 - calcPE15(a, x);
373  }
374  else if ( x < x0)
375  {
376  // Eqn. (11)
377  const scalar R = (exp(-x)*pow(x, a))/tgamma(a);
378 
379  return R*calcQE11(a, x);
380  }
381  else
382  {
383  // Eqn. (16)
384  return calcQE16(a, x);
385  }
386  }
387  else
388  {
389  if (floor(2*a) == 2*a)
390  {
391  // Eqn. (14)
392  if (floor(a) == a)
393  {
394  scalar sum = 0;
395 
396  for (int n = 0; n <= (a - 1); n++)
397  {
398  sum += pow(x, n)/factorial(n);
399  }
400 
401  return exp(-x)*sum;
402  }
403  else
404  {
405  int i = a - 0.5;
406  scalar prod = 1;
407  scalar sum = 0;
408 
409  for (int n = 1; n <= i; n++)
410  {
411  prod *= (n - 0.5);
412  sum += pow(x, n)/prod;
413  }
414 
415  return erfc(sqrt(x)) + exp(-x)/sqrt(pi*x)*sum;
416  }
417  }
418  else if (x <= max(a, log(10.0)))
419  {
420  // Eqn. (15)
421  return 1 - calcPE15(a, x);
422  }
423  else if ( x < x0)
424  {
425  // Eqn. (11)
426  const scalar R = (exp(-x)*pow(x, a))/tgamma(a);
427 
428  return R*calcQE11(a, x);
429  }
430  else
431  {
432  // Eqn. (16)
433  return calcQE16(a, x);
434  }
435  }
436  }
437 }
438 
439 
440 Foam::scalar Foam::incGammaRatio_P(const scalar a, const scalar x)
441 {
442  return 1 - incGammaRatio_Q(a, x);
443 }
444 
445 
446 Foam::scalar Foam::incGamma_Q(const scalar a, const scalar x)
447 {
448  return incGammaRatio_Q(a, x)*tgamma(a);
449 }
450 
451 
452 Foam::scalar Foam::incGamma_P(const scalar a, const scalar x)
453 {
454  return incGammaRatio_P(a, x)*tgamma(a);
455 }
456 
457 
458 // ************************************************************************* //
L
const vector L(dict.get< vector >("L"))
mathematicalConstants.H
Foam::incGamma_Q
scalar incGamma_Q(const scalar a, const scalar x)
Upper incomplete gamma function.
Definition: incGamma.C:443
Foam::calcTE18
static scalar calcTE18(const scalar a, const scalar e0, const scalar x, const scalar lambda, const scalar sigma, const scalar phi)
Definition: incGamma.C:126
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
H
volScalarField H(IOobject("H", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
Foam::constant::mathematical::e
constexpr scalar e(M_E)
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::calcQE16
static scalar calcQE16(const scalar a, const scalar x, const int N=20)
Definition: incGamma.C:107
Foam::erf
dimensionedScalar erf(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:276
Foam::incGammaRatio_Q
scalar incGammaRatio_Q(const scalar a, const scalar x)
Normalized upper incomplete gamma function.
Definition: incGamma.C:223
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
R
#define R(A, B, C, D, E, F, K, M)
Foam::pow6
dimensionedScalar pow6(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:122
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::erfc
dimensionedScalar erfc(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:277
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
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::incGamma_P
scalar incGamma_P(const scalar a, const scalar x)
Lower incomplete gamma function.
Definition: incGamma.C:449
Foam::pow5
dimensionedScalar pow5(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:111
lambda
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::calcPE15
static scalar calcPE15(const scalar a, const scalar x, const int nmax=20)
Definition: incGamma.C:89
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::constant::mathematical
Mathematical constants.
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::incGammaRatio_P
scalar incGammaRatio_P(const scalar a, const scalar x)
Normalized lower incomplete gamma function.
Definition: incGamma.C:437
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::factorial
label factorial(label n)
Evaluate n! : 0 < n <= 12.
Definition: label.C:145
sigma
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
N
const Vector< label > N(dict.get< Vector< label >>("N"))
Foam::calcQE11
static scalar calcQE11(const scalar a, const scalar x, const int e=30)
Definition: incGamma.C:56
y
scalar y
Definition: LISASMDCalcMethod1.H:14