ibeta_inverse.hpp 36 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093
  1. // Copyright John Maddock 2006.
  2. // Copyright Paul A. Bristow 2007
  3. // Use, modification and distribution are subject to the
  4. // Boost Software License, Version 1.0. (See accompanying file
  5. // LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
  6. #ifndef BOOST_MATH_SPECIAL_FUNCTIONS_IBETA_INVERSE_HPP
  7. #define BOOST_MATH_SPECIAL_FUNCTIONS_IBETA_INVERSE_HPP
  8. #ifdef _MSC_VER
  9. #pragma once
  10. #endif
  11. #include <boost/math/special_functions/beta.hpp>
  12. #include <boost/math/special_functions/erf.hpp>
  13. #include <boost/math/tools/roots.hpp>
  14. #include <boost/math/special_functions/detail/t_distribution_inv.hpp>
  15. #include <boost/math/special_functions/fpclassify.hpp>
  16. #include <boost/math/tools/precision.hpp>
  17. namespace boost{ namespace math{ namespace detail{
  18. //
  19. // Helper object used by root finding
  20. // code to convert eta to x.
  21. //
  22. template <class T>
  23. struct temme_root_finder
  24. {
  25. temme_root_finder(const T t_, const T a_) : t(t_), a(a_) {
  26. BOOST_MATH_ASSERT(
  27. math::tools::epsilon<T>() <= a && !(boost::math::isinf)(a));
  28. }
  29. boost::math::tuple<T, T> operator()(T x)
  30. {
  31. BOOST_MATH_STD_USING // ADL of std names
  32. T y = 1 - x;
  33. T f = log(x) + a * log(y) + t;
  34. T f1 = (1 / x) - (a / (y));
  35. return boost::math::make_tuple(f, f1);
  36. }
  37. private:
  38. T t, a;
  39. };
  40. //
  41. // See:
  42. // "Asymptotic Inversion of the Incomplete Beta Function"
  43. // N.M. Temme
  44. // Journal of Computation and Applied Mathematics 41 (1992) 145-157.
  45. // Section 2.
  46. //
  47. template <class T, class Policy>
  48. T temme_method_1_ibeta_inverse(T a, T b, T z, const Policy& pol)
  49. {
  50. BOOST_MATH_STD_USING // ADL of std names
  51. const T r2 = sqrt(T(2));
  52. //
  53. // get the first approximation for eta from the inverse
  54. // error function (Eq: 2.9 and 2.10).
  55. //
  56. T eta0 = boost::math::erfc_inv(2 * z, pol);
  57. eta0 /= -sqrt(a / 2);
  58. T terms[4] = { eta0 };
  59. T workspace[7];
  60. //
  61. // calculate powers:
  62. //
  63. T B = b - a;
  64. T B_2 = B * B;
  65. T B_3 = B_2 * B;
  66. //
  67. // Calculate correction terms:
  68. //
  69. // See eq following 2.15:
  70. workspace[0] = -B * r2 / 2;
  71. workspace[1] = (1 - 2 * B) / 8;
  72. workspace[2] = -(B * r2 / 48);
  73. workspace[3] = T(-1) / 192;
  74. workspace[4] = -B * r2 / 3840;
  75. terms[1] = tools::evaluate_polynomial(workspace, eta0, 5);
  76. // Eq Following 2.17:
  77. workspace[0] = B * r2 * (3 * B - 2) / 12;
  78. workspace[1] = (20 * B_2 - 12 * B + 1) / 128;
  79. workspace[2] = B * r2 * (20 * B - 1) / 960;
  80. workspace[3] = (16 * B_2 + 30 * B - 15) / 4608;
  81. workspace[4] = B * r2 * (21 * B + 32) / 53760;
  82. workspace[5] = (-32 * B_2 + 63) / 368640;
  83. workspace[6] = -B * r2 * (120 * B + 17) / 25804480;
  84. terms[2] = tools::evaluate_polynomial(workspace, eta0, 7);
  85. // Eq Following 2.17:
  86. workspace[0] = B * r2 * (-75 * B_2 + 80 * B - 16) / 480;
  87. workspace[1] = (-1080 * B_3 + 868 * B_2 - 90 * B - 45) / 9216;
  88. workspace[2] = B * r2 * (-1190 * B_2 + 84 * B + 373) / 53760;
  89. workspace[3] = (-2240 * B_3 - 2508 * B_2 + 2100 * B - 165) / 368640;
  90. terms[3] = tools::evaluate_polynomial(workspace, eta0, 4);
  91. //
  92. // Bring them together to get a final estimate for eta:
  93. //
  94. T eta = tools::evaluate_polynomial(terms, T(1/a), 4);
  95. //
  96. // now we need to convert eta to x, by solving the appropriate
  97. // quadratic equation:
  98. //
  99. T eta_2 = eta * eta;
  100. T c = -exp(-eta_2 / 2);
  101. T x;
  102. if(eta_2 == 0)
  103. x = static_cast<T>(0.5f);
  104. else
  105. x = (1 + eta * sqrt((1 + c) / eta_2)) / 2;
  106. //
  107. // These are post-conditions of the method, but the addition above
  108. // may result in us being out by 1ulp either side of the boundary,
  109. // so just check that we're in bounds and adjust as needed.
  110. // See https://github.com/boostorg/math/issues/961
  111. //
  112. if (x < 0)
  113. x = 0;
  114. else if (x > 1)
  115. x = 1;
  116. BOOST_MATH_ASSERT(eta * (x - 0.5) >= 0);
  117. #ifdef BOOST_INSTRUMENT
  118. std::cout << "Estimating x with Temme method 1: " << x << std::endl;
  119. #endif
  120. return x;
  121. }
  122. //
  123. // See:
  124. // "Asymptotic Inversion of the Incomplete Beta Function"
  125. // N.M. Temme
  126. // Journal of Computation and Applied Mathematics 41 (1992) 145-157.
  127. // Section 3.
  128. //
  129. template <class T, class Policy>
  130. T temme_method_2_ibeta_inverse(T /*a*/, T /*b*/, T z, T r, T theta, const Policy& pol)
  131. {
  132. BOOST_MATH_STD_USING // ADL of std names
  133. //
  134. // Get first estimate for eta, see Eq 3.9 and 3.10,
  135. // but note there is a typo in Eq 3.10:
  136. //
  137. T eta0 = boost::math::erfc_inv(2 * z, pol);
  138. eta0 /= -sqrt(r / 2);
  139. T s = sin(theta);
  140. T c = cos(theta);
  141. //
  142. // Now we need to perturb eta0 to get eta, which we do by
  143. // evaluating the polynomial in 1/r at the bottom of page 151,
  144. // to do this we first need the error terms e1, e2 e3
  145. // which we'll fill into the array "terms". Since these
  146. // terms are themselves polynomials, we'll need another
  147. // array "workspace" to calculate those...
  148. //
  149. T terms[4] = { eta0 };
  150. T workspace[6];
  151. //
  152. // some powers of sin(theta)cos(theta) that we'll need later:
  153. //
  154. T sc = s * c;
  155. T sc_2 = sc * sc;
  156. T sc_3 = sc_2 * sc;
  157. T sc_4 = sc_2 * sc_2;
  158. T sc_5 = sc_2 * sc_3;
  159. T sc_6 = sc_3 * sc_3;
  160. T sc_7 = sc_4 * sc_3;
  161. //
  162. // Calculate e1 and put it in terms[1], see the middle of page 151:
  163. //
  164. workspace[0] = (2 * s * s - 1) / (3 * s * c);
  165. static const BOOST_MATH_INT_TABLE_TYPE(T, int) co1[] = { -1, -5, 5 };
  166. workspace[1] = -tools::evaluate_even_polynomial(co1, s, 3) / (36 * sc_2);
  167. static const BOOST_MATH_INT_TABLE_TYPE(T, int) co2[] = { 1, 21, -69, 46 };
  168. workspace[2] = tools::evaluate_even_polynomial(co2, s, 4) / (1620 * sc_3);
  169. static const BOOST_MATH_INT_TABLE_TYPE(T, int) co3[] = { 7, -2, 33, -62, 31 };
  170. workspace[3] = -tools::evaluate_even_polynomial(co3, s, 5) / (6480 * sc_4);
  171. static const BOOST_MATH_INT_TABLE_TYPE(T, int) co4[] = { 25, -52, -17, 88, -115, 46 };
  172. workspace[4] = tools::evaluate_even_polynomial(co4, s, 6) / (90720 * sc_5);
  173. terms[1] = tools::evaluate_polynomial(workspace, eta0, 5);
  174. //
  175. // Now evaluate e2 and put it in terms[2]:
  176. //
  177. static const BOOST_MATH_INT_TABLE_TYPE(T, int) co5[] = { 7, 12, -78, 52 };
  178. workspace[0] = -tools::evaluate_even_polynomial(co5, s, 4) / (405 * sc_3);
  179. static const BOOST_MATH_INT_TABLE_TYPE(T, int) co6[] = { -7, 2, 183, -370, 185 };
  180. workspace[1] = tools::evaluate_even_polynomial(co6, s, 5) / (2592 * sc_4);
  181. static const BOOST_MATH_INT_TABLE_TYPE(T, int) co7[] = { -533, 776, -1835, 10240, -13525, 5410 };
  182. workspace[2] = -tools::evaluate_even_polynomial(co7, s, 6) / (204120 * sc_5);
  183. static const BOOST_MATH_INT_TABLE_TYPE(T, int) co8[] = { -1579, 3747, -3372, -15821, 45588, -45213, 15071 };
  184. workspace[3] = -tools::evaluate_even_polynomial(co8, s, 7) / (2099520 * sc_6);
  185. terms[2] = tools::evaluate_polynomial(workspace, eta0, 4);
  186. //
  187. // And e3, and put it in terms[3]:
  188. //
  189. static const BOOST_MATH_INT_TABLE_TYPE(T, int) co9[] = {449, -1259, -769, 6686, -9260, 3704 };
  190. workspace[0] = tools::evaluate_even_polynomial(co9, s, 6) / (102060 * sc_5);
  191. static const BOOST_MATH_INT_TABLE_TYPE(T, int) co10[] = { 63149, -151557, 140052, -727469, 2239932, -2251437, 750479 };
  192. workspace[1] = -tools::evaluate_even_polynomial(co10, s, 7) / (20995200 * sc_6);
  193. static const BOOST_MATH_INT_TABLE_TYPE(T, int) co11[] = { 29233, -78755, 105222, 146879, -1602610, 3195183, -2554139, 729754 };
  194. workspace[2] = tools::evaluate_even_polynomial(co11, s, 8) / (36741600 * sc_7);
  195. terms[3] = tools::evaluate_polynomial(workspace, eta0, 3);
  196. //
  197. // Bring the correction terms together to evaluate eta,
  198. // this is the last equation on page 151:
  199. //
  200. T eta = tools::evaluate_polynomial(terms, T(1/r), 4);
  201. //
  202. // Now that we have eta we need to back solve for x,
  203. // we seek the value of x that gives eta in Eq 3.2.
  204. // The two methods used are described in section 5.
  205. //
  206. // Begin by defining a few variables we'll need later:
  207. //
  208. T x;
  209. T s_2 = s * s;
  210. T c_2 = c * c;
  211. T alpha = c / s;
  212. alpha *= alpha;
  213. T lu = (-(eta * eta) / (2 * s_2) + log(s_2) + c_2 * log(c_2) / s_2);
  214. //
  215. // Temme doesn't specify what value to switch on here,
  216. // but this seems to work pretty well:
  217. //
  218. if(fabs(eta) < 0.7)
  219. {
  220. //
  221. // Small eta use the expansion Temme gives in the second equation
  222. // of section 5, it's a polynomial in eta:
  223. //
  224. workspace[0] = s * s;
  225. workspace[1] = s * c;
  226. workspace[2] = (1 - 2 * workspace[0]) / 3;
  227. static const BOOST_MATH_INT_TABLE_TYPE(T, int) co12[] = { 1, -13, 13 };
  228. workspace[3] = tools::evaluate_polynomial(co12, workspace[0], 3) / (36 * s * c);
  229. static const BOOST_MATH_INT_TABLE_TYPE(T, int) co13[] = { 1, 21, -69, 46 };
  230. workspace[4] = tools::evaluate_polynomial(co13, workspace[0], 4) / (270 * workspace[0] * c * c);
  231. x = tools::evaluate_polynomial(workspace, eta, 5);
  232. #ifdef BOOST_INSTRUMENT
  233. std::cout << "Estimating x with Temme method 2 (small eta): " << x << std::endl;
  234. #endif
  235. }
  236. else
  237. {
  238. //
  239. // If eta is large we need to solve Eq 3.2 more directly,
  240. // begin by getting an initial approximation for x from
  241. // the last equation on page 155, this is a polynomial in u:
  242. //
  243. T u = exp(lu);
  244. workspace[0] = u;
  245. workspace[1] = alpha;
  246. workspace[2] = 0;
  247. workspace[3] = 3 * alpha * (3 * alpha + 1) / 6;
  248. workspace[4] = 4 * alpha * (4 * alpha + 1) * (4 * alpha + 2) / 24;
  249. workspace[5] = 5 * alpha * (5 * alpha + 1) * (5 * alpha + 2) * (5 * alpha + 3) / 120;
  250. x = tools::evaluate_polynomial(workspace, u, 6);
  251. //
  252. // At this point we may or may not have the right answer, Eq-3.2 has
  253. // two solutions for x for any given eta, however the mapping in 3.2
  254. // is 1:1 with the sign of eta and x-sin^2(theta) being the same.
  255. // So we can check if we have the right root of 3.2, and if not
  256. // switch x for 1-x. This transformation is motivated by the fact
  257. // that the distribution is *almost* symmetric so 1-x will be in the right
  258. // ball park for the solution:
  259. //
  260. if((x - s_2) * eta < 0)
  261. x = 1 - x;
  262. #ifdef BOOST_INSTRUMENT
  263. std::cout << "Estimating x with Temme method 2 (large eta): " << x << std::endl;
  264. #endif
  265. }
  266. //
  267. // The final step is a few Newton-Raphson iterations to
  268. // clean up our approximation for x, this is pretty cheap
  269. // in general, and very cheap compared to an incomplete beta
  270. // evaluation. The limits set on x come from the observation
  271. // that the sign of eta and x-sin^2(theta) are the same.
  272. //
  273. T lower, upper;
  274. if(eta < 0)
  275. {
  276. lower = 0;
  277. upper = s_2;
  278. }
  279. else
  280. {
  281. lower = s_2;
  282. upper = 1;
  283. }
  284. //
  285. // If our initial approximation is out of bounds then bisect:
  286. //
  287. if((x < lower) || (x > upper))
  288. x = (lower+upper) / 2;
  289. //
  290. // And iterate:
  291. //
  292. x = tools::newton_raphson_iterate(
  293. temme_root_finder<T>(-lu, alpha), x, lower, upper, policies::digits<T, Policy>() / 2);
  294. return x;
  295. }
  296. //
  297. // See:
  298. // "Asymptotic Inversion of the Incomplete Beta Function"
  299. // N.M. Temme
  300. // Journal of Computation and Applied Mathematics 41 (1992) 145-157.
  301. // Section 4.
  302. //
  303. template <class T, class Policy>
  304. T temme_method_3_ibeta_inverse(T a, T b, T p, T q, const Policy& pol)
  305. {
  306. BOOST_MATH_STD_USING // ADL of std names
  307. //
  308. // Begin by getting an initial approximation for the quantity
  309. // eta from the dominant part of the incomplete beta:
  310. //
  311. T eta0;
  312. if(p < q)
  313. eta0 = boost::math::gamma_q_inv(b, p, pol);
  314. else
  315. eta0 = boost::math::gamma_p_inv(b, q, pol);
  316. eta0 /= a;
  317. //
  318. // Define the variables and powers we'll need later on:
  319. //
  320. T mu = b / a;
  321. T w = sqrt(1 + mu);
  322. T w_2 = w * w;
  323. T w_3 = w_2 * w;
  324. T w_4 = w_2 * w_2;
  325. T w_5 = w_3 * w_2;
  326. T w_6 = w_3 * w_3;
  327. T w_7 = w_4 * w_3;
  328. T w_8 = w_4 * w_4;
  329. T w_9 = w_5 * w_4;
  330. T w_10 = w_5 * w_5;
  331. T d = eta0 - mu;
  332. T d_2 = d * d;
  333. T d_3 = d_2 * d;
  334. T d_4 = d_2 * d_2;
  335. T w1 = w + 1;
  336. T w1_2 = w1 * w1;
  337. T w1_3 = w1 * w1_2;
  338. T w1_4 = w1_2 * w1_2;
  339. //
  340. // Now we need to compute the perturbation error terms that
  341. // convert eta0 to eta, these are all polynomials of polynomials.
  342. // Probably these should be re-written to use tabulated data
  343. // (see examples above), but it's less of a win in this case as we
  344. // need to calculate the individual powers for the denominator terms
  345. // anyway, so we might as well use them for the numerator-polynomials
  346. // as well....
  347. //
  348. // Refer to p154-p155 for the details of these expansions:
  349. //
  350. T e1 = (w + 2) * (w - 1) / (3 * w);
  351. e1 += (w_3 + 9 * w_2 + 21 * w + 5) * d / (36 * w_2 * w1);
  352. e1 -= (w_4 - 13 * w_3 + 69 * w_2 + 167 * w + 46) * d_2 / (1620 * w1_2 * w_3);
  353. e1 -= (7 * w_5 + 21 * w_4 + 70 * w_3 + 26 * w_2 - 93 * w - 31) * d_3 / (6480 * w1_3 * w_4);
  354. e1 -= (75 * w_6 + 202 * w_5 + 188 * w_4 - 888 * w_3 - 1345 * w_2 + 118 * w + 138) * d_4 / (272160 * w1_4 * w_5);
  355. T e2 = (28 * w_4 + 131 * w_3 + 402 * w_2 + 581 * w + 208) * (w - 1) / (1620 * w1 * w_3);
  356. e2 -= (35 * w_6 - 154 * w_5 - 623 * w_4 - 1636 * w_3 - 3983 * w_2 - 3514 * w - 925) * d / (12960 * w1_2 * w_4);
  357. e2 -= (2132 * w_7 + 7915 * w_6 + 16821 * w_5 + 35066 * w_4 + 87490 * w_3 + 141183 * w_2 + 95993 * w + 21640) * d_2 / (816480 * w_5 * w1_3);
  358. e2 -= (11053 * w_8 + 53308 * w_7 + 117010 * w_6 + 163924 * w_5 + 116188 * w_4 - 258428 * w_3 - 677042 * w_2 - 481940 * w - 105497) * d_3 / (T(14696640) * w1_4 * w_6);
  359. T e3 = -((3592 * w_7 + 8375 * w_6 - 1323 * w_5 - 29198 * w_4 - 89578 * w_3 - 154413 * w_2 - 116063 * w - 29632) * (w - 1)) / (816480 * w_5 * w1_2);
  360. e3 -= (442043 * w_9 + T(2054169) * w_8 + T(3803094) * w_7 + T(3470754) * w_6 + T(2141568) * w_5 - T(2393568) * w_4 - T(19904934) * w_3 - T(34714674) * w_2 - T(23128299) * w - T(5253353)) * d / (T(146966400) * w_6 * w1_3);
  361. e3 -= (116932 * w_10 + 819281 * w_9 + T(2378172) * w_8 + T(4341330) * w_7 + T(6806004) * w_6 + T(10622748) * w_5 + T(18739500) * w_4 + T(30651894) * w_3 + T(30869976) * w_2 + T(15431867) * w + T(2919016)) * d_2 / (T(146966400) * w1_4 * w_7);
  362. //
  363. // Combine eta0 and the error terms to compute eta (Second equation p155):
  364. //
  365. T eta = eta0 + e1 / a + e2 / (a * a) + e3 / (a * a * a);
  366. //
  367. // Now we need to solve Eq 4.2 to obtain x. For any given value of
  368. // eta there are two solutions to this equation, and since the distribution
  369. // may be very skewed, these are not related by x ~ 1-x we used when
  370. // implementing section 3 above. However we know that:
  371. //
  372. // cross < x <= 1 ; iff eta < mu
  373. // x == cross ; iff eta == mu
  374. // 0 <= x < cross ; iff eta > mu
  375. //
  376. // Where cross == 1 / (1 + mu)
  377. // Many thanks to Prof Temme for clarifying this point.
  378. //
  379. // Therefore we'll just jump straight into Newton iterations
  380. // to solve Eq 4.2 using these bounds, and simple bisection
  381. // as the first guess, in practice this converges pretty quickly
  382. // and we only need a few digits correct anyway:
  383. //
  384. if(eta <= 0)
  385. eta = tools::min_value<T>();
  386. T u = eta - mu * log(eta) + (1 + mu) * log(1 + mu) - mu;
  387. T cross = 1 / (1 + mu);
  388. T lower = eta < mu ? cross : 0;
  389. T upper = eta < mu ? 1 : cross;
  390. T x = (lower + upper) / 2;
  391. // Early exit for cases with numerical precision issues.
  392. if (cross == 0 || cross == 1) { return cross; }
  393. x = tools::newton_raphson_iterate(
  394. temme_root_finder<T>(u, mu), x, lower, upper, policies::digits<T, Policy>() / 2);
  395. #ifdef BOOST_INSTRUMENT
  396. std::cout << "Estimating x with Temme method 3: " << x << std::endl;
  397. #endif
  398. return x;
  399. }
  400. template <class T, class Policy>
  401. struct ibeta_roots
  402. {
  403. ibeta_roots(T _a, T _b, T t, bool inv = false)
  404. : a(_a), b(_b), target(t), invert(inv) {}
  405. boost::math::tuple<T, T, T> operator()(T x)
  406. {
  407. BOOST_MATH_STD_USING // ADL of std names
  408. BOOST_FPU_EXCEPTION_GUARD
  409. T f1;
  410. T y = 1 - x;
  411. T f = ibeta_imp(a, b, x, Policy(), invert, true, &f1) - target;
  412. if(invert)
  413. f1 = -f1;
  414. if(y == 0)
  415. y = tools::min_value<T>() * 64;
  416. if(x == 0)
  417. x = tools::min_value<T>() * 64;
  418. T f2 = f1 * (-y * a + (b - 2) * x + 1);
  419. if(fabs(f2) < y * x * tools::max_value<T>())
  420. f2 /= (y * x);
  421. if(invert)
  422. f2 = -f2;
  423. // make sure we don't have a zero derivative:
  424. if(f1 == 0)
  425. f1 = (invert ? -1 : 1) * tools::min_value<T>() * 64;
  426. return boost::math::make_tuple(f, f1, f2);
  427. }
  428. private:
  429. T a, b, target;
  430. bool invert;
  431. };
  432. template <class T, class Policy>
  433. T ibeta_inv_imp(T a, T b, T p, T q, const Policy& pol, T* py)
  434. {
  435. BOOST_MATH_STD_USING // For ADL of math functions.
  436. //
  437. // The flag invert is set to true if we swap a for b and p for q,
  438. // in which case the result has to be subtracted from 1:
  439. //
  440. bool invert = false;
  441. //
  442. // Handle trivial cases first:
  443. //
  444. if(q == 0)
  445. {
  446. if(py) *py = 0;
  447. return 1;
  448. }
  449. else if(p == 0)
  450. {
  451. if(py) *py = 1;
  452. return 0;
  453. }
  454. else if(a == 1)
  455. {
  456. if(b == 1)
  457. {
  458. if(py) *py = 1 - p;
  459. return p;
  460. }
  461. // Change things around so we can handle as b == 1 special case below:
  462. std::swap(a, b);
  463. std::swap(p, q);
  464. invert = true;
  465. }
  466. //
  467. // Depending upon which approximation method we use, we may end up
  468. // calculating either x or y initially (where y = 1-x):
  469. //
  470. T x = 0; // Set to a safe zero to avoid a
  471. // MSVC 2005 warning C4701: potentially uninitialized local variable 'x' used
  472. // But code inspection appears to ensure that x IS assigned whatever the code path.
  473. T y;
  474. // For some of the methods we can put tighter bounds
  475. // on the result than simply [0,1]:
  476. //
  477. T lower = 0;
  478. T upper = 1;
  479. //
  480. // Student's T with b = 0.5 gets handled as a special case, swap
  481. // around if the arguments are in the "wrong" order:
  482. //
  483. if(a == 0.5f)
  484. {
  485. if(b == 0.5f)
  486. {
  487. x = sin(p * constants::half_pi<T>());
  488. x *= x;
  489. if(py)
  490. {
  491. *py = sin(q * constants::half_pi<T>());
  492. *py *= *py;
  493. }
  494. return x;
  495. }
  496. else if(b > 0.5f)
  497. {
  498. std::swap(a, b);
  499. std::swap(p, q);
  500. invert = !invert;
  501. }
  502. }
  503. //
  504. // Select calculation method for the initial estimate:
  505. //
  506. if((b == 0.5f) && (a >= 0.5f) && (p != 1))
  507. {
  508. //
  509. // We have a Student's T distribution:
  510. x = find_ibeta_inv_from_t_dist(a, p, q, &y, pol);
  511. }
  512. else if(b == 1)
  513. {
  514. if(p < q)
  515. {
  516. if(a > 1)
  517. {
  518. x = pow(p, 1 / a);
  519. y = -boost::math::expm1(log(p) / a, pol);
  520. }
  521. else
  522. {
  523. x = pow(p, 1 / a);
  524. y = 1 - x;
  525. }
  526. }
  527. else
  528. {
  529. x = exp(boost::math::log1p(-q, pol) / a);
  530. y = -boost::math::expm1(boost::math::log1p(-q, pol) / a, pol);
  531. }
  532. if(invert)
  533. std::swap(x, y);
  534. if(py)
  535. *py = y;
  536. return x;
  537. }
  538. else if(a + b > 5)
  539. {
  540. //
  541. // When a+b is large then we can use one of Prof Temme's
  542. // asymptotic expansions, begin by swapping things around
  543. // so that p < 0.5, we do this to avoid cancellations errors
  544. // when p is large.
  545. //
  546. if(p > 0.5)
  547. {
  548. std::swap(a, b);
  549. std::swap(p, q);
  550. invert = !invert;
  551. }
  552. T minv = (std::min)(a, b);
  553. T maxv = (std::max)(a, b);
  554. if((sqrt(minv) > (maxv - minv)) && (minv > 5))
  555. {
  556. //
  557. // When a and b differ by a small amount
  558. // the curve is quite symmetrical and we can use an error
  559. // function to approximate the inverse. This is the cheapest
  560. // of the three Temme expansions, and the calculated value
  561. // for x will never be much larger than p, so we don't have
  562. // to worry about cancellation as long as p is small.
  563. //
  564. x = temme_method_1_ibeta_inverse(a, b, p, pol);
  565. y = 1 - x;
  566. }
  567. else
  568. {
  569. T r = a + b;
  570. T theta = asin(sqrt(a / r));
  571. T lambda = minv / r;
  572. if((lambda >= 0.2) && (lambda <= 0.8) && (r >= 10))
  573. {
  574. //
  575. // The second error function case is the next cheapest
  576. // to use, it brakes down when the result is likely to be
  577. // very small, if a+b is also small, but we can use a
  578. // cheaper expansion there in any case. As before x won't
  579. // be much larger than p, so as long as p is small we should
  580. // be free of cancellation error.
  581. //
  582. T ppa = pow(p, 1/a);
  583. if((ppa < 0.0025) && (a + b < 200))
  584. {
  585. x = ppa * pow(a * boost::math::beta(a, b, pol), 1/a);
  586. }
  587. else
  588. x = temme_method_2_ibeta_inverse(a, b, p, r, theta, pol);
  589. y = 1 - x;
  590. }
  591. else
  592. {
  593. //
  594. // If we get here then a and b are very different in magnitude
  595. // and we need to use the third of Temme's methods which
  596. // involves inverting the incomplete gamma. This is much more
  597. // expensive than the other methods. We also can only use this
  598. // method when a > b, which can lead to cancellation errors
  599. // if we really want y (as we will when x is close to 1), so
  600. // a different expansion is used in that case.
  601. //
  602. if(a < b)
  603. {
  604. std::swap(a, b);
  605. std::swap(p, q);
  606. invert = !invert;
  607. }
  608. //
  609. // Try and compute the easy way first:
  610. //
  611. T bet = 0;
  612. if (b < 2)
  613. {
  614. #ifndef BOOST_MATH_NO_EXCEPTIONS
  615. try
  616. #endif
  617. {
  618. bet = boost::math::beta(a, b, pol);
  619. typedef typename Policy::overflow_error_type overflow_type;
  620. BOOST_MATH_IF_CONSTEXPR(overflow_type::value != boost::math::policies::throw_on_error)
  621. if(bet > tools::max_value<T>())
  622. bet = tools::max_value<T>();
  623. }
  624. #ifndef BOOST_MATH_NO_EXCEPTIONS
  625. catch (const std::overflow_error&)
  626. {
  627. bet = tools::max_value<T>();
  628. }
  629. #endif
  630. }
  631. if(bet != 0)
  632. {
  633. y = pow(b * q * bet, 1/b);
  634. x = 1 - y;
  635. }
  636. else
  637. y = 1;
  638. if(y > 1e-5)
  639. {
  640. x = temme_method_3_ibeta_inverse(a, b, p, q, pol);
  641. y = 1 - x;
  642. }
  643. }
  644. }
  645. }
  646. else if((a < 1) && (b < 1))
  647. {
  648. //
  649. // Both a and b less than 1,
  650. // there is a point of inflection at xs:
  651. //
  652. T xs = (1 - a) / (2 - a - b);
  653. //
  654. // Now we need to ensure that we start our iteration from the
  655. // right side of the inflection point:
  656. //
  657. T fs = boost::math::ibeta(a, b, xs, pol) - p;
  658. if(fabs(fs) / p < tools::epsilon<T>() * 3)
  659. {
  660. // The result is at the point of inflection, best just return it:
  661. *py = invert ? xs : 1 - xs;
  662. return invert ? 1-xs : xs;
  663. }
  664. if(fs < 0)
  665. {
  666. std::swap(a, b);
  667. std::swap(p, q);
  668. invert = !invert;
  669. xs = 1 - xs;
  670. }
  671. if ((a < tools::min_value<T>()) && (b > tools::min_value<T>()))
  672. {
  673. if (py)
  674. {
  675. *py = invert ? 0 : 1;
  676. }
  677. return invert ? 1 : 0; // nothing interesting going on here.
  678. }
  679. //
  680. // The call to beta may overflow, plus the alternative using lgamma may do the same
  681. // if T is a type where 1/T is infinite for small values (denorms for example).
  682. //
  683. T bet = 0;
  684. T xg;
  685. bool overflow = false;
  686. #ifndef BOOST_MATH_NO_EXCEPTIONS
  687. try {
  688. #endif
  689. bet = boost::math::beta(a, b, pol);
  690. #ifndef BOOST_MATH_NO_EXCEPTIONS
  691. }
  692. catch (const std::runtime_error&)
  693. {
  694. overflow = true;
  695. }
  696. #endif
  697. if (overflow || !(boost::math::isfinite)(bet))
  698. {
  699. xg = exp((boost::math::lgamma(a + 1, pol) + boost::math::lgamma(b, pol) - boost::math::lgamma(a + b, pol) + log(p)) / a);
  700. if (xg > 2 / tools::epsilon<T>())
  701. xg = 2 / tools::epsilon<T>();
  702. }
  703. else
  704. xg = pow(a * p * bet, 1/a);
  705. x = xg / (1 + xg);
  706. y = 1 / (1 + xg);
  707. //
  708. // And finally we know that our result is below the inflection
  709. // point, so set an upper limit on our search:
  710. //
  711. if(x > xs)
  712. x = xs;
  713. upper = xs;
  714. }
  715. else if((a > 1) && (b > 1))
  716. {
  717. //
  718. // Small a and b, both greater than 1,
  719. // there is a point of inflection at xs,
  720. // and it's complement is xs2, we must always
  721. // start our iteration from the right side of the
  722. // point of inflection.
  723. //
  724. T xs = (a - 1) / (a + b - 2);
  725. T xs2 = (b - 1) / (a + b - 2);
  726. T ps = boost::math::ibeta(a, b, xs, pol) - p;
  727. if(ps < 0)
  728. {
  729. std::swap(a, b);
  730. std::swap(p, q);
  731. std::swap(xs, xs2);
  732. invert = !invert;
  733. }
  734. //
  735. // Estimate x and y, using expm1 to get a good estimate
  736. // for y when it's very small:
  737. //
  738. T lx = log(p * a * boost::math::beta(a, b, pol)) / a;
  739. x = exp(lx);
  740. y = x < 0.9 ? T(1 - x) : (T)(-boost::math::expm1(lx, pol));
  741. if((b < a) && (x < 0.2))
  742. {
  743. //
  744. // Under a limited range of circumstances we can improve
  745. // our estimate for x, frankly it's clear if this has much effect!
  746. //
  747. T ap1 = a - 1;
  748. T bm1 = b - 1;
  749. T a_2 = a * a;
  750. T a_3 = a * a_2;
  751. T b_2 = b * b;
  752. T terms[5] = { 0, 1 };
  753. terms[2] = bm1 / ap1;
  754. ap1 *= ap1;
  755. terms[3] = bm1 * (3 * a * b + 5 * b + a_2 - a - 4) / (2 * (a + 2) * ap1);
  756. ap1 *= (a + 1);
  757. terms[4] = bm1 * (33 * a * b_2 + 31 * b_2 + 8 * a_2 * b_2 - 30 * a * b - 47 * b + 11 * a_2 * b + 6 * a_3 * b + 18 + 4 * a - a_3 + a_2 * a_2 - 10 * a_2)
  758. / (3 * (a + 3) * (a + 2) * ap1);
  759. x = tools::evaluate_polynomial(terms, x, 5);
  760. }
  761. //
  762. // And finally we know that our result is below the inflection
  763. // point, so set an upper limit on our search:
  764. //
  765. if(x > xs)
  766. x = xs;
  767. upper = xs;
  768. }
  769. else /*if((a <= 1) != (b <= 1))*/
  770. {
  771. //
  772. // If all else fails we get here, only one of a and b
  773. // is above 1, and a+b is small. Start by swapping
  774. // things around so that we have a concave curve with b > a
  775. // and no points of inflection in [0,1]. As long as we expect
  776. // x to be small then we can use the simple (and cheap) power
  777. // term to estimate x, but when we expect x to be large then
  778. // this greatly underestimates x and leaves us trying to
  779. // iterate "round the corner" which may take almost forever...
  780. //
  781. // We could use Temme's inverse gamma function case in that case,
  782. // this works really rather well (albeit expensively) even though
  783. // strictly speaking we're outside it's defined range.
  784. //
  785. // However it's expensive to compute, and an alternative approach
  786. // which models the curve as a distorted quarter circle is much
  787. // cheaper to compute, and still keeps the number of iterations
  788. // required down to a reasonable level. With thanks to Prof Temme
  789. // for this suggestion.
  790. //
  791. if(b < a)
  792. {
  793. std::swap(a, b);
  794. std::swap(p, q);
  795. invert = !invert;
  796. }
  797. if (a < tools::min_value<T>())
  798. {
  799. // Avoid spurious overflows for denorms:
  800. if (p < 1)
  801. {
  802. x = 1;
  803. y = 0;
  804. }
  805. else
  806. {
  807. x = 0;
  808. y = 1;
  809. }
  810. }
  811. else if(pow(p, 1/a) < 0.5)
  812. {
  813. #ifndef BOOST_MATH_NO_EXCEPTIONS
  814. try
  815. {
  816. #endif
  817. x = pow(p * a * boost::math::beta(a, b, pol), 1 / a);
  818. if ((x > 1) || !(boost::math::isfinite)(x))
  819. x = 1;
  820. #ifndef BOOST_MATH_NO_EXCEPTIONS
  821. }
  822. catch (const std::overflow_error&)
  823. {
  824. x = 1;
  825. }
  826. #endif
  827. if(x == 0)
  828. x = boost::math::tools::min_value<T>();
  829. y = 1 - x;
  830. }
  831. else /*if(pow(q, 1/b) < 0.1)*/
  832. {
  833. // model a distorted quarter circle:
  834. #ifndef BOOST_MATH_NO_EXCEPTIONS
  835. try
  836. {
  837. #endif
  838. y = pow(1 - pow(p, b * boost::math::beta(a, b, pol)), 1/b);
  839. if ((y > 1) || !(boost::math::isfinite)(y))
  840. y = 1;
  841. #ifndef BOOST_MATH_NO_EXCEPTIONS
  842. }
  843. catch (const std::overflow_error&)
  844. {
  845. y = 1;
  846. }
  847. #endif
  848. if(y == 0)
  849. y = boost::math::tools::min_value<T>();
  850. x = 1 - y;
  851. }
  852. }
  853. //
  854. // Now we have a guess for x (and for y) we can set things up for
  855. // iteration. If x > 0.5 it pays to swap things round:
  856. //
  857. if(x > 0.5)
  858. {
  859. std::swap(a, b);
  860. std::swap(p, q);
  861. std::swap(x, y);
  862. invert = !invert;
  863. T l = 1 - upper;
  864. T u = 1 - lower;
  865. lower = l;
  866. upper = u;
  867. }
  868. //
  869. // lower bound for our search:
  870. //
  871. // We're not interested in denormalised answers as these tend to
  872. // these tend to take up lots of iterations, given that we can't get
  873. // accurate derivatives in this area (they tend to be infinite).
  874. //
  875. if(lower == 0)
  876. {
  877. if(invert && (py == 0))
  878. {
  879. //
  880. // We're not interested in answers smaller than machine epsilon:
  881. //
  882. lower = boost::math::tools::epsilon<T>();
  883. if(x < lower)
  884. x = lower;
  885. }
  886. else
  887. lower = boost::math::tools::min_value<T>();
  888. if(x < lower)
  889. x = lower;
  890. }
  891. std::uintmax_t max_iter = policies::get_max_root_iterations<Policy>();
  892. std::uintmax_t max_iter_used = 0;
  893. //
  894. // Figure out how many digits to iterate towards:
  895. //
  896. int digits = boost::math::policies::digits<T, Policy>() / 2;
  897. if((x < 1e-50) && ((a < 1) || (b < 1)))
  898. {
  899. //
  900. // If we're in a region where the first derivative is very
  901. // large, then we have to take care that the root-finder
  902. // doesn't terminate prematurely. We'll bump the precision
  903. // up to avoid this, but we have to take care not to set the
  904. // precision too high or the last few iterations will just
  905. // thrash around and convergence may be slow in this case.
  906. // Try 3/4 of machine epsilon:
  907. //
  908. digits *= 3;
  909. digits /= 2;
  910. }
  911. //
  912. // Now iterate, we can use either p or q as the target here
  913. // depending on which is smaller:
  914. //
  915. x = boost::math::tools::halley_iterate(
  916. boost::math::detail::ibeta_roots<T, Policy>(a, b, (p < q ? p : q), (p < q ? false : true)), x, lower, upper, digits, max_iter);
  917. policies::check_root_iterations<T>("boost::math::ibeta<%1%>(%1%, %1%, %1%)", max_iter + max_iter_used, pol);
  918. //
  919. // We don't really want these asserts here, but they are useful for sanity
  920. // checking that we have the limits right, uncomment if you suspect bugs *only*.
  921. //
  922. //BOOST_MATH_ASSERT(x != upper);
  923. //BOOST_MATH_ASSERT((x != lower) || (x == boost::math::tools::min_value<T>()) || (x == boost::math::tools::epsilon<T>()));
  924. //
  925. // Tidy up, if we "lower" was too high then zero is the best answer we have:
  926. //
  927. if(x == lower)
  928. x = 0;
  929. if(py)
  930. *py = invert ? x : 1 - x;
  931. return invert ? 1-x : x;
  932. }
  933. } // namespace detail
  934. template <class T1, class T2, class T3, class T4, class Policy>
  935. inline typename tools::promote_args<T1, T2, T3, T4>::type
  936. ibeta_inv(T1 a, T2 b, T3 p, T4* py, const Policy& pol)
  937. {
  938. static const char* function = "boost::math::ibeta_inv<%1%>(%1%,%1%,%1%)";
  939. BOOST_FPU_EXCEPTION_GUARD
  940. typedef typename tools::promote_args<T1, T2, T3, T4>::type result_type;
  941. typedef typename policies::evaluation<result_type, Policy>::type value_type;
  942. typedef typename policies::normalise<
  943. Policy,
  944. policies::promote_float<false>,
  945. policies::promote_double<false>,
  946. policies::discrete_quantile<>,
  947. policies::assert_undefined<> >::type forwarding_policy;
  948. if(a <= 0)
  949. return policies::raise_domain_error<result_type>(function, "The argument a to the incomplete beta function inverse must be greater than zero (got a=%1%).", a, pol);
  950. if(b <= 0)
  951. return policies::raise_domain_error<result_type>(function, "The argument b to the incomplete beta function inverse must be greater than zero (got b=%1%).", b, pol);
  952. if((p < 0) || (p > 1))
  953. return policies::raise_domain_error<result_type>(function, "Argument p outside the range [0,1] in the incomplete beta function inverse (got p=%1%).", p, pol);
  954. value_type rx, ry;
  955. rx = detail::ibeta_inv_imp(
  956. static_cast<value_type>(a),
  957. static_cast<value_type>(b),
  958. static_cast<value_type>(p),
  959. static_cast<value_type>(1 - p),
  960. forwarding_policy(), &ry);
  961. if(py) *py = policies::checked_narrowing_cast<T4, forwarding_policy>(ry, function);
  962. return policies::checked_narrowing_cast<result_type, forwarding_policy>(rx, function);
  963. }
  964. template <class T1, class T2, class T3, class T4>
  965. inline typename tools::promote_args<T1, T2, T3, T4>::type
  966. ibeta_inv(T1 a, T2 b, T3 p, T4* py)
  967. {
  968. return ibeta_inv(a, b, p, py, policies::policy<>());
  969. }
  970. template <class T1, class T2, class T3>
  971. inline typename tools::promote_args<T1, T2, T3>::type
  972. ibeta_inv(T1 a, T2 b, T3 p)
  973. {
  974. typedef typename tools::promote_args<T1, T2, T3>::type result_type;
  975. return ibeta_inv(a, b, p, static_cast<result_type*>(nullptr), policies::policy<>());
  976. }
  977. template <class T1, class T2, class T3, class Policy>
  978. inline typename tools::promote_args<T1, T2, T3>::type
  979. ibeta_inv(T1 a, T2 b, T3 p, const Policy& pol)
  980. {
  981. typedef typename tools::promote_args<T1, T2, T3>::type result_type;
  982. return ibeta_inv(a, b, p, static_cast<result_type*>(nullptr), pol);
  983. }
  984. template <class T1, class T2, class T3, class T4, class Policy>
  985. inline typename tools::promote_args<T1, T2, T3, T4>::type
  986. ibetac_inv(T1 a, T2 b, T3 q, T4* py, const Policy& pol)
  987. {
  988. static const char* function = "boost::math::ibetac_inv<%1%>(%1%,%1%,%1%)";
  989. BOOST_FPU_EXCEPTION_GUARD
  990. typedef typename tools::promote_args<T1, T2, T3, T4>::type result_type;
  991. typedef typename policies::evaluation<result_type, Policy>::type value_type;
  992. typedef typename policies::normalise<
  993. Policy,
  994. policies::promote_float<false>,
  995. policies::promote_double<false>,
  996. policies::discrete_quantile<>,
  997. policies::assert_undefined<> >::type forwarding_policy;
  998. if(a <= 0)
  999. return policies::raise_domain_error<result_type>(function, "The argument a to the incomplete beta function inverse must be greater than zero (got a=%1%).", a, pol);
  1000. if(b <= 0)
  1001. return policies::raise_domain_error<result_type>(function, "The argument b to the incomplete beta function inverse must be greater than zero (got b=%1%).", b, pol);
  1002. if((q < 0) || (q > 1))
  1003. return policies::raise_domain_error<result_type>(function, "Argument q outside the range [0,1] in the incomplete beta function inverse (got q=%1%).", q, pol);
  1004. value_type rx, ry;
  1005. rx = detail::ibeta_inv_imp(
  1006. static_cast<value_type>(a),
  1007. static_cast<value_type>(b),
  1008. static_cast<value_type>(1 - q),
  1009. static_cast<value_type>(q),
  1010. forwarding_policy(), &ry);
  1011. if(py) *py = policies::checked_narrowing_cast<T4, forwarding_policy>(ry, function);
  1012. return policies::checked_narrowing_cast<result_type, forwarding_policy>(rx, function);
  1013. }
  1014. template <class T1, class T2, class T3, class T4>
  1015. inline typename tools::promote_args<T1, T2, T3, T4>::type
  1016. ibetac_inv(T1 a, T2 b, T3 q, T4* py)
  1017. {
  1018. return ibetac_inv(a, b, q, py, policies::policy<>());
  1019. }
  1020. template <class RT1, class RT2, class RT3>
  1021. inline typename tools::promote_args<RT1, RT2, RT3>::type
  1022. ibetac_inv(RT1 a, RT2 b, RT3 q)
  1023. {
  1024. typedef typename tools::promote_args<RT1, RT2, RT3>::type result_type;
  1025. return ibetac_inv(a, b, q, static_cast<result_type*>(nullptr), policies::policy<>());
  1026. }
  1027. template <class RT1, class RT2, class RT3, class Policy>
  1028. inline typename tools::promote_args<RT1, RT2, RT3>::type
  1029. ibetac_inv(RT1 a, RT2 b, RT3 q, const Policy& pol)
  1030. {
  1031. typedef typename tools::promote_args<RT1, RT2, RT3>::type result_type;
  1032. return ibetac_inv(a, b, q, static_cast<result_type*>(nullptr), pol);
  1033. }
  1034. } // namespace math
  1035. } // namespace boost
  1036. #endif // BOOST_MATH_SPECIAL_FUNCTIONS_IGAMMA_INVERSE_HPP