specfun.c 19 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863
  1. /* This file is part of the GNU plotutils package. */
  2. /* Copyright (C) 1997, 1998, 1998, 2005, 2008, 2009, Free Software
  3. Foundation, Inc. */
  4. /* This file contains the five functions
  5. ibeta, igamma, norm, invnorm, inverf
  6. and versions of lgamma, erf, erfc for machines without them. */
  7. /* The inspiration for this file was the file specfun.c that has long been
  8. a part of the gnuplot distribution. However, the current version of the
  9. present file has been rewritten, largely from scratch, on the basis of
  10. published algorithms. It is not dependent on the gnuplot specfun.c
  11. (which has included, and still includes, copyrighted code). */
  12. #include "sys-defines.h"
  13. #include "ode.h"
  14. #include "extern.h"
  15. #include <errno.h>
  16. #define ITERMAX 200
  17. #ifdef FLT_EPSILON
  18. #define MACHEPS FLT_EPSILON /* 1.0E-08 */
  19. #else
  20. #define MACHEPS 1.0E-08
  21. #endif
  22. #ifdef FLT_MIN_EXP
  23. #define MINEXP FLT_MIN_EXP /* -88.0 */
  24. #else
  25. #define MINEXP -88.0
  26. #endif
  27. #ifdef FLT_MAX_EXP
  28. #define MAXEXP FLT_MAX_EXP /* +88.0 */
  29. #else
  30. #define MAXEXP 88.0
  31. #endif
  32. #ifdef FLT_MAX
  33. #define OFLOW FLT_MAX /* 1.0E+37 */
  34. #else
  35. #define OFLOW 1.0E+37
  36. #endif
  37. #ifdef FLT_MAX_10_EXP
  38. #define XBIG FLT_MAX_10_EXP /* 2.55E+305 */
  39. #else
  40. #define XBIG 2.55E+305
  41. #endif
  42. #ifndef HUGE_VAL
  43. #ifdef HUGE
  44. #define HUGE_VAL HUGE
  45. #else
  46. #ifdef INF
  47. #define HUGE_VAL INF
  48. #else
  49. #define HUGE_VAL OFLOW
  50. #endif
  51. #endif
  52. #endif
  53. /*
  54. * Mathematical constants
  55. */
  56. #ifdef M_PI
  57. #undef M_PI
  58. #endif
  59. #define M_PI 3.14159265358979323846264338327950288
  60. #ifdef M_SQRT2
  61. #undef M_SQRT2
  62. #endif
  63. #define M_SQRT2 1.41421356237309504880168872420969809
  64. #define M_LNSQRT2PI 0.9189385332046727
  65. /* Forward references */
  66. /* The following gamma-related nonsense is necessary because (1) some
  67. vendors have lgamma(), some have gamma(), and some have neither [see
  68. include/sys-defines.h for further comments], (2) some vendors do not
  69. declare whichever function they have [e.g. Irix 5.3 requires an
  70. auxiliary preprocessing symbol to be defined for the declaration in
  71. math.h to be visible], and (3) some vendors supply broken versions which
  72. we can't use [e.g. AIX's libm.a gamma support is conspicuously broken],
  73. so we need to link in a replacement, but we can't use the same name for
  74. the external symbol `signgam'. What a mess! -- rsm */
  75. #ifdef NO_SYSTEM_GAMMA
  76. #define SIGNGAM our_signgam
  77. static int SIGNGAM;
  78. double f_lgamma (double x);
  79. static double lgamma_neg (double x);
  80. static double lgamma_pos (double x);
  81. #else /* not NO_SYSTEM_GAMMA, we link in vendor code */
  82. #define SIGNGAM signgam
  83. extern int SIGNGAM;
  84. #endif
  85. double f_gamma (double x);
  86. #ifndef HAVE_ERF
  87. double erf (double x);
  88. double erfc (double x);
  89. #endif
  90. double ibeta (double a, double b, double x);
  91. double igamma (double a, double x);
  92. double inverf (double x);
  93. double invnorm (double x);
  94. double norm (double x);
  95. static double ibeta_internal (double a, double b, double x);
  96. /*****************************************************************/
  97. /************ Functions related to gamma function ****************/
  98. /*****************************************************************/
  99. /* Our gamma function. F_LGAMMA(), which this calls, computes the log of
  100. the gamma function, with the sign being returned in SIGNGAM. F_LGAMMA()
  101. is defined in include/sys-defines.h. It may be a vendor-supplied
  102. lgamma(), a vendor-supplied gamma(), or our own f_lgamma (see below). */
  103. double
  104. f_gamma (double x)
  105. {
  106. #ifdef HAVE_MATHERR
  107. #ifdef __cplusplus
  108. struct __exception exc;
  109. #else
  110. struct exception exc;
  111. #endif
  112. #endif
  113. double y = F_LGAMMA(x);
  114. if (y > MAXEXP)
  115. {
  116. #ifdef HAVE_MATHERR
  117. exc.name = (char *)"gamma";
  118. exc.arg1 = x;
  119. exc.retval = HUGE_VAL;
  120. exc.type = OVERFLOW;
  121. if (!matherr (&exc))
  122. {
  123. fprintf (stderr, "gamma: OVERFLOW error\n");
  124. errno = ERANGE;
  125. }
  126. return exc.retval;
  127. #else
  128. errno = ERANGE;
  129. return HUGE_VAL;
  130. #endif
  131. }
  132. else
  133. return SIGNGAM * exp (y);
  134. }
  135. #ifdef NO_SYSTEM_GAMMA
  136. /*
  137. Define our own lgamma(): compute log(Gamma(z)) for positive z as the sum
  138. of a suitably generated and truncated Lanczos series. Lanczos series
  139. resemble Stirling (i.e. DeMoivre) asymptotic approximations, but unlike
  140. them are not asymptotic series; rather, they are convergent in the entire
  141. right half plane, on which a uniform bound on the error can be obtained.
  142. See C. Lanczos, "A Precision Approximation of the Gamma Function", SIAM
  143. J. Numerical Analysis 1B (1964), 86--96. */
  144. double
  145. f_lgamma (double z)
  146. {
  147. SIGNGAM = 1; /* will return sign of Gamma(z) in SIGNGAM */
  148. if (z <= 0.0)
  149. return lgamma_neg (z);
  150. else
  151. return lgamma_pos (z);
  152. }
  153. /* Case I. z<=0 (if z < 0, reducible to z>0 case by reflection formula) */
  154. static double
  155. lgamma_neg (double z)
  156. {
  157. double intpart, trigfac, retval;
  158. #ifdef HAVE_MATHERR
  159. #ifdef __cplusplus
  160. struct __exception exc;
  161. #else
  162. struct exception exc;
  163. #endif
  164. #endif
  165. if (modf (-z, &intpart) == 0.0)
  166. /* z is nonpositive integer, so SING error */
  167. {
  168. #ifdef HAVE_MATHERR
  169. exc.name = "lgamma";
  170. exc.arg1 = z;
  171. exc.retval = HUGE_VAL;
  172. exc.type = SING;
  173. if (!matherr (&exc))
  174. {
  175. fprintf (stderr, "lgamma: SING error\n");
  176. errno = EDOM;
  177. }
  178. return (exc.retval);
  179. #else
  180. errno = EDOM;
  181. return HUGE_VAL;
  182. #endif /* HAVE_MATHERR */
  183. }
  184. /* use Euler's reflection formula, and call lgamma_pos() */
  185. trigfac = sin (M_PI * z) / M_PI;
  186. if (trigfac < 0.0)
  187. {
  188. trigfac = - trigfac;
  189. SIGNGAM = -1;
  190. }
  191. retval = - lgamma_pos (1.0 - z) - log (trigfac);
  192. if (fabs (retval) == HUGE_VAL)
  193. {
  194. #ifdef HAVE_MATHERR
  195. exc.name = "lgamma";
  196. exc.arg1 = z;
  197. exc.retval = HUGE_VAL;
  198. exc.type = OVERFLOW;
  199. if (!matherr(&exc))
  200. errno = ERANGE;
  201. return (exc.retval);
  202. #else
  203. errno = ERANGE;
  204. return HUGE_VAL;
  205. #endif
  206. }
  207. return retval;
  208. }
  209. /* Case II. z>0, the primary case */
  210. /* Lanczos parameter G (called lower-case gamma by him). "[A] large value
  211. of G is advocated if very high accuracy is demanded, but then the
  212. required number of terms will also be larger." */
  213. #define LANCZOS_G 6
  214. /* Values for coeffs of as many terms in the Lanczos expansion as are
  215. needed for this value of G, computed by Ray Toy <toy@rtp.ericsson.se>.
  216. (In his 1964 paper, Lanczos only went up to G=5.) It is claimed (see
  217. gnuplot's specfun.c) that this value of G (i.e., 6) and number of terms
  218. will yield 14-digit accuracy everywhere except near z=1 and z=2. */
  219. #define NUM_LANCZOS_TERMS 9
  220. static const double lanczos[NUM_LANCZOS_TERMS] =
  221. {
  222. .99999999999980993227684700473478296744476168282198,
  223. 676.52036812188509856700919044401903816411251975244084,
  224. -1259.13921672240287047156078755282840836424300664868028,
  225. 771.32342877765307884865282588943070775227268469602500,
  226. -176.61502916214059906584551353999392943274507608117860,
  227. 12.50734327868690481445893685327104972970563021816420,
  228. -.13857109526572011689554706984971501358032683492780,
  229. .00000998436957801957085956266828104544089848531228,
  230. .00000015056327351493115583383579667028994545044040
  231. };
  232. static double
  233. lgamma_pos (double z)
  234. {
  235. double accum, retval;
  236. int i;
  237. #ifdef HAVE_MATHERR
  238. #ifdef __cplusplus
  239. struct __exception exc;
  240. #else
  241. struct exception exc;
  242. #endif
  243. #endif
  244. accum = lanczos[0];
  245. for (i = 1; i < NUM_LANCZOS_TERMS; i++)
  246. accum += lanczos[i] / (z + i - 1);
  247. retval = (log (accum) + M_LNSQRT2PI - z - LANCZOS_G - 0.5
  248. + (z - 0.5) * log (z + LANCZOS_G + 0.5));
  249. if (retval == HUGE_VAL)
  250. {
  251. #ifdef HAVE_MATHERR
  252. exc.name = "lgamma";
  253. exc.arg1 = z;
  254. exc.retval = HUGE_VAL;
  255. exc.type = OVERFLOW;
  256. if (!matherr (&exc))
  257. {
  258. fprintf (stderr, "lgamma: OVERFLOW error\n");
  259. errno = ERANGE;
  260. }
  261. return exc.retval;
  262. #else
  263. errno = ERANGE;
  264. return HUGE_VAL;
  265. #endif
  266. }
  267. return retval;
  268. }
  269. #endif /* NO_SYSTEM_GAMMA */
  270. /*****************************************************************/
  271. /************ Functions related to inverse beta function *********/
  272. /*****************************************************************/
  273. /* Our incomplete beta function, I_x(a,b). Here a,b>0 and x is in [0,1].
  274. Returned value is in [0,1]. Note: this normalization convention is not
  275. universal.
  276. The formula given in Abramowitz & Stegun (Eq. 26.5.8) is used. It
  277. includes a continued fraction expansion. They say, "Best results are
  278. obtained when x < (a-1)/(a+b-2)." We use it when x <= a/(a+b).
  279. This calls F_LGAMMA (e.g., lgamma()) to compute the prefactor in the
  280. formula. */
  281. double
  282. ibeta (double a, double b, double x)
  283. {
  284. double retval;
  285. #ifdef HAVE_MATHERR
  286. #ifdef __cplusplus
  287. struct __exception exc;
  288. #else
  289. struct exception exc;
  290. #endif
  291. #endif
  292. if (x < 0.0 || x > 1.0 || a <= 0.0 || b <= 0.0) /* DOMAIN error */
  293. {
  294. #ifdef HAVE_MATHERR
  295. exc.name = (char *)"ibeta";
  296. exc.arg1 = a;
  297. exc.arg2 = b; /* have no arg3, can't return x (!) */
  298. exc.retval = HUGE_VAL;
  299. exc.type = DOMAIN;
  300. if (!matherr (&exc))
  301. {
  302. fprintf (stderr, "ibeta: DOMAIN error\n");
  303. errno = EDOM;
  304. }
  305. return exc.retval;
  306. #else
  307. errno = EDOM;
  308. return HUGE_VAL;
  309. #endif
  310. }
  311. if (x == 0.0 || x == 1.0)
  312. return x;
  313. if (a < x * (a + b))
  314. /* interchange */
  315. retval = 1.0 - ibeta_internal (b, a, 1.0 - x);
  316. else
  317. retval = ibeta_internal (a, b, x);
  318. if (retval < 0.0) /* error: failure of convergence */
  319. {
  320. #ifdef HAVE_MATHERR
  321. exc.name = (char *)"ibeta";
  322. exc.arg1 = a;
  323. exc.arg2 = b; /* have no arg3, can't return x (!) */
  324. exc.retval = HUGE_VAL;
  325. exc.type = TLOSS;
  326. if (!matherr (&exc))
  327. {
  328. fprintf (stderr, "ibeta: TLOSS error\n");
  329. errno = EDOM;
  330. }
  331. return exc.retval;
  332. #else
  333. errno = EDOM;
  334. return HUGE_VAL;
  335. #endif
  336. }
  337. return retval;
  338. }
  339. /* Evaluate convergents of the continued fraction by Wallis's method;
  340. return value will be positive, except that -1.0 is returned if there is
  341. no convergence. */
  342. static double
  343. ibeta_internal (double a, double b, double x)
  344. {
  345. double A0, B0;
  346. double A2 = 1.0;
  347. double B2 = 0.0;
  348. double A1 = 1.0;
  349. double B1 = 1.0;
  350. double prefactor;
  351. double f0 = 0.0, f1 = 1.0; /* f0 initted to quiet compiler */
  352. int goodf0, goodf1 = 1;
  353. int j;
  354. prefactor = exp (a * log (x) + b * log (1.0 - x)
  355. + F_LGAMMA(a + b) - F_LGAMMA(a + 1.0) - F_LGAMMA(b));
  356. for (j = 1; j <= ITERMAX; j++)
  357. {
  358. double aj;
  359. int m;
  360. if (j % 2) /* j odd, j = 2m + 1 */
  361. {
  362. m = (j - 1)/2;
  363. aj = - (a + m) * (a + b + m) * x / ((a + 2 * m) * (a + 2 * m + 1));
  364. }
  365. else /* j even, j = 2m */
  366. {
  367. m = j/2;
  368. aj = m * (b - m) * x / ((a + 2 * m - 1) * (a + 2 * m));
  369. }
  370. A0 = 1.0 * A1 + aj * A2;
  371. B0 = 1.0 * B1 + aj * B2;
  372. if (B0 != 0.0)
  373. {
  374. double ren;
  375. /* renormalize; don't really need to do this on each pass */
  376. ren = 1.0 / B0;
  377. A0 *= ren;
  378. B0 = 1.0;
  379. A1 *= ren;
  380. B1 *= ren;
  381. f0 = A0;
  382. goodf0 = 1;
  383. /* test f0 = A0/B0 = A0 for exit */
  384. if (goodf1 && fabs (f0 - f1) <= DMIN(MACHEPS, fabs (f0) * MACHEPS))
  385. return (prefactor / f0);
  386. }
  387. else
  388. goodf0 = 0;
  389. /* shift down */
  390. A2 = A1;
  391. B2 = B1;
  392. A1 = A0;
  393. B1 = B0;
  394. f1 = f0;
  395. goodf1 = goodf0;
  396. }
  397. /* if we reached here, convergence failed */
  398. return -1.0;
  399. }
  400. /*****************************************************************/
  401. /************ Functions related to incomplete gamma function *****/
  402. /*****************************************************************/
  403. /* Our incomplete gamma function, igamma(a,x) with a>0.0, x>=0.0. Return
  404. value is in [0,1]. The algorithm is AS 239, documented in B. L. Shea,
  405. "Chi-Squared and Incomplete Gamma Integral", Applied Statistics 37
  406. (1988), 466-473.
  407. There have been claims that if 0<=x<=1, in which case Shea's algorithm
  408. uses Pearson's series rather than a continued fraction representation,
  409. an inaccurate value may result. This has not been verified. There have
  410. also been claims that the continued fraction representation is reliable
  411. only if x >= a+2, rather than x >= a (the latter being Shea's
  412. condition). For safety, we use it only if x >= a+2. */
  413. double
  414. igamma (double a, double x)
  415. {
  416. double arg, prefactor;
  417. int i;
  418. #ifdef HAVE_MATHERR
  419. #ifdef __cplusplus
  420. struct __exception exc;
  421. #else
  422. struct exception exc;
  423. #endif
  424. #endif
  425. if (x < 0.0 || a <= 0.0) /* DOMAIN error */
  426. {
  427. #ifdef HAVE_MATHERR
  428. exc.name = (char *)"igamma";
  429. exc.arg1 = a;
  430. exc.arg2 = x;
  431. exc.retval = HUGE_VAL;
  432. exc.type = DOMAIN;
  433. if (!matherr (&exc))
  434. {
  435. fprintf (stderr, "igamma: DOMAIN error\n");
  436. errno = EDOM;
  437. }
  438. return exc.retval;
  439. #else
  440. errno = EDOM;
  441. return HUGE_VAL;
  442. #endif
  443. }
  444. if (x > XBIG) /* TLOSS error */
  445. {
  446. #ifdef HAVE_MATHERR
  447. exc.name = (char *)"igamma";
  448. exc.arg1 = a;
  449. exc.arg2 = x;
  450. exc.retval = 1.0;
  451. exc.type = TLOSS;
  452. if (!matherr (&exc))
  453. {
  454. fprintf (stderr, "igamma: TLOSS error\n");
  455. errno = EDOM;
  456. }
  457. return exc.retval;
  458. #else
  459. errno = EDOM;
  460. return 1.0;
  461. #endif
  462. }
  463. if (x == 0.0)
  464. return 0.0;
  465. /* check exponentiation in prefactor */
  466. arg = a * log (x) - x - F_LGAMMA(a + 1.0);
  467. if (arg < MINEXP)
  468. {
  469. #ifdef HAVE_MATHERR
  470. exc.name = (char *)"igamma";
  471. exc.arg1 = a;
  472. exc.arg2 = x;
  473. exc.retval = 0.0;
  474. exc.type = TLOSS;
  475. if (!matherr (&exc))
  476. {
  477. fprintf (stderr, "ibeta: TLOSS error\n");
  478. errno = EDOM;
  479. }
  480. return exc.retval;
  481. #else
  482. errno = EDOM;
  483. return 0.0;
  484. #endif
  485. }
  486. prefactor = exp (arg);
  487. if ((x > 1.0) && (x >= a + 2.0))
  488. /* use the continued fraction, not Pearson's series; generate its
  489. convergents by Wallis's method */
  490. {
  491. double A0, B0, A1, B1, A2, B2;
  492. double f0 = 0.0, f1; /* f0 initted to quiet compiler */
  493. double aa, bb;
  494. int goodf0, goodf1 = 1;
  495. aa = 1.0 - a;
  496. bb = aa + x + 1.0;
  497. A2 = 1.0;
  498. B2 = x;
  499. A1 = x + 1.0;
  500. B1 = x * bb;
  501. f1 = A1 / B1;
  502. for (i = 1; i <= ITERMAX; i++)
  503. {
  504. aa++;
  505. bb += 2.0;
  506. A0 = bb * A1 - i * aa * A2;
  507. B0 = bb * B1 - i * aa * B2;
  508. if (B0 != 0.0)
  509. {
  510. f0 = A0 / B0;
  511. if (goodf1 &&
  512. fabs (f0 - f1) <= DMIN(MACHEPS, fabs (f0) * MACHEPS))
  513. return (1.0 - prefactor * a * f0);
  514. goodf0 = 1;
  515. }
  516. else
  517. goodf0 = 0;
  518. /* shift down */
  519. A2 = A1;
  520. B2 = B1;
  521. A1 = A0;
  522. B1 = B0;
  523. f1 = f0;
  524. goodf1 = goodf0;
  525. if (fabs(A0) >= OFLOW)
  526. /* renormalize */
  527. {
  528. A2 /= OFLOW;
  529. B2 /= OFLOW;
  530. A1 /= OFLOW;
  531. B1 /= OFLOW;
  532. }
  533. }
  534. }
  535. else
  536. /* use Pearson's series, not the continued fraction */
  537. {
  538. double aa, bb, cc;
  539. aa = a;
  540. bb = 1.0;
  541. cc = 1.0;
  542. for (i = 0; i <= ITERMAX; i++)
  543. {
  544. aa++;
  545. cc *= (x / aa);
  546. bb += cc;
  547. if (cc < bb * MACHEPS)
  548. return prefactor * bb;
  549. }
  550. }
  551. /* if we reached here, convergence failed */
  552. #ifdef HAVE_MATHERR
  553. exc.name = (char *)"igamma";
  554. exc.arg1 = a;
  555. exc.arg2 = x;
  556. exc.retval = HUGE_VAL;
  557. exc.type = TLOSS;
  558. if (!matherr (&exc))
  559. {
  560. fprintf (stderr, "ibeta: TLOSS error\n");
  561. errno = EDOM;
  562. }
  563. return exc.retval;
  564. #else
  565. errno = EDOM;
  566. return HUGE_VAL;
  567. #endif
  568. }
  569. #ifndef HAVE_ERF
  570. double
  571. erf (double x)
  572. {
  573. return x < 0.0 ? -igamma (0.5, x * x) : igamma (0.5, x * x);
  574. }
  575. double
  576. erfc (double x)
  577. {
  578. return x < 0.0 ? 1.0 + igamma (0.5, x * x) : 1.0 - igamma (0.5, x * x);
  579. }
  580. #endif /* not HAVE_ERF */
  581. double
  582. norm (double x)
  583. {
  584. return 0.5 * (1.0 + erf (0.5 * M_SQRT2 * x));
  585. }
  586. /*****************************************************************/
  587. /************ Functions related to inverse error function ********/
  588. /*****************************************************************/
  589. /* Our inverse error function, inverf(x) with -1.0<x<1.0. It approximates
  590. this (odd) function by four distinct degree-3 rational functions, each
  591. applying in a sub-interval of the right half-interval 0.0<x<1.0. */
  592. /* rational function R0(x) */
  593. static const double n0[4] =
  594. {
  595. -6.075593, 9.577584, -4.026908, 0.3110567
  596. };
  597. static const double d0[4] =
  598. {
  599. -6.855572, 12.601905, -6.855985, 1.0
  600. };
  601. /* rational function R1(w) */
  602. static const double n1[4] =
  603. {
  604. -39.202359, 19.332985, -6.953050, 0.9360732
  605. };
  606. static const double d1[4] =
  607. {
  608. -44.27977, 21.98546, -7.586103, 1.0
  609. };
  610. /* rational function R2(w) */
  611. static const double n2[4] =
  612. {
  613. -5.911558, 4.795462, -3.111584, 1.005405
  614. };
  615. static const double d2[4] =
  616. {
  617. -6.266786, 4.666263, -2.962883, 1.0
  618. };
  619. /* rational function R3(1/w) */
  620. static const double n3[4] =
  621. {
  622. 0.09952975, 0.51914515, -0.2187214, 1.0107864
  623. };
  624. static const double d3[4] =
  625. {
  626. 0.09952975, 0.5211733, -0.06888301, 1.0
  627. };
  628. double
  629. inverf (double x)
  630. {
  631. static double num, den, retval;
  632. static int xsign;
  633. #ifdef HAVE_MATHERR
  634. #ifdef __cplusplus
  635. struct __exception exc;
  636. #else
  637. struct exception exc;
  638. #endif
  639. #endif
  640. if (x <= -1.0 || x >= 1.0) /* DOMAIN error */
  641. {
  642. #ifdef HAVE_MATHERR
  643. exc.name = (char *)"inverf";
  644. exc.arg1 = x;
  645. exc.retval = (x < 0.0 ? -HUGE_VAL : HUGE_VAL);
  646. exc.type = DOMAIN;
  647. if (!matherr (&exc))
  648. {
  649. fprintf (stderr, "inverf: DOMAIN error\n");
  650. errno = EDOM;
  651. }
  652. return exc.retval;
  653. #else
  654. errno = EDOM;
  655. return (x < 0.0 ? -HUGE_VAL : HUGE_VAL);
  656. #endif
  657. }
  658. /* exploit oddness in x */
  659. xsign = (x >= 0.0 ? 1 : -1);
  660. x = (xsign > 0 ? x : -x);
  661. /* N.B. The numerator and denominator of each of these four rational
  662. approximants should really be written in nested polynomial form. */
  663. if (x <= 0.85)
  664. /* 0.0 <= x <= 0.85; use f = x R0(x**2), where R0 is degree-3 rational */
  665. {
  666. double y;
  667. y = x * x;
  668. num = n0[0] + n0[1]*y + n0[2]*y*y + n0[3]*y*y*y;
  669. den = d0[0] + d0[1]*y + d0[2]*y*y + d0[3]*y*y*y;
  670. retval = x * num / den;
  671. }
  672. else /* x > 0.85 */
  673. {
  674. double w;
  675. w = sqrt (- log (1 - x * x)); /* w > 1.132 */
  676. /* note that as x->1-, i.e., w->infinity, retval is asymptotic to w,
  677. to leading order */
  678. if (w <= 2.5)
  679. /* 1.132 < w <= 2.5; use f = w R1(w), where R1 is degree-3 rational */
  680. {
  681. num = n1[0] + n1[1]*w + n1[2]*w*w + n1[3]*w*w*w;
  682. den = d1[0] + d1[1]*w + d1[2]*w*w + d1[3]*w*w*w;
  683. retval = w * num / den;
  684. }
  685. else if (w <= 4.0)
  686. /* 2.5 < w <= 4.0; use f = w R2(w), where R2 is degree-3 rational */
  687. {
  688. num = n2[0] + n2[1]*w + n2[2]*w*w + n2[3]*w*w*w;
  689. den = d2[0] + d2[1]*w + d2[2]*w*w + d2[3]*w*w*w;
  690. retval = w * num / den;
  691. }
  692. else
  693. /* w > 4.0; use f = w R3(1/w), where R3 is degree-3 rational
  694. with equal constant terms in numerator and denominator */
  695. {
  696. double w1;
  697. w1 = 1.0 / w;
  698. num = n3[0] + n3[1]*w1 + n3[2]*w1*w1 + n3[3]*w1*w1*w1;
  699. den = d3[0] + d3[1]*w1 + d3[2]*w1*w1 + d3[3]*w1*w1*w1;
  700. retval = w * num / den;
  701. }
  702. }
  703. return (xsign > 0 ? retval : -retval);
  704. }
  705. /* Our inverse normal function (i.e., inverse Gaussian probability
  706. function), invnorm(x) with 0.0<x<1.0. Trivially expressed in terms of
  707. inverf(), just as erf() can be expressed in terms of erfc(). */
  708. double
  709. invnorm (double x)
  710. {
  711. #ifdef HAVE_MATHERR
  712. #ifdef __cplusplus
  713. struct __exception exc;
  714. #else
  715. struct exception exc;
  716. #endif
  717. #endif
  718. if (x <= 0.0 || x >= 1.0) /* DOMAIN error */
  719. {
  720. #ifdef HAVE_MATHERR
  721. exc.name = (char *)"invnorm";
  722. exc.arg1 = x;
  723. exc.retval = HUGE_VAL;
  724. exc.type = DOMAIN;
  725. if (!matherr (&exc))
  726. {
  727. fprintf (stderr, "invnorm: DOMAIN error\n");
  728. errno = EDOM;
  729. }
  730. return exc.retval;
  731. #else
  732. errno = EDOM;
  733. return HUGE_VAL;
  734. #endif
  735. }
  736. return -M_SQRT2 * inverf (1.0 - 2 * x);
  737. }