ibeta_inverse.hpp 38 KB

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