File size: 12,877 Bytes
7885a28
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
/* Translated into C++ by SciPy developers in 2024.
 * Original header with Copyright information appears below.
 */

/*                                                     igam.c
 *
 *     Incomplete Gamma integral
 *
 *
 *
 * SYNOPSIS:
 *
 * double a, x, y, igam();
 *
 * y = igam( a, x );
 *
 * DESCRIPTION:
 *
 * The function is defined by
 *
 *                           x
 *                            -
 *                   1       | |  -t  a-1
 *  igam(a,x)  =   -----     |   e   t   dt.
 *                  -      | |
 *                 | (a)    -
 *                           0
 *
 *
 * In this implementation both arguments must be positive.
 * The integral is evaluated by either a power series or
 * continued fraction expansion, depending on the relative
 * values of a and x.
 *
 * ACCURACY:
 *
 *                      Relative error:
 * arithmetic   domain     # trials      peak         rms
 *    IEEE      0,30       200000       3.6e-14     2.9e-15
 *    IEEE      0,100      300000       9.9e-14     1.5e-14
 */
/*							igamc()
 *
 *	Complemented incomplete Gamma integral
 *
 *
 *
 * SYNOPSIS:
 *
 * double a, x, y, igamc();
 *
 * y = igamc( a, x );
 *
 * DESCRIPTION:
 *
 * The function is defined by
 *
 *
 *  igamc(a,x)   =   1 - igam(a,x)
 *
 *                            inf.
 *                              -
 *                     1       | |  -t  a-1
 *               =   -----     |   e   t   dt.
 *                    -      | |
 *                   | (a)    -
 *                             x
 *
 *
 * In this implementation both arguments must be positive.
 * The integral is evaluated by either a power series or
 * continued fraction expansion, depending on the relative
 * values of a and x.
 *
 * ACCURACY:
 *
 * Tested at random a, x.
 *                a         x                      Relative error:
 * arithmetic   domain   domain     # trials      peak         rms
 *    IEEE     0.5,100   0,100      200000       1.9e-14     1.7e-15
 *    IEEE     0.01,0.5  0,100      200000       1.4e-13     1.6e-15
 */

/*
 * Cephes Math Library Release 2.0:  April, 1987
 * Copyright 1985, 1987 by Stephen L. Moshier
 * Direct inquiries to 30 Frost Street, Cambridge, MA 02140
 */

/* Sources
 * [1] "The Digital Library of Mathematical Functions", dlmf.nist.gov
 * [2] Maddock et. al., "Incomplete Gamma Functions",
 *     https://www.boost.org/doc/libs/1_61_0/libs/math/doc/html/math_toolkit/sf_gamma/igamma.html
 */

/* Scipy changes:
 * - 05-01-2016: added asymptotic expansion for igam to improve the
 *   a ~ x regime.
 * - 06-19-2016: additional series expansion added for igamc to
 *   improve accuracy at small arguments.
 * - 06-24-2016: better choice of domain for the asymptotic series;
 *   improvements in accuracy for the asymptotic series when a and x
 *   are very close.
 */
#pragma once

#include "../config.h"
#include "../error.h"

#include "const.h"
#include "gamma.h"
#include "igam_asymp_coeff.h"
#include "lanczos.h"
#include "ndtr.h"
#include "unity.h"

namespace xsf {
namespace cephes {

    namespace detail {

        constexpr int igam_MAXITER = 2000;
        constexpr int IGAM = 1;
        constexpr int IGAMC = 0;
        constexpr double igam_SMALL = 20;
        constexpr double igam_LARGE = 200;
        constexpr double igam_SMALLRATIO = 0.3;
        constexpr double igam_LARGERATIO = 4.5;

        constexpr double igam_big = 4.503599627370496e15;
        constexpr double igam_biginv = 2.22044604925031308085e-16;

