igam.c 3.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197
  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 "mconf.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) double a, x;
  98. {
  99. double ans, ax, c, yc, r, t, y, z;
  100. double pk, pkm1, pkm2, qk, qkm1, qkm2;
  101. if ((x <= 0) || (a <= 0))
  102. return (1.0);
  103. if ((x < 1.0) || (x < a))
  104. return (1.0 - igam(a, x));
  105. ax = a * log(x) - x - lgam(a);
  106. if (ax < -MAXLOG) {
  107. mtherr("igamc", UNDERFLOW);
  108. return (0.0);
  109. }
  110. ax = exp(ax);
  111. /* continued fraction */
  112. y = 1.0 - a;
  113. z = x + y + 1.0;
  114. c = 0.0;
  115. pkm2 = 1.0;
  116. qkm2 = x;
  117. pkm1 = x + 1.0;
  118. qkm1 = z * x;
  119. ans = pkm1 / qkm1;
  120. do {
  121. c += 1.0;
  122. y += 1.0;
  123. z += 2.0;
  124. yc = y * c;
  125. pk = pkm1 * z - pkm2 * yc;
  126. qk = qkm1 * z - qkm2 * yc;
  127. if (qk != 0) {
  128. r = pk / qk;
  129. t = fabs((ans - r) / r);
  130. ans = r;
  131. } else
  132. t = 1.0;
  133. pkm2 = pkm1;
  134. pkm1 = pk;
  135. qkm2 = qkm1;
  136. qkm1 = qk;
  137. if (fabs(pk) > big) {
  138. pkm2 *= biginv;
  139. pkm1 *= biginv;
  140. qkm2 *= biginv;
  141. qkm1 *= biginv;
  142. }
  143. } while (t > MACHEP);
  144. return (ans * ax);
  145. }
  146. /* left tail of incomplete gamma function:
  147. *
  148. * inf. k
  149. * a -x - x
  150. * x e > ----------
  151. * - -
  152. * k=0 | (a+k+1)
  153. *
  154. */
  155. double igam(a, x) double a, x;
  156. {
  157. double ans, ax, c, r;
  158. if ((x <= 0) || (a <= 0))
  159. return (0.0);
  160. if ((x > 1.0) && (x > a))
  161. return (1.0 - igamc(a, x));
  162. /* Compute x**a * exp(-x) / gamma(a) */
  163. ax = a * log(x) - x - lgam(a);
  164. if (ax < -MAXLOG) {
  165. mtherr("igam", UNDERFLOW);
  166. return (0.0);
  167. }
  168. ax = exp(ax);
  169. /* power series */
  170. r = a;
  171. c = 1.0;
  172. ans = 1.0;
  173. do {
  174. r += 1.0;
  175. c *= x / r;
  176. ans += c;
  177. } while (c / ans > MACHEP);
  178. return (ans * ax / a);
  179. }