hyperg.c 7.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357
  1. /* hyperg.c
  2. *
  3. * Confluent hypergeometric function
  4. *
  5. *
  6. *
  7. * SYNOPSIS:
  8. *
  9. * double a, b, x, y, hyperg();
  10. *
  11. * y = hyperg( a, b, x );
  12. *
  13. *
  14. *
  15. * DESCRIPTION:
  16. *
  17. * Computes the confluent hypergeometric function
  18. *
  19. * 1 2
  20. * a x a(a+1) x
  21. * F ( a,b;x ) = 1 + ---- + --------- + ...
  22. * 1 1 b 1! b(b+1) 2!
  23. *
  24. * Many higher transcendental functions are special cases of
  25. * this power series.
  26. *
  27. * As is evident from the formula, b must not be a negative
  28. * integer or zero unless a is an integer with 0 >= a > b.
  29. *
  30. * The routine attempts both a direct summation of the series
  31. * and an asymptotic expansion. In each case error due to
  32. * roundoff, cancellation, and nonconvergence is estimated.
  33. * The result with smaller estimated error is returned.
  34. *
  35. *
  36. *
  37. * ACCURACY:
  38. *
  39. * Tested at random points (a, b, x), all three variables
  40. * ranging from 0 to 30.
  41. * Relative error:
  42. * arithmetic domain # trials peak rms
  43. * DEC 0,30 2000 1.2e-15 1.3e-16
  44. qtst1:
  45. 21800 max = 1.4200E-14 rms = 1.0841E-15 ave = -5.3640E-17
  46. ltstd:
  47. 25500 max = 1.2759e-14 rms = 3.7155e-16 ave = 1.5384e-18
  48. * IEEE 0,30 30000 1.8e-14 1.1e-15
  49. *
  50. * Larger errors can be observed when b is near a negative
  51. * integer or zero. Certain combinations of arguments yield
  52. * serious cancellation error in the power series summation
  53. * and also are not in the region of near convergence of the
  54. * asymptotic series. An error message is printed if the
  55. * self-estimated relative error is greater than 1.0e-12.
  56. *
  57. */
  58. /* hyperg.c */
  59. /*
  60. Cephes Math Library Release 2.8: June, 2000
  61. Copyright 1984, 1987, 1988, 2000 by Stephen L. Moshier
  62. */
  63. #include "mconf.h"
  64. #ifdef ANSIPROT
  65. extern double exp(double);
  66. extern double log(double);
  67. extern double gamma(double);
  68. extern double lgam(double);
  69. extern double fabs(double);
  70. double hyp2f0(double, double, double, int, double *);
  71. static double hy1f1p(double, double, double, double *);
  72. static double hy1f1a(double, double, double, double *);
  73. double hyperg(double, double, double);
  74. #else
  75. double exp(), log(), gamma(), lgam(), fabs(), hyp2f0();
  76. static double hy1f1p();
  77. static double hy1f1a();
  78. double hyperg();
  79. #endif
  80. extern double MAXNUM, MACHEP;
  81. double hyperg(a, b, x) double a, b, x;
  82. {
  83. double asum, psum, acanc, pcanc, temp;
  84. /* See if a Kummer transformation will help */
  85. temp = b - a;
  86. if (fabs(temp) < 0.001 * fabs(a))
  87. return (exp(x) * hyperg(temp, b, -x));
  88. psum = hy1f1p(a, b, x, &pcanc);
  89. if (pcanc < 1.0e-15)
  90. goto done;
  91. /* try asymptotic series */
  92. asum = hy1f1a(a, b, x, &acanc);
  93. /* Pick the result with less estimated error */
  94. if (acanc < pcanc) {
  95. pcanc = acanc;
  96. psum = asum;
  97. }
  98. done:
  99. if (pcanc > 1.0e-12)
  100. mtherr("hyperg", PLOSS);
  101. return (psum);
  102. }
  103. /* Power series summation for confluent hypergeometric function */
  104. static double hy1f1p(a, b, x, err) double a, b, x;
  105. double *err;
  106. {
  107. double n, a0, sum, t, u, temp;
  108. double an, bn, maxt, pcanc;
  109. /* set up for power series summation */
  110. an = a;
  111. bn = b;
  112. a0 = 1.0;
  113. sum = 1.0;
  114. n = 1.0;
  115. t = 1.0;
  116. maxt = 0.0;
  117. while (t > MACHEP) {
  118. if (bn == 0) /* check bn first since if both */
  119. {
  120. mtherr("hyperg", SING);
  121. return (MAXNUM); /* an and bn are zero it is */
  122. }
  123. if (an == 0) /* a singularity */
  124. return (sum);
  125. if (n > 200)
  126. goto pdone;
  127. u = x * (an / (bn * n));
  128. /* check for blowup */
  129. temp = fabs(u);
  130. if ((temp > 1.0) && (maxt > (MAXNUM / temp))) {
  131. pcanc = 1.0; /* estimate 100% error */
  132. goto blowup;
  133. }
  134. a0 *= u;
  135. sum += a0;
  136. t = fabs(a0);
  137. if (t > maxt)
  138. maxt = t;
  139. /*
  140. if( (maxt/fabs(sum)) > 1.0e17 )
  141. {
  142. pcanc = 1.0;
  143. goto blowup;
  144. }
  145. */
  146. an += 1.0;
  147. bn += 1.0;
  148. n += 1.0;
  149. }
  150. pdone:
  151. /* estimate error due to roundoff and cancellation */
  152. if (sum != 0.0)
  153. maxt /= fabs(sum);
  154. maxt *= MACHEP; /* this way avoids multiply overflow */
  155. pcanc = fabs(MACHEP * n + maxt);
  156. blowup:
  157. *err = pcanc;
  158. return (sum);
  159. }
  160. /* hy1f1a() */
  161. /* asymptotic formula for hypergeometric function:
  162. *
  163. * ( -a
  164. * -- ( |z|
  165. * | (b) ( -------- 2f0( a, 1+a-b, -1/x )
  166. * ( --
  167. * ( | (b-a)
  168. *
  169. *
  170. * x a-b )
  171. * e |x| )
  172. * + -------- 2f0( b-a, 1-a, 1/x ) )
  173. * -- )
  174. * | (a) )
  175. */
  176. static double hy1f1a(a, b, x, err) double a, b, x;
  177. double *err;
  178. {
  179. double h1, h2, t, u, temp, acanc, asum, err1, err2;
  180. if (x == 0) {
  181. acanc = 1.0;
  182. asum = MAXNUM;
  183. goto adone;
  184. }
  185. temp = log(fabs(x));
  186. t = x + temp * (a - b);
  187. u = -temp * a;
  188. if (b > 0) {
  189. temp = lgam(b);
  190. t += temp;
  191. u += temp;
  192. }
  193. h1 = hyp2f0(a, a - b + 1, -1.0 / x, 1, &err1);
  194. temp = exp(u) / gamma(b - a);
  195. h1 *= temp;
  196. err1 *= temp;
  197. h2 = hyp2f0(b - a, 1.0 - a, 1.0 / x, 2, &err2);
  198. if (a < 0)
  199. temp = exp(t) / gamma(a);
  200. else
  201. temp = exp(t - lgam(a));
  202. h2 *= temp;
  203. err2 *= temp;
  204. if (x < 0.0)
  205. asum = h1;
  206. else
  207. asum = h2;
  208. acanc = fabs(err1) + fabs(err2);
  209. if (b < 0) {
  210. temp = gamma(b);
  211. asum *= temp;
  212. acanc *= fabs(temp);
  213. }
  214. if (asum != 0.0)
  215. acanc /= fabs(asum);
  216. acanc *= 30.0; /* fudge factor, since error of asymptotic formula
  217. * often seems this much larger than advertised */
  218. adone:
  219. *err = acanc;
  220. return (asum);
  221. }
  222. /* hyp2f0() */
  223. double hyp2f0(a, b, x, type, err) double a, b, x;
  224. int type; /* determines what converging factor to use */
  225. double *err;
  226. {
  227. double a0, alast, t, tlast, maxt;
  228. double n, an, bn, u, sum, temp;
  229. an = a;
  230. bn = b;
  231. a0 = 1.0e0;
  232. alast = 1.0e0;
  233. sum = 0.0;
  234. n = 1.0e0;
  235. t = 1.0e0;
  236. tlast = 1.0e9;
  237. maxt = 0.0;
  238. do {
  239. if (an == 0)
  240. goto pdone;
  241. if (bn == 0)
  242. goto pdone;
  243. u = an * (bn * x / n);
  244. /* check for blowup */
  245. temp = fabs(u);
  246. if ((temp > 1.0) && (maxt > (MAXNUM / temp)))
  247. goto error;
  248. a0 *= u;
  249. t = fabs(a0);
  250. /* terminating condition for asymptotic series */
  251. if (t > tlast)
  252. goto ndone;
  253. tlast = t;
  254. sum += alast; /* the sum is one term behind */
  255. alast = a0;
  256. if (n > 200)
  257. goto ndone;
  258. an += 1.0e0;
  259. bn += 1.0e0;
  260. n += 1.0e0;
  261. if (t > maxt)
  262. maxt = t;
  263. } while (t > MACHEP);
  264. pdone: /* series converged! */
  265. /* estimate error due to roundoff and cancellation */
  266. *err = fabs(MACHEP * (n + maxt));
  267. alast = a0;
  268. goto done;
  269. ndone: /* series did not converge */
  270. /* The following "Converging factors" are supposed to improve accuracy,
  271. * but do not actually seem to accomplish very much. */
  272. n -= 1.0;
  273. x = 1.0 / x;
  274. switch (type) /* "type" given as subroutine argument */
  275. {
  276. case 1:
  277. alast *= (0.5 + (0.125 + 0.25 * b - 0.5 * a + 0.25 * x - 0.25 * n) / x);
  278. break;
  279. case 2:
  280. alast *= 2.0 / 3.0 - b + 2.0 * a + x - n;
  281. break;
  282. default:;
  283. }
  284. /* estimate error due to roundoff, cancellation, and nonconvergence */
  285. *err = MACHEP * (n + maxt) + fabs(a0);
  286. done:
  287. sum += alast;
  288. return (sum);
  289. /* series blew up: */
  290. error:
  291. *err = MAXNUM;
  292. mtherr("hyperg", TLOSS);
  293. return (sum);
  294. }