gammainc.js (9278B)
1 /** 2 * @license Apache-2.0 3 * 4 * Copyright (c) 2018 The Stdlib Authors. 5 * 6 * Licensed under the Apache License, Version 2.0 (the "License"); 7 * you may not use this file except in compliance with the License. 8 * You may obtain a copy of the License at 9 * 10 * http://www.apache.org/licenses/LICENSE-2.0 11 * 12 * Unless required by applicable law or agreed to in writing, software 13 * distributed under the License is distributed on an "AS IS" BASIS, 14 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. 15 * See the License for the specific language governing permissions and 16 * limitations under the License. 17 * 18 * 19 * ## Notice 20 * 21 * The original C++ code and copyright notice are from the [Boost library]{@link http://www.boost.org/doc/libs/1_62_0/boost/math/special_functions/gamma.hpp}. The implementation has been modified for JavaScript. 22 * 23 * ```text 24 * (C) Copyright John Maddock 2006-7, 2013-14. 25 * (C) Copyright Paul A. Bristow 2007, 2013-14. 26 * (C) Copyright Nikhar Agrawal 2013-14. 27 * (C) Christopher Kormanyos 2013-14. 28 * 29 * Use, modification and distribution are subject to the 30 * Boost Software License, Version 1.0. (See accompanying file 31 * LICENSE or copy at http://www.boost.org/LICENSE_1_0.txt) 32 * ``` 33 */ 34 35 'use strict'; 36 37 // MODULES // 38 39 var gammaln = require( './../../../../base/special/gammaln' ); 40 var floor = require( './../../../../base/special/floor' ); 41 var gamma = require( './../../../../base/special/gamma' ); 42 var abs = require( './../../../../base/special/abs' ); 43 var exp = require( './../../../../base/special/exp' ); 44 var pow = require( './../../../../base/special/pow' ); 45 var ln = require( './../../../../base/special/ln' ); 46 var SQRT_EPSILON = require( '@stdlib/constants/float64/sqrt-eps' ); 47 var FLOAT64_MAX = require( '@stdlib/constants/float64/max' ); 48 var SQRT_TWO_PI = require( '@stdlib/constants/float64/sqrt-two-pi' ); 49 var MAX_LN = require( '@stdlib/constants/float64/max-ln' ); 50 var PINF = require( '@stdlib/constants/float64/pinf' ); 51 var finiteGammaQ = require( './finite_gamma_q.js' ); 52 var finiteHalfGammaQ = require( './finite_half_gamma_q.js' ); 53 var fullIGammaPrefix = require( './full_igamma_prefix.js' ); 54 var igammaTemmeLarge = require( './igamma_temme_large.js' ); 55 var lowerGammaSeries = require( './lower_gamma_series.js' ); 56 var regularisedGammaPrefix = require( './regularised_gamma_prefix.js' ); 57 var tgammaSmallUpperPart = require( './tgamma_small_upper_part.js' ); 58 var upperGammaFraction = require( './upper_gamma_fraction.js' ); 59 60 61 // VARIABLES // 62 63 var MAX_FACTORIAL = 170; // TODO: consider extracting as a constant 64 65 66 // MAIN // 67 68 /** 69 * Computes the regularized incomplete gamma function. The upper tail is calculated via the modified Lentz's method for computing continued fractions, the lower tail using a power expansion. 70 * 71 * 72 * ## Notes 73 * 74 * - When `a >= MAX_FACTORIAL` and computing the non-normalized incomplete gamma, result is rather hard to compute unless we use logs. There are really two options a) if `x` is a long way from `a` in value then we can reliably use methods 2 and 4 below in logarithmic form and go straight to the result. Otherwise we let the regularized gamma take the strain (the result is unlikely to underflow in the central region anyway) and combine with `lgamma` in the hopes that we get a finite result. 75 * 76 * @param {NonNegativeNumber} x - function parameter 77 * @param {PositiveNumber} a - function parameter 78 * @param {boolean} [regularized=true] - boolean indicating if the function should evaluate the regularized or non-regularized incomplete gamma functions 79 * @param {boolean} [upper=false] - boolean indicating if the function should return the upper tail of the incomplete gamma function 80 * @returns {number} function value 81 */ 82 function gammainc( x, a, regularized, upper ) { 83 var optimisedInvert; 84 var normalized; 85 var evalMethod; 86 var initValue; 87 var isHalfInt; 88 var useTemme; 89 var isSmallA; 90 var invert; 91 var result; 92 var isInt; 93 var sigma; 94 var gam; 95 var res; 96 var fa; 97 var g; 98 99 if ( x < 0.0 || a <= 0.0 ) { 100 return NaN; 101 } 102 normalized = ( regularized === void 0 ) ? true : regularized; 103 invert = upper; 104 result = 0.0; 105 if ( a >= MAX_FACTORIAL && !normalized ) { 106 if ( invert && ( a * 4.0 < x ) ) { 107 // This is method 4 below, done in logs: 108 result = ( a * ln(x) ) - x; 109 result += ln( upperGammaFraction( a, x ) ); 110 } 111 else if ( !invert && ( a > 4.0 * x ) ) { 112 // This is method 2 below, done in logs: 113 result = ( a * ln(x) ) - x; 114 initValue = 0; 115 result += ln( lowerGammaSeries( a, x, initValue ) / a ); 116 } 117 else { 118 result = gammainc( a, x, true, invert ); 119 if ( result === 0.0 ) { 120 if ( invert ) { 121 // Try http://functions.wolfram.com/06.06.06.0039.01 122 result = 1.0 + ( 1.0 / (12.0*a) ) + ( 1.0 / (288.0*a*a) ); 123 result = ln( result ) - a + ( ( a-0.5 ) * ln(a) ); 124 result += ln( SQRT_TWO_PI ); 125 } else { 126 // This is method 2 below, done in logs, we're really outside the range of this method, but since the result is almost certainly infinite, we should probably be OK: 127 result = ( a * ln( x ) ) - x; 128 initValue = 0.0; 129 result += ln( lowerGammaSeries( a, x, initValue ) / a); 130 } 131 } 132 else { 133 result = ln( result ) + gammaln( a ); 134 } 135 } 136 if ( result > MAX_LN ) { 137 return PINF; 138 } 139 return exp( result ); 140 } 141 isSmallA = ( a < 30 ) && ( a <= x + 1.0 ) && ( x < MAX_LN ); 142 if ( isSmallA ) { 143 fa = floor( a ); 144 isInt = ( fa === a ); 145 isHalfInt = ( isInt ) ? false : ( abs( fa - a ) === 0.5 ); 146 } else { 147 isInt = isHalfInt = false; 148 } 149 if ( isInt && x > 0.6 ) { 150 // Calculate Q via finite sum: 151 invert = !invert; 152 evalMethod = 0; 153 } 154 else if ( isHalfInt && x > 0.2 ) { 155 // Calculate Q via finite sum for half integer a: 156 invert = !invert; 157 evalMethod = 1; 158 } 159 else if ( x < SQRT_EPSILON && a > 1.0 ) { 160 evalMethod = 6; 161 } 162 else if ( x < 0.5 ) { 163 // Changeover criterion chosen to give a changeover at Q ~ 0.33: 164 if ( -0.4 / ln( x ) < a ) { 165 evalMethod = 2; 166 } else { 167 evalMethod = 3; 168 } 169 } 170 else if ( x < 1.1 ) { 171 // Changeover here occurs when P ~ 0.75 or Q ~ 0.25: 172 if ( x * 0.75 < a ) { 173 evalMethod = 2; 174 } else { 175 evalMethod = 3; 176 } 177 } 178 else { 179 // Begin by testing whether we're in the "bad" zone where the result will be near 0.5 and the usual series and continued fractions are slow to converge: 180 useTemme = false; 181 if ( normalized && a > 20 ) { 182 sigma = abs( (x-a)/a ); 183 if ( a > 200 ) { 184 // Limit chosen so that we use Temme's expansion only if the result would be larger than about 10^-6. Below that the regular series and continued fractions converge OK, and if we use Temme's method we get increasing errors from the dominant erfc term as it's (inexact) argument increases in magnitude. 185 if ( 20 / a > sigma * sigma ) { 186 useTemme = true; 187 } 188 } else if ( sigma < 0.4 ) { 189 useTemme = true; 190 } 191 } 192 if ( useTemme ) { 193 evalMethod = 5; 194 } 195 // Regular case where the result will not be too close to 0.5: Changeover occurs at P ~ Q ~ 0.5. Note that series computation of P is about x2 faster than continued fraction calculation of Q, so try and use the CF only when really necessary, especially for small x. 196 else if ( x - ( 1.0 / (3.0 * x) ) < a ) { 197 evalMethod = 2; 198 } else { 199 evalMethod = 4; 200 invert = !invert; 201 } 202 } 203 204 /* eslint-disable default-case */ 205 switch ( evalMethod ) { 206 case 0: 207 result = finiteGammaQ( a, x ); 208 if (normalized === false ) { 209 result *= gamma( a ); 210 } 211 break; 212 case 1: 213 result = finiteHalfGammaQ( a, x ); 214 if ( normalized === false ) { 215 result *= gamma( a ); 216 } 217 break; 218 case 2: 219 // Compute P: 220 result = ( normalized ) ? 221 regularisedGammaPrefix( a, x ) : 222 fullIGammaPrefix( a, x ); 223 if ( result !== 0.0 ) { 224 initValue = 0.0; 225 optimisedInvert = false; 226 if ( invert ) { 227 initValue = ( normalized ) ? 1.0 : gamma(a); 228 if ( 229 normalized || 230 result >= 1.0 || 231 FLOAT64_MAX * result > initValue 232 ) { 233 initValue /= result; 234 if ( 235 normalized || 236 a < 1.0 || 237 ( FLOAT64_MAX / a > initValue ) 238 ) { 239 initValue *= -a; 240 optimisedInvert = true; 241 } 242 else { 243 initValue = 0.0; 244 } 245 } 246 else { 247 initValue = 0.0; 248 } 249 } 250 } 251 result *= lowerGammaSeries( a, x, initValue ) / a; 252 if ( optimisedInvert ) { 253 invert = false; 254 result = -result; 255 } 256 break; 257 case 3: 258 // Compute Q: 259 invert = !invert; 260 res = tgammaSmallUpperPart( a, x, invert ); 261 result = res[ 0 ]; 262 g = res[ 1 ]; 263 invert = false; 264 if ( normalized ) { 265 result /= g; 266 } 267 break; 268 case 4: 269 // Compute Q: 270 result = ( normalized ) ? 271 regularisedGammaPrefix( a, x ) : 272 fullIGammaPrefix( a, x ); 273 if ( result !== 0 ) { 274 result *= upperGammaFraction( a, x ); 275 } 276 break; 277 case 5: 278 result = igammaTemmeLarge( a, x ); 279 if ( x >= a ) { 280 invert = !invert; 281 } 282 break; 283 case 6: 284 // Since x is so small that P is necessarily very small too, use http://functions.wolfram.com/GammaBetaErf/GammaRegularized/06/01/05/01/01/ 285 result = ( normalized ) ? 286 pow(x, a) / gamma( a + 1.0 ) : 287 pow( x, a ) / a; 288 result *= 1.0 - ( a * x / ( a + 1.0 ) ); 289 break; 290 } 291 if ( normalized && result > 1.0 ) { 292 result = 1.0; 293 } 294 if ( invert ) { 295 gam = ( normalized ) ? 1.0 : gamma( a ); 296 result = gam - result; 297 } 298 return result; 299 } 300 301 302 // EXPORTS // 303 304 module.exports = gammainc;