igam.c 3.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210
  1. /* igam.c
  2. *
  3. * Incomplete gamma integral
  4. *
  5. *
  6. *
  7. * SYNOPSIS:
  8. *
  9. * double a, x, y, igam();
  10. *
  11. * y = igam( a, x );
  12. *
  13. * DESCRIPTION:
  14. *
  15. * The function is defined by
  16. *
  17. * x
  18. * -
  19. * 1 | | -t a-1
  20. * igam(a,x) = ----- | e t dt.
  21. * - | |
  22. * | (a) -
  23. * 0
  24. *
  25. *
  26. * In this implementation both arguments must be positive.
  27. * The integral is evaluated by either a power series or
  28. * continued fraction expansion, depending on the relative
  29. * values of a and x.
  30. *
  31. * ACCURACY:
  32. *
  33. * Relative error:
  34. * arithmetic domain # trials peak rms
  35. * IEEE 0,30 200000 3.6e-14 2.9e-15
  36. * IEEE 0,100 300000 9.9e-14 1.5e-14
  37. */
  38. /* igamc()
  39. *
  40. * Complemented incomplete gamma integral
  41. *
  42. *
  43. *
  44. * SYNOPSIS:
  45. *
  46. * double a, x, y, igamc();
  47. *
  48. * y = igamc( a, x );
  49. *
  50. * DESCRIPTION:
  51. *
  52. * The function is defined by
  53. *
  54. *
  55. * igamc(a,x) = 1 - igam(a,x)
  56. *
  57. * inf.
  58. * -
  59. * 1 | | -t a-1
  60. * = ----- | e t dt.
  61. * - | |
  62. * | (a) -
  63. * x
  64. *
  65. *
  66. * In this implementation both arguments must be positive.
  67. * The integral is evaluated by either a power series or
  68. * continued fraction expansion, depending on the relative
  69. * values of a and x.
  70. *
  71. * ACCURACY:
  72. *
  73. * Tested at random a, x.
  74. * a x Relative error:
  75. * arithmetic domain domain # trials peak rms
  76. * IEEE 0.5,100 0,100 200000 1.9e-14 1.7e-15
  77. * IEEE 0.01,0.5 0,100 200000 1.4e-13 1.6e-15
  78. */
  79. /*
  80. Cephes Math Library Release 2.8: June, 2000
  81. Copyright 1985, 1987, 2000 by Stephen L. Moshier
  82. */
  83. #include <math.h>
  84. #ifdef ANSIPROT
  85. extern double lgam ( double );
  86. extern double exp ( double );
  87. extern double log ( double );
  88. extern double fabs ( double );
  89. extern double igam ( double, double );
  90. extern double igamc ( double, double );
  91. #else
  92. double lgam(), exp(), log(), fabs(), igam(), igamc();
  93. #endif
  94. extern double MACHEP, MAXLOG;
  95. static double big = 4.503599627370496e15;
  96. static double biginv = 2.22044604925031308085e-16;
  97. double igamc( a, x )
  98. double a, x;
  99. {
  100. double ans, ax, c, yc, r, t, y, z;
  101. double pk, pkm1, pkm2, qk, qkm1, qkm2;
  102. if( (x <= 0) || ( a <= 0) )
  103. return( 1.0 );
  104. if( (x < 1.0) || (x < a) )
  105. return( 1.0 - igam(a,x) );
  106. ax = a * log(x) - x - lgam(a);
  107. if( ax < -MAXLOG )
  108. {
  109. mtherr( "igamc", UNDERFLOW );
  110. return( 0.0 );
  111. }
  112. ax = exp(ax);
  113. /* continued fraction */
  114. y = 1.0 - a;
  115. z = x + y + 1.0;
  116. c = 0.0;
  117. pkm2 = 1.0;
  118. qkm2 = x;
  119. pkm1 = x + 1.0;
  120. qkm1 = z * x;
  121. ans = pkm1/qkm1;
  122. do
  123. {
  124. c += 1.0;
  125. y += 1.0;
  126. z += 2.0;
  127. yc = y * c;
  128. pk = pkm1 * z - pkm2 * yc;
  129. qk = qkm1 * z - qkm2 * yc;
  130. if( qk != 0 )
  131. {
  132. r = pk/qk;
  133. t = fabs( (ans - r)/r );
  134. ans = r;
  135. }
  136. else
  137. t = 1.0;
  138. pkm2 = pkm1;
  139. pkm1 = pk;
  140. qkm2 = qkm1;
  141. qkm1 = qk;
  142. if( fabs(pk) > big )
  143. {
  144. pkm2 *= biginv;
  145. pkm1 *= biginv;
  146. qkm2 *= biginv;
  147. qkm1 *= biginv;
  148. }
  149. }
  150. while( t > MACHEP );
  151. return( ans * ax );
  152. }
  153. /* left tail of incomplete gamma function:
  154. *
  155. * inf. k
  156. * a -x - x
  157. * x e > ----------
  158. * - -
  159. * k=0 | (a+k+1)
  160. *
  161. */
  162. double igam( a, x )
  163. double a, x;
  164. {
  165. double ans, ax, c, r;
  166. if( (x <= 0) || ( a <= 0) )
  167. return( 0.0 );
  168. if( (x > 1.0) && (x > a ) )
  169. return( 1.0 - igamc(a,x) );
  170. /* Compute x**a * exp(-x) / gamma(a) */
  171. ax = a * log(x) - x - lgam(a);
  172. if( ax < -MAXLOG )
  173. {
  174. mtherr( "igam", UNDERFLOW );
  175. return( 0.0 );
  176. }
  177. ax = exp(ax);
  178. /* power series */
  179. r = a;
  180. c = 1.0;
  181. ans = 1.0;
  182. do
  183. {
  184. r += 1.0;
  185. c *= x/r;
  186. ans += c;
  187. }
  188. while( c/ans > MACHEP );
  189. return( ans * ax/a );
  190. }