ibeta_inverse.hpp 33 KB

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