hyperg.c 7.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386
  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 <math.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)
  82. double a, b, x;
  83. {
  84. double asum, psum, acanc, pcanc, temp;
  85. /* See if a Kummer transformation will help */
  86. temp = b - a;
  87. if( fabs(temp) < 0.001 * fabs(a) )
  88. return( exp(x) * hyperg( temp, b, -x ) );
  89. psum = hy1f1p( a, b, x, &pcanc );
  90. if( pcanc < 1.0e-15 )
  91. goto done;
  92. /* try asymptotic series */
  93. asum = hy1f1a( a, b, x, &acanc );
  94. /* Pick the result with less estimated error */
  95. if( acanc < pcanc )
  96. {
  97. pcanc = acanc;
  98. psum = asum;
  99. }
  100. done:
  101. if( pcanc > 1.0e-12 )
  102. mtherr( "hyperg", PLOSS );
  103. return( psum );
  104. }
  105. /* Power series summation for confluent hypergeometric function */
  106. static double hy1f1p( a, b, x, err )
  107. double a, b, x;
  108. double *err;
  109. {
  110. double n, a0, sum, t, u, temp;
  111. double an, bn, maxt, pcanc;
  112. /* set up for power series summation */
  113. an = a;
  114. bn = b;
  115. a0 = 1.0;
  116. sum = 1.0;
  117. n = 1.0;
  118. t = 1.0;
  119. maxt = 0.0;
  120. while( t > MACHEP )
  121. {
  122. if( bn == 0 ) /* check bn first since if both */
  123. {
  124. mtherr( "hyperg", SING );
  125. return( MAXNUM ); /* an and bn are zero it is */
  126. }
  127. if( an == 0 ) /* a singularity */
  128. return( sum );
  129. if( n > 200 )
  130. goto pdone;
  131. u = x * ( an / (bn * n) );
  132. /* check for blowup */
  133. temp = fabs(u);
  134. if( (temp > 1.0 ) && (maxt > (MAXNUM/temp)) )
  135. {
  136. pcanc = 1.0; /* estimate 100% error */
  137. goto blowup;
  138. }
  139. a0 *= u;
  140. sum += a0;
  141. t = fabs(a0);
  142. if( t > maxt )
  143. maxt = t;
  144. /*
  145. if( (maxt/fabs(sum)) > 1.0e17 )
  146. {
  147. pcanc = 1.0;
  148. goto blowup;
  149. }
  150. */
  151. an += 1.0;
  152. bn += 1.0;
  153. n += 1.0;
  154. }
  155. pdone:
  156. /* estimate error due to roundoff and cancellation */
  157. if( sum != 0.0 )
  158. maxt /= fabs(sum);
  159. maxt *= MACHEP; /* this way avoids multiply overflow */
  160. pcanc = fabs( MACHEP * n + maxt );
  161. blowup:
  162. *err = pcanc;
  163. return( sum );
  164. }
  165. /* hy1f1a() */
  166. /* asymptotic formula for hypergeometric function:
  167. *
  168. * ( -a
  169. * -- ( |z|
  170. * | (b) ( -------- 2f0( a, 1+a-b, -1/x )
  171. * ( --
  172. * ( | (b-a)
  173. *
  174. *
  175. * x a-b )
  176. * e |x| )
  177. * + -------- 2f0( b-a, 1-a, 1/x ) )
  178. * -- )
  179. * | (a) )
  180. */
  181. static double hy1f1a( a, b, x, err )
  182. double a, b, x;
  183. double *err;
  184. {
  185. double h1, h2, t, u, temp, acanc, asum, err1, err2;
  186. if( x == 0 )
  187. {
  188. acanc = 1.0;
  189. asum = MAXNUM;
  190. goto adone;
  191. }
  192. temp = log( fabs(x) );
  193. t = x + temp * (a-b);
  194. u = -temp * a;
  195. if( b > 0 )
  196. {
  197. temp = lgam(b);
  198. t += temp;
  199. u += temp;
  200. }
  201. h1 = hyp2f0( a, a-b+1, -1.0/x, 1, &err1 );
  202. temp = exp(u) / gamma(b-a);
  203. h1 *= temp;
  204. err1 *= temp;
  205. h2 = hyp2f0( b-a, 1.0-a, 1.0/x, 2, &err2 );
  206. if( a < 0 )
  207. temp = exp(t) / gamma(a);
  208. else
  209. temp = exp( t - lgam(a) );
  210. h2 *= temp;
  211. err2 *= temp;
  212. if( x < 0.0 )
  213. asum = h1;
  214. else
  215. asum = h2;
  216. acanc = fabs(err1) + fabs(err2);
  217. if( b < 0 )
  218. {
  219. temp = gamma(b);
  220. asum *= temp;
  221. acanc *= fabs(temp);
  222. }
  223. if( asum != 0.0 )
  224. acanc /= fabs(asum);
  225. acanc *= 30.0; /* fudge factor, since error of asymptotic formula
  226. * often seems this much larger than advertised */
  227. adone:
  228. *err = acanc;
  229. return( asum );
  230. }
  231. /* hyp2f0() */
  232. double hyp2f0( a, b, x, type, err )
  233. double a, b, x;
  234. int type; /* determines what converging factor to use */
  235. double *err;
  236. {
  237. double a0, alast, t, tlast, maxt;
  238. double n, an, bn, u, sum, temp;
  239. an = a;
  240. bn = b;
  241. a0 = 1.0e0;
  242. alast = 1.0e0;
  243. sum = 0.0;
  244. n = 1.0e0;
  245. t = 1.0e0;
  246. tlast = 1.0e9;
  247. maxt = 0.0;
  248. do
  249. {
  250. if( an == 0 )
  251. goto pdone;
  252. if( bn == 0 )
  253. goto pdone;
  254. u = an * (bn * x / n);
  255. /* check for blowup */
  256. temp = fabs(u);
  257. if( (temp > 1.0 ) && (maxt > (MAXNUM/temp)) )
  258. goto error;
  259. a0 *= u;
  260. t = fabs(a0);
  261. /* terminating condition for asymptotic series */
  262. if( t > tlast )
  263. goto ndone;
  264. tlast = t;
  265. sum += alast; /* the sum is one term behind */
  266. alast = a0;
  267. if( n > 200 )
  268. goto ndone;
  269. an += 1.0e0;
  270. bn += 1.0e0;
  271. n += 1.0e0;
  272. if( t > maxt )
  273. maxt = t;
  274. }
  275. while( t > MACHEP );
  276. pdone: /* series converged! */
  277. /* estimate error due to roundoff and cancellation */
  278. *err = fabs( MACHEP * (n + maxt) );
  279. alast = a0;
  280. goto done;
  281. ndone: /* series did not converge */
  282. /* The following "Converging factors" are supposed to improve accuracy,
  283. * but do not actually seem to accomplish very much. */
  284. n -= 1.0;
  285. x = 1.0/x;
  286. switch( type ) /* "type" given as subroutine argument */
  287. {
  288. case 1:
  289. alast *= ( 0.5 + (0.125 + 0.25*b - 0.5*a + 0.25*x - 0.25*n)/x );
  290. break;
  291. case 2:
  292. alast *= 2.0/3.0 - b + 2.0*a + x - n;
  293. break;
  294. default:
  295. ;
  296. }
  297. /* estimate error due to roundoff, cancellation, and nonconvergence */
  298. *err = MACHEP * (n + maxt) + fabs ( a0 );
  299. done:
  300. sum += alast;
  301. return( sum );
  302. /* series blew up: */
  303. error:
  304. *err = MAXNUM;
  305. mtherr( "hyperg", TLOSS );
  306. return( sum );
  307. }