        /* Compute
         *
         * x^a * exp(-x) / gamma(a)
         *
         * corrected from (15) and (16) in [2] by replacing exp(x - a) with
         * exp(a - x).
         */
        XSF_HOST_DEVICE inline double igam_fac(double a, double x) {
            double ax, fac, res, num;

            if (std::abs(a - x) > 0.4 * std::abs(a)) {
                ax = a * std::log(x) - x - xsf::cephes::lgam(a);
                if (ax < -MAXLOG) {
                    set_error("igam", SF_ERROR_UNDERFLOW, NULL);
                    return 0.0;
                }
                return std::exp(ax);
            }

            fac = a + xsf::cephes::lanczos_g - 0.5;
            res = std::sqrt(fac / std::exp(1)) / xsf::cephes::lanczos_sum_expg_scaled(a);

            if ((a < 200) && (x < 200)) {
                res *= std::exp(a - x) * std::pow(x / fac, a);
            } else {
                num = x - a - xsf::cephes::lanczos_g + 0.5;
                res *= std::exp(a * xsf::cephes::log1pmx(num / fac) + x * (0.5 - xsf::cephes::lanczos_g) / fac);
            }

            return res;
        }

        /* Compute igamc using DLMF 8.9.2. */
        XSF_HOST_DEVICE inline double igamc_continued_fraction(double a, double x) {
            int i;
            double ans, ax, c, yc, r, t, y, z;
            double pk, pkm1, pkm2, qk, qkm1, qkm2;

            ax = igam_fac(a, x);
            if (ax == 0.0) {
                return 0.0;
            }

            /* continued fraction */
            y = 1.0 - a;
            z = x + y + 1.0;
            c = 0.0;
            pkm2 = 1.0;
            qkm2 = x;
            pkm1 = x + 1.0;
            qkm1 = z * x;
            ans = pkm1 / qkm1;

            for (i = 0; i < igam_MAXITER; i++) {
                c += 1.0;
                y += 1.0;
                z += 2.0;
                yc = y * c;
                pk = pkm1 * z - pkm2 * yc;
                qk = qkm1 * z - qkm2 * yc;
                if (qk != 0) {
                    r = pk / qk;
                    t = std::abs((ans - r) / r);
                    ans = r;
                } else
                    t = 1.0;
                pkm2 = pkm1;
                pkm1 = pk;
                qkm2 = qkm1;
                qkm1 = qk;
                if (std::abs(pk) > igam_big) {
                    pkm2 *= igam_biginv;
                    pkm1 *= igam_biginv;
                    qkm2 *= igam_biginv;
                    qkm1 *= igam_biginv;
                }
                if (t <= MACHEP) {
                    break;
                }
            }

            return (ans * ax);
        }

        /* Compute igam using DLMF 8.11.4. */
        XSF_HOST_DEVICE inline double igam_series(double a, double x) {
            int i;
            double ans, ax, c, r;

            ax = igam_fac(a, x);
            if (ax == 0.0) {
                return 0.0;
            }

            /* power series */
            r = a;
            c = 1.0;
            ans = 1.0;

            for (i = 0; i < igam_MAXITER; i++) {
                r += 1.0;
                c *= x / r;
                ans += c;
                if (c <= MACHEP * ans) {
                    break;
                }
            }

            return (ans * ax / a);
        }

        /* Compute igamc using DLMF 8.7.3. This is related to the series in
         * igam_series but extra care is taken to avoid cancellation.
         */
        XSF_HOST_DEVICE inline double igamc_series(double a, double x) {
            int n;
            double fac = 1;
            double sum = 0;
            double term, logx;

            for (n = 1; n < igam_MAXITER; n++) {
                fac *= -x / n;
                term = fac / (a + n);
                sum += term;
                if (std::abs(term) <= MACHEP * std::abs(sum)) {
                    break;
                }
            }

            logx = std::log(x);
            term = -xsf::cephes::expm1(a * logx - xsf::cephes::lgam1p(a));
            return term - std::exp(a * logx - xsf::cephes::lgam(a)) * sum;
        }

