ibeta_inverse.hpp 37 KB

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