        /* Compute igam/igamc using DLMF 8.12.3/8.12.4. */
        XSF_HOST_DEVICE inline double asymptotic_series(double a, double x, int func) {
            int k, n, sgn;
            int maxpow = 0;
            double lambda = x / a;
            double sigma = (x - a) / a;
            double eta, res, ck, ckterm, term, absterm;
            double absoldterm = std::numeric_limits<double>::infinity();
            double etapow[detail::igam_asymp_coeff_N] = {1};
            double sum = 0;
            double afac = 1;

            if (func == detail::IGAM) {
                sgn = -1;
            } else {
                sgn = 1;
            }

            if (lambda > 1) {
                eta = std::sqrt(-2 * xsf::cephes::log1pmx(sigma));
            } else if (lambda < 1) {
                eta = -std::sqrt(-2 * xsf::cephes::log1pmx(sigma));
            } else {
                eta = 0;
            }
            res = 0.5 * xsf::cephes::erfc(sgn * eta * std::sqrt(a / 2));

            for (k = 0; k < igam_asymp_coeff_K; k++) {
                ck = igam_asymp_coeff_d[k][0];
                for (n = 1; n < igam_asymp_coeff_N; n++) {
                    if (n > maxpow) {
                        etapow[n] = eta * etapow[n - 1];
                        maxpow += 1;
                    }
                    ckterm = igam_asymp_coeff_d[k][n] * etapow[n];
                    ck += ckterm;
                    if (std::abs(ckterm) < MACHEP * std::abs(ck)) {
                        break;
                    }
                }
                term = ck * afac;
                absterm = std::abs(term);
                if (absterm > absoldterm) {
                    break;
                }
                sum += term;
                if (absterm < MACHEP * std::abs(sum)) {
                    break;
                }
                absoldterm = absterm;
                afac /= a;
            }
            res += sgn * std::exp(-0.5 * a * eta * eta) * sum / std::sqrt(2 * M_PI * a);

            return res;
        }

    } // namespace detail

    XSF_HOST_DEVICE inline double igamc(double a, double x);

    XSF_HOST_DEVICE inline double igam(double a, double x) {
        double absxma_a;

        if (x < 0 || a < 0) {
            set_error("gammainc", SF_ERROR_DOMAIN, NULL);
            return std::numeric_limits<double>::quiet_NaN();
        } else if (a == 0) {
            if (x > 0) {
                return 1;
            } else {
                return std::numeric_limits<double>::quiet_NaN();
            }
        } else if (x == 0) {
            /* Zero integration limit */
            return 0;
        } else if (std::isinf(a)) {
            if (std::isinf(x)) {
                return std::numeric_limits<double>::quiet_NaN();
            }
            return 0;
        } else if (std::isinf(x)) {
            return 1;
        }

        /* Asymptotic regime where a ~ x; see [2]. */
        absxma_a = std::abs(x - a) / a;
        if ((a > detail::igam_SMALL) && (a < detail::igam_LARGE) && (absxma_a < detail::igam_SMALLRATIO)) {
            return detail::asymptotic_series(a, x, detail::IGAM);
        } else if ((a > detail::igam_LARGE) && (absxma_a < detail::igam_LARGERATIO / std::sqrt(a))) {
            return detail::asymptotic_series(a, x, detail::IGAM);
        }

        if ((x > 1.0) && (x > a)) {
            return (1.0 - igamc(a, x));
        }

        return detail::igam_series(a, x);
    }

    XSF_HOST_DEVICE double igamc(double a, double x) {
        double absxma_a;

        if (x < 0 || a < 0) {
            set_error("gammaincc", SF_ERROR_DOMAIN, NULL);
            return std::numeric_limits<double>::quiet_NaN();
        } else if (a == 0) {
            if (x > 0) {
                return 0;
            } else {
                return std::numeric_limits<double>::quiet_NaN();
            }
        } else if (x == 0) {
            return 1;
        } else if (std::isinf(a)) {
            if (std::isinf(x)) {
                return std::numeric_limits<double>::quiet_NaN();
            }
            return 1;
        } else if (std::isinf(x)) {
            return 0;
        }

        /* Asymptotic regime where a ~ x; see [2]. */
        absxma_a = std::abs(x - a) / a;
        if ((a > detail::igam_SMALL) && (a < detail::igam_LARGE) && (absxma_a < detail::igam_SMALLRATIO)) {
            return detail::asymptotic_series(a, x, detail::IGAMC);
        } else if ((a > detail::igam_LARGE) && (absxma_a < detail::igam_LARGERATIO / std::sqrt(a))) {
            return detail::asymptotic_series(a, x, detail::IGAMC);
        }

        /* Everywhere else; see [2]. */
        if (x > 1.1) {
            if (x < a) {
                return 1.0 - detail::igam_series(a, x);
            } else {
                return detail::igamc_continued_fraction(a, x);
            }
        } else if (x <= 0.5) {
            if (-0.4 / std::log(x) < a) {
                return 1.0 - detail::igam_series(a, x);
            } else {
                return detail::igamc_series(a, x);
            }
        } else {
            if (x * 1.1 < a) {
                return 1.0 - detail::igam_series(a, x);
            } else {
                return detail::igamc_series(a, x);
            }
        }
    }

} // namespace cephes
} // namespace xsf