beta.hpp 58 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536153715381539154015411542154315441545154615471548154915501551155215531554155515561557155815591560156115621563156415651566156715681569157015711572157315741575157615771578157915801581158215831584158515861587158815891590159115921593159415951596159715981599160016011602160316041605160616071608160916101611161216131614161516161617161816191620162116221623162416251626162716281629163016311632163316341635163616371638163916401641164216431644164516461647164816491650165116521653165416551656165716581659166016611662166316641665166616671668166916701671167216731674167516761677167816791680168116821683168416851686168716881689169016911692169316941695169616971698169917001701170217031704170517061707170817091710171117121713171417151716171717181719172017211722
  1. // (C) Copyright John Maddock 2006.
  2. // Use, modification and distribution are subject to the
  3. // Boost Software License, Version 1.0. (See accompanying file
  4. // LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
  5. #ifndef BOOST_MATH_SPECIAL_BETA_HPP
  6. #define BOOST_MATH_SPECIAL_BETA_HPP
  7. #ifdef _MSC_VER
  8. #pragma once
  9. #endif
  10. #include <boost/math/special_functions/math_fwd.hpp>
  11. #include <boost/math/tools/config.hpp>
  12. #include <boost/math/special_functions/gamma.hpp>
  13. #include <boost/math/special_functions/binomial.hpp>
  14. #include <boost/math/special_functions/factorials.hpp>
  15. #include <boost/math/special_functions/erf.hpp>
  16. #include <boost/math/special_functions/log1p.hpp>
  17. #include <boost/math/special_functions/expm1.hpp>
  18. #include <boost/math/special_functions/trunc.hpp>
  19. #include <boost/math/tools/roots.hpp>
  20. #include <boost/math/tools/assert.hpp>
  21. #include <cmath>
  22. namespace boost{ namespace math{
  23. namespace detail{
  24. //
  25. // Implementation of Beta(a,b) using the Lanczos approximation:
  26. //
  27. template <class T, class Lanczos, class Policy>
  28. T beta_imp(T a, T b, const Lanczos&, const Policy& pol)
  29. {
  30. BOOST_MATH_STD_USING // for ADL of std names
  31. if(a <= 0)
  32. return policies::raise_domain_error<T>("boost::math::beta<%1%>(%1%,%1%)", "The arguments to the beta function must be greater than zero (got a=%1%).", a, pol);
  33. if(b <= 0)
  34. return policies::raise_domain_error<T>("boost::math::beta<%1%>(%1%,%1%)", "The arguments to the beta function must be greater than zero (got b=%1%).", b, pol);
  35. T result; // LCOV_EXCL_LINE
  36. T prefix = 1;
  37. T c = a + b;
  38. // Special cases:
  39. if((c == a) && (b < tools::epsilon<T>()))
  40. return 1 / b;
  41. else if((c == b) && (a < tools::epsilon<T>()))
  42. return 1 / a;
  43. if(b == 1)
  44. return 1/a;
  45. else if(a == 1)
  46. return 1/b;
  47. else if(c < tools::epsilon<T>())
  48. {
  49. result = c / a;
  50. result /= b;
  51. return result;
  52. }
  53. /*
  54. //
  55. // This code appears to be no longer necessary: it was
  56. // used to offset errors introduced from the Lanczos
  57. // approximation, but the current Lanczos approximations
  58. // are sufficiently accurate for all z that we can ditch
  59. // this. It remains in the file for future reference...
  60. //
  61. // If a or b are less than 1, shift to greater than 1:
  62. if(a < 1)
  63. {
  64. prefix *= c / a;
  65. c += 1;
  66. a += 1;
  67. }
  68. if(b < 1)
  69. {
  70. prefix *= c / b;
  71. c += 1;
  72. b += 1;
  73. }
  74. */
  75. if(a < b)
  76. std::swap(a, b);
  77. // Lanczos calculation:
  78. T agh = static_cast<T>(a + Lanczos::g() - 0.5f);
  79. T bgh = static_cast<T>(b + Lanczos::g() - 0.5f);
  80. T cgh = static_cast<T>(c + Lanczos::g() - 0.5f);
  81. result = Lanczos::lanczos_sum_expG_scaled(a) * (Lanczos::lanczos_sum_expG_scaled(b) / Lanczos::lanczos_sum_expG_scaled(c));
  82. T ambh = a - 0.5f - b;
  83. if((fabs(b * ambh) < (cgh * 100)) && (a > 100))
  84. {
  85. // Special case where the base of the power term is close to 1
  86. // compute (1+x)^y instead:
  87. result *= exp(ambh * boost::math::log1p(-b / cgh, pol));
  88. }
  89. else
  90. {
  91. result *= pow(agh / cgh, a - T(0.5) - b);
  92. }
  93. if(cgh > 1e10f)
  94. // this avoids possible overflow, but appears to be marginally less accurate:
  95. result *= pow((agh / cgh) * (bgh / cgh), b);
  96. else
  97. result *= pow((agh * bgh) / (cgh * cgh), b);
  98. result *= sqrt(boost::math::constants::e<T>() / bgh);
  99. // If a and b were originally less than 1 we need to scale the result:
  100. result *= prefix;
  101. return result;
  102. } // template <class T, class Lanczos> beta_imp(T a, T b, const Lanczos&)
  103. //
  104. // Generic implementation of Beta(a,b) without Lanczos approximation support
  105. // (Caution this is slow!!!):
  106. //
  107. template <class T, class Policy>
  108. T beta_imp(T a, T b, const lanczos::undefined_lanczos& l, const Policy& pol)
  109. {
  110. BOOST_MATH_STD_USING
  111. if(a <= 0)
  112. return policies::raise_domain_error<T>("boost::math::beta<%1%>(%1%,%1%)", "The arguments to the beta function must be greater than zero (got a=%1%).", a, pol);
  113. if(b <= 0)
  114. return policies::raise_domain_error<T>("boost::math::beta<%1%>(%1%,%1%)", "The arguments to the beta function must be greater than zero (got b=%1%).", b, pol);
  115. const T c = a + b;
  116. // Special cases:
  117. if ((c == a) && (b < tools::epsilon<T>()))
  118. return 1 / b;
  119. else if ((c == b) && (a < tools::epsilon<T>()))
  120. return 1 / a;
  121. if (b == 1)
  122. return 1 / a;
  123. else if (a == 1)
  124. return 1 / b;
  125. else if (c < tools::epsilon<T>())
  126. {
  127. T result = c / a;
  128. result /= b;
  129. return result;
  130. }
  131. // Regular cases start here:
  132. const T min_sterling = minimum_argument_for_bernoulli_recursion<T>();
  133. long shift_a = 0;
  134. long shift_b = 0;
  135. if(a < min_sterling)
  136. shift_a = 1 + ltrunc(min_sterling - a);
  137. if(b < min_sterling)
  138. shift_b = 1 + ltrunc(min_sterling - b);
  139. long shift_c = shift_a + shift_b;
  140. if ((shift_a == 0) && (shift_b == 0))
  141. {
  142. return pow(a / c, a) * pow(b / c, b) * scaled_tgamma_no_lanczos(a, pol) * scaled_tgamma_no_lanczos(b, pol) / scaled_tgamma_no_lanczos(c, pol);
  143. }
  144. else if ((a < 1) && (b < 1))
  145. {
  146. return boost::math::tgamma(a, pol) * (boost::math::tgamma(b, pol) / boost::math::tgamma(c));
  147. }
  148. else if(a < 1)
  149. return boost::math::tgamma(a, pol) * boost::math::tgamma_delta_ratio(b, a, pol);
  150. else if(b < 1)
  151. return boost::math::tgamma(b, pol) * boost::math::tgamma_delta_ratio(a, b, pol);
  152. else
  153. {
  154. T result = beta_imp(T(a + shift_a), T(b + shift_b), l, pol);
  155. //
  156. // Recursion:
  157. //
  158. for (long i = 0; i < shift_c; ++i)
  159. {
  160. result *= c + i;
  161. if (i < shift_a)
  162. result /= a + i;
  163. if (i < shift_b)
  164. result /= b + i;
  165. }
  166. return result;
  167. }
  168. } // template <class T>T beta_imp(T a, T b, const lanczos::undefined_lanczos& l)
  169. //
  170. // Compute the leading power terms in the incomplete Beta:
  171. //
  172. // (x^a)(y^b)/Beta(a,b) when normalised, and
  173. // (x^a)(y^b) otherwise.
  174. //
  175. // Almost all of the error in the incomplete beta comes from this
  176. // function: particularly when a and b are large. Computing large
  177. // powers are *hard* though, and using logarithms just leads to
  178. // horrendous cancellation errors.
  179. //
  180. template <class T, class Lanczos, class Policy>
  181. T ibeta_power_terms(T a,
  182. T b,
  183. T x,
  184. T y,
  185. const Lanczos&,
  186. bool normalised,
  187. const Policy& pol,
  188. T prefix = 1,
  189. const char* function = "boost::math::ibeta<%1%>(%1%, %1%, %1%)")
  190. {
  191. BOOST_MATH_STD_USING
  192. if(!normalised)
  193. {
  194. // can we do better here?
  195. return pow(x, a) * pow(y, b);
  196. }
  197. T result; // LCOV_EXCL_LINE
  198. T c = a + b;
  199. // combine power terms with Lanczos approximation:
  200. T agh = static_cast<T>(a + Lanczos::g() - 0.5f);
  201. T bgh = static_cast<T>(b + Lanczos::g() - 0.5f);
  202. T cgh = static_cast<T>(c + Lanczos::g() - 0.5f);
  203. if ((a < tools::min_value<T>()) || (b < tools::min_value<T>()))
  204. result = 0; // denominator overflows in this case
  205. else
  206. result = Lanczos::lanczos_sum_expG_scaled(c) / (Lanczos::lanczos_sum_expG_scaled(a) * Lanczos::lanczos_sum_expG_scaled(b));
  207. result *= prefix;
  208. // combine with the leftover terms from the Lanczos approximation:
  209. result *= sqrt(bgh / boost::math::constants::e<T>());
  210. result *= sqrt(agh / cgh);
  211. // l1 and l2 are the base of the exponents minus one:
  212. T l1 = (x * b - y * agh) / agh;
  213. T l2 = (y * a - x * bgh) / bgh;
  214. if(((std::min)(fabs(l1), fabs(l2)) < 0.2))
  215. {
  216. // when the base of the exponent is very near 1 we get really
  217. // gross errors unless extra care is taken:
  218. if((l1 * l2 > 0) || ((std::min)(a, b) < 1))
  219. {
  220. //
  221. // This first branch handles the simple cases where either:
  222. //
  223. // * The two power terms both go in the same direction
  224. // (towards zero or towards infinity). In this case if either
  225. // term overflows or underflows, then the product of the two must
  226. // do so also.
  227. // *Alternatively if one exponent is less than one, then we
  228. // can't productively use it to eliminate overflow or underflow
  229. // from the other term. Problems with spurious overflow/underflow
  230. // can't be ruled out in this case, but it is *very* unlikely
  231. // since one of the power terms will evaluate to a number close to 1.
  232. //
  233. if(fabs(l1) < 0.1)
  234. {
  235. result *= exp(a * boost::math::log1p(l1, pol));
  236. BOOST_MATH_INSTRUMENT_VARIABLE(result);
  237. }
  238. else
  239. {
  240. result *= pow((x * cgh) / agh, a);
  241. BOOST_MATH_INSTRUMENT_VARIABLE(result);
  242. }
  243. if(fabs(l2) < 0.1)
  244. {
  245. result *= exp(b * boost::math::log1p(l2, pol));
  246. BOOST_MATH_INSTRUMENT_VARIABLE(result);
  247. }
  248. else
  249. {
  250. result *= pow((y * cgh) / bgh, b);
  251. BOOST_MATH_INSTRUMENT_VARIABLE(result);
  252. }
  253. }
  254. else if((std::max)(fabs(l1), fabs(l2)) < 0.5)
  255. {
  256. //
  257. // Both exponents are near one and both the exponents are
  258. // greater than one and further these two
  259. // power terms tend in opposite directions (one towards zero,
  260. // the other towards infinity), so we have to combine the terms
  261. // to avoid any risk of overflow or underflow.
  262. //
  263. // We do this by moving one power term inside the other, we have:
  264. //
  265. // (1 + l1)^a * (1 + l2)^b
  266. // = ((1 + l1)*(1 + l2)^(b/a))^a
  267. // = (1 + l1 + l3 + l1*l3)^a ; l3 = (1 + l2)^(b/a) - 1
  268. // = exp((b/a) * log(1 + l2)) - 1
  269. //
  270. // The tricky bit is deciding which term to move inside :-)
  271. // By preference we move the larger term inside, so that the
  272. // size of the largest exponent is reduced. However, that can
  273. // only be done as long as l3 (see above) is also small.
  274. //
  275. bool small_a = a < b;
  276. T ratio = b / a;
  277. if((small_a && (ratio * l2 < 0.1)) || (!small_a && (l1 / ratio > 0.1)))
  278. {
  279. T l3 = boost::math::expm1(ratio * boost::math::log1p(l2, pol), pol);
  280. l3 = l1 + l3 + l3 * l1;
  281. l3 = a * boost::math::log1p(l3, pol);
  282. result *= exp(l3);
  283. BOOST_MATH_INSTRUMENT_VARIABLE(result);
  284. }
  285. else
  286. {
  287. T l3 = boost::math::expm1(boost::math::log1p(l1, pol) / ratio, pol);
  288. l3 = l2 + l3 + l3 * l2;
  289. l3 = b * boost::math::log1p(l3, pol);
  290. result *= exp(l3);
  291. BOOST_MATH_INSTRUMENT_VARIABLE(result);
  292. }
  293. }
  294. else if(fabs(l1) < fabs(l2))
  295. {
  296. // First base near 1 only:
  297. T l = a * boost::math::log1p(l1, pol)
  298. + b * log((y * cgh) / bgh);
  299. if((l <= tools::log_min_value<T>()) || (l >= tools::log_max_value<T>()))
  300. {
  301. l += log(result);
  302. if(l >= tools::log_max_value<T>())
  303. return policies::raise_overflow_error<T>(function, nullptr, pol); // LCOV_EXCL_LINE we can probably never get here, probably!
  304. result = exp(l);
  305. }
  306. else
  307. result *= exp(l);
  308. BOOST_MATH_INSTRUMENT_VARIABLE(result);
  309. }
  310. else
  311. {
  312. // Second base near 1 only:
  313. T l = b * boost::math::log1p(l2, pol)
  314. + a * log((x * cgh) / agh);
  315. if((l <= tools::log_min_value<T>()) || (l >= tools::log_max_value<T>()))
  316. {
  317. l += log(result);
  318. if(l >= tools::log_max_value<T>())
  319. return policies::raise_overflow_error<T>(function, nullptr, pol); // LCOV_EXCL_LINE we can probably never get here, probably!
  320. result = exp(l);
  321. }
  322. else
  323. result *= exp(l);
  324. BOOST_MATH_INSTRUMENT_VARIABLE(result);
  325. }
  326. }
  327. else
  328. {
  329. // general case:
  330. T b1 = (x * cgh) / agh;
  331. T b2 = (y * cgh) / bgh;
  332. l1 = a * log(b1);
  333. l2 = b * log(b2);
  334. BOOST_MATH_INSTRUMENT_VARIABLE(b1);
  335. BOOST_MATH_INSTRUMENT_VARIABLE(b2);
  336. BOOST_MATH_INSTRUMENT_VARIABLE(l1);
  337. BOOST_MATH_INSTRUMENT_VARIABLE(l2);
  338. if((l1 >= tools::log_max_value<T>())
  339. || (l1 <= tools::log_min_value<T>())
  340. || (l2 >= tools::log_max_value<T>())
  341. || (l2 <= tools::log_min_value<T>())
  342. )
  343. {
  344. // Oops, under/overflow, sidestep if we can:
  345. if(a < b)
  346. {
  347. T p1 = pow(b2, b / a);
  348. T l3 = (b1 != 0) && (p1 != 0) ? (a * (log(b1) + log(p1))) : tools::max_value<T>(); // arbitrary large value if the logs would fail!
  349. if((l3 < tools::log_max_value<T>())
  350. && (l3 > tools::log_min_value<T>()))
  351. {
  352. result *= pow(p1 * b1, a);
  353. }
  354. else
  355. {
  356. l2 += l1 + log(result);
  357. if(l2 >= tools::log_max_value<T>())
  358. return policies::raise_overflow_error<T>(function, nullptr, pol); // LCOV_EXCL_LINE we can probably never get here, probably!
  359. result = exp(l2);
  360. }
  361. }
  362. else
  363. {
  364. // This protects against spurious overflow in a/b:
  365. T p1 = (b1 < 1) && (b < 1) && (tools::max_value<T>() * b < a) ? static_cast<T>(0) : static_cast<T>(pow(b1, a / b));
  366. T l3 = (p1 != 0) && (b2 != 0) ? (log(p1) + log(b2)) * b : tools::max_value<T>(); // arbitrary large value if the logs would fail!
  367. if((l3 < tools::log_max_value<T>())
  368. && (l3 > tools::log_min_value<T>()))
  369. {
  370. result *= pow(p1 * b2, b);
  371. }
  372. else if(result != 0) // we can elude the calculation below if we're already going to be zero
  373. {
  374. l2 += l1 + log(result);
  375. if(l2 >= tools::log_max_value<T>())
  376. return policies::raise_overflow_error<T>(function, nullptr, pol); // LCOV_EXCL_LINE we can probably never get here, probably!
  377. result = exp(l2);
  378. }
  379. }
  380. BOOST_MATH_INSTRUMENT_VARIABLE(result);
  381. }
  382. else
  383. {
  384. // finally the normal case:
  385. result *= pow(b1, a) * pow(b2, b);
  386. BOOST_MATH_INSTRUMENT_VARIABLE(result);
  387. }
  388. }
  389. BOOST_MATH_INSTRUMENT_VARIABLE(result);
  390. if (0 == result)
  391. {
  392. if ((a > 1) && (x == 0))
  393. return result; // true zero LCOV_EXCL_LINE we can probably never get here
  394. if ((b > 1) && (y == 0))
  395. return result; // true zero LCOV_EXCL_LINE we can probably never get here
  396. return boost::math::policies::raise_underflow_error<T>(function, nullptr, pol);
  397. }
  398. return result;
  399. }
  400. //
  401. // Compute the leading power terms in the incomplete Beta:
  402. //
  403. // (x^a)(y^b)/Beta(a,b) when normalised, and
  404. // (x^a)(y^b) otherwise.
  405. //
  406. // Almost all of the error in the incomplete beta comes from this
  407. // function: particularly when a and b are large. Computing large
  408. // powers are *hard* though, and using logarithms just leads to
  409. // horrendous cancellation errors.
  410. //
  411. // This version is generic, slow, and does not use the Lanczos approximation.
  412. //
  413. template <class T, class Policy>
  414. T ibeta_power_terms(T a,
  415. T b,
  416. T x,
  417. T y,
  418. const boost::math::lanczos::undefined_lanczos& l,
  419. bool normalised,
  420. const Policy& pol,
  421. T prefix = 1,
  422. const char* = "boost::math::ibeta<%1%>(%1%, %1%, %1%)")
  423. {
  424. BOOST_MATH_STD_USING
  425. if(!normalised)
  426. {
  427. return prefix * pow(x, a) * pow(y, b);
  428. }
  429. T c = a + b;
  430. const T min_sterling = minimum_argument_for_bernoulli_recursion<T>();
  431. long shift_a = 0;
  432. long shift_b = 0;
  433. if (a < min_sterling)
  434. shift_a = 1 + ltrunc(min_sterling - a);
  435. if (b < min_sterling)
  436. shift_b = 1 + ltrunc(min_sterling - b);
  437. if ((shift_a == 0) && (shift_b == 0))
  438. {
  439. T power1, power2;
  440. bool need_logs = false;
  441. if (a < b)
  442. {
  443. BOOST_MATH_IF_CONSTEXPR(std::numeric_limits<T>::has_infinity)
  444. {
  445. power1 = pow((x * y * c * c) / (a * b), a);
  446. power2 = pow((y * c) / b, b - a);
  447. }
  448. else
  449. {
  450. // We calculate these logs purely so we can check for overflow in the power functions
  451. T l1 = log((x * y * c * c) / (a * b));
  452. T l2 = log((y * c) / b);
  453. if ((l1 * a > tools::log_min_value<T>()) && (l1 * a < tools::log_max_value<T>()) && (l2 * (b - a) < tools::log_max_value<T>()) && (l2 * (b - a) > tools::log_min_value<T>()))
  454. {
  455. power1 = pow((x * y * c * c) / (a * b), a);
  456. power2 = pow((y * c) / b, b - a);
  457. }
  458. else
  459. {
  460. need_logs = true;
  461. }
  462. }
  463. }
  464. else
  465. {
  466. BOOST_MATH_IF_CONSTEXPR(std::numeric_limits<T>::has_infinity)
  467. {
  468. power1 = pow((x * y * c * c) / (a * b), b);
  469. power2 = pow((x * c) / a, a - b);
  470. }
  471. else
  472. {
  473. // We calculate these logs purely so we can check for overflow in the power functions
  474. T l1 = log((x * y * c * c) / (a * b)) * b;
  475. T l2 = log((x * c) / a) * (a - b);
  476. if ((l1 * a > tools::log_min_value<T>()) && (l1 * a < tools::log_max_value<T>()) && (l2 * (b - a) < tools::log_max_value<T>()) && (l2 * (b - a) > tools::log_min_value<T>()))
  477. {
  478. power1 = pow((x * y * c * c) / (a * b), b);
  479. power2 = pow((x * c) / a, a - b);
  480. }
  481. else
  482. need_logs = true;
  483. }
  484. }
  485. BOOST_MATH_IF_CONSTEXPR(std::numeric_limits<T>::has_infinity)
  486. {
  487. if (!(boost::math::isnormal)(power1) || !(boost::math::isnormal)(power2))
  488. {
  489. need_logs = true;
  490. }
  491. }
  492. if (need_logs)
  493. {
  494. //
  495. // We want:
  496. //
  497. // (xc / a)^a (yc / b)^b
  498. //
  499. // But we know that one or other term will over / underflow and combining the logs will be next to useless as that will cause significant cancellation.
  500. // If we assume b > a and express z ^ b as(z ^ b / a) ^ a with z = (yc / b) then we can move one power term inside the other :
  501. //
  502. // ((xc / a) * (yc / b)^(b / a))^a
  503. //
  504. // However, we're not quite there yet, as the term being exponentiated is quite likely to be close to unity, so let:
  505. //
  506. // xc / a = 1 + (xb - ya) / a
  507. //
  508. // analogously let :
  509. //
  510. // 1 + p = (yc / b) ^ (b / a) = 1 + expm1((b / a) * log1p((ya - xb) / b))
  511. //
  512. // so putting the two together we have :
  513. //
  514. // exp(a * log1p((xb - ya) / a + p + p(xb - ya) / a))
  515. //
  516. // Analogously, when a > b we can just swap all the terms around.
  517. //
  518. // Finally, there are a few cases (x or y is unity) when the above logic can't be used
  519. // or where there is no logarithmic cancellation and accuracy is better just using
  520. // the regular formula:
  521. //
  522. T xc_a = x * c / a;
  523. T yc_b = y * c / b;
  524. if ((x == 1) || (y == 1) || (fabs(xc_a - 1) > 0.25) || (fabs(yc_b - 1) > 0.25))
  525. {
  526. // The above logic fails, the result is almost certainly zero:
  527. power1 = exp(log(xc_a) * a + log(yc_b) * b);
  528. power2 = 1;
  529. }
  530. else if (b > a)
  531. {
  532. T p = boost::math::expm1((b / a) * boost::math::log1p((y * a - x * b) / b));
  533. power1 = exp(a * boost::math::log1p((x * b - y * a) / a + p * (x * c / a)));
  534. power2 = 1;
  535. }
  536. else
  537. {
  538. T p = boost::math::expm1((a / b) * boost::math::log1p((x * b - y * a) / a));
  539. power1 = exp(b * boost::math::log1p((y * a - x * b) / b + p * (y * c / b)));
  540. power2 = 1;
  541. }
  542. }
  543. return prefix * power1 * power2 * scaled_tgamma_no_lanczos(c, pol) / (scaled_tgamma_no_lanczos(a, pol) * scaled_tgamma_no_lanczos(b, pol));
  544. }
  545. T power1 = pow(x, a);
  546. T power2 = pow(y, b);
  547. T bet = beta_imp(a, b, l, pol);
  548. if(!(boost::math::isnormal)(power1) || !(boost::math::isnormal)(power2) || !(boost::math::isnormal)(bet))
  549. {
  550. int shift_c = shift_a + shift_b;
  551. T result = ibeta_power_terms(T(a + shift_a), T(b + shift_b), x, y, l, normalised, pol, prefix);
  552. if ((boost::math::isnormal)(result))
  553. {
  554. for (int i = 0; i < shift_c; ++i)
  555. {
  556. result /= c + i;
  557. if (i < shift_a)
  558. {
  559. result *= a + i;
  560. result /= x;
  561. }
  562. if (i < shift_b)
  563. {
  564. result *= b + i;
  565. result /= y;
  566. }
  567. }
  568. return prefix * result;
  569. }
  570. else
  571. {
  572. T log_result = log(x) * a + log(y) * b + log(prefix);
  573. if ((boost::math::isnormal)(bet))
  574. log_result -= log(bet);
  575. else
  576. log_result += boost::math::lgamma(c, pol) - boost::math::lgamma(a, pol) - boost::math::lgamma(b, pol);
  577. return exp(log_result);
  578. }
  579. }
  580. return prefix * power1 * (power2 / bet);
  581. }
  582. //
  583. // Series approximation to the incomplete beta:
  584. //
  585. template <class T>
  586. struct ibeta_series_t
  587. {
  588. typedef T result_type;
  589. ibeta_series_t(T a_, T b_, T x_, T mult) : result(mult), x(x_), apn(a_), poch(1-b_), n(1) {}
  590. T operator()()
  591. {
  592. T r = result / apn;
  593. apn += 1;
  594. result *= poch * x / n;
  595. ++n;
  596. poch += 1;
  597. return r;
  598. }
  599. private:
  600. T result, x, apn, poch;
  601. int n;
  602. };
  603. template <class T, class Lanczos, class Policy>
  604. T ibeta_series(T a, T b, T x, T s0, const Lanczos&, bool normalised, T* p_derivative, T y, const Policy& pol)
  605. {
  606. BOOST_MATH_STD_USING
  607. T result;
  608. BOOST_MATH_ASSERT((p_derivative == 0) || normalised);
  609. if(normalised)
  610. {
  611. T c = a + b;
  612. // incomplete beta power term, combined with the Lanczos approximation:
  613. T agh = static_cast<T>(a + Lanczos::g() - 0.5f);
  614. T bgh = static_cast<T>(b + Lanczos::g() - 0.5f);
  615. T cgh = static_cast<T>(c + Lanczos::g() - 0.5f);
  616. if ((a < tools::min_value<T>()) || (b < tools::min_value<T>()))
  617. result = 0; // denorms cause overflow in the Lanzos series, result will be zero anyway
  618. else
  619. result = Lanczos::lanczos_sum_expG_scaled(c) / (Lanczos::lanczos_sum_expG_scaled(a) * Lanczos::lanczos_sum_expG_scaled(b));
  620. if (!(boost::math::isfinite)(result))
  621. result = 0; // LCOV_EXCL_LINE we can probably never get here, covered already above?
  622. T l1 = log(cgh / bgh) * (b - 0.5f);
  623. T l2 = log(x * cgh / agh) * a;
  624. //
  625. // Check for over/underflow in the power terms:
  626. //
  627. if((l1 > tools::log_min_value<T>())
  628. && (l1 < tools::log_max_value<T>())
  629. && (l2 > tools::log_min_value<T>())
  630. && (l2 < tools::log_max_value<T>()))
  631. {
  632. if(a * b < bgh * 10)
  633. result *= exp((b - 0.5f) * boost::math::log1p(a / bgh, pol));
  634. else
  635. result *= pow(cgh / bgh, T(b - T(0.5)));
  636. result *= pow(x * cgh / agh, a);
  637. result *= sqrt(agh / boost::math::constants::e<T>());
  638. if(p_derivative)
  639. {
  640. *p_derivative = result * pow(y, b);
  641. BOOST_MATH_ASSERT(*p_derivative >= 0);
  642. }
  643. }
  644. else
  645. {
  646. //
  647. // Oh dear, we need logs, and this *will* cancel:
  648. //
  649. if (result != 0) // elude calculation when result will be zero.
  650. {
  651. result = log(result) + l1 + l2 + (log(agh) - 1) / 2;
  652. if (p_derivative)
  653. *p_derivative = exp(result + b * log(y));
  654. result = exp(result);
  655. }
  656. }
  657. }
  658. else
  659. {
  660. // Non-normalised, just compute the power:
  661. result = pow(x, a);
  662. }
  663. if(result < tools::min_value<T>())
  664. return s0; // Safeguard: series can't cope with denorms.
  665. ibeta_series_t<T> s(a, b, x, result);
  666. std::uintmax_t max_iter = policies::get_max_series_iterations<Policy>();
  667. result = boost::math::tools::sum_series(s, boost::math::policies::get_epsilon<T, Policy>(), max_iter, s0);
  668. policies::check_series_iterations<T>("boost::math::ibeta<%1%>(%1%, %1%, %1%) in ibeta_series (with lanczos)", max_iter, pol);
  669. return result;
  670. }
  671. //
  672. // Incomplete Beta series again, this time without Lanczos support:
  673. //
  674. template <class T, class Policy>
  675. T ibeta_series(T a, T b, T x, T s0, const boost::math::lanczos::undefined_lanczos& l, bool normalised, T* p_derivative, T y, const Policy& pol)
  676. {
  677. BOOST_MATH_STD_USING
  678. T result;
  679. BOOST_MATH_ASSERT((p_derivative == 0) || normalised);
  680. if(normalised)
  681. {
  682. const T min_sterling = minimum_argument_for_bernoulli_recursion<T>();
  683. long shift_a = 0;
  684. long shift_b = 0;
  685. if (a < min_sterling)
  686. shift_a = 1 + ltrunc(min_sterling - a);
  687. if (b < min_sterling)
  688. shift_b = 1 + ltrunc(min_sterling - b);
  689. T c = a + b;
  690. if ((shift_a == 0) && (shift_b == 0))
  691. {
  692. result = pow(x * c / a, a) * pow(c / b, b) * scaled_tgamma_no_lanczos(c, pol) / (scaled_tgamma_no_lanczos(a, pol) * scaled_tgamma_no_lanczos(b, pol));
  693. }
  694. else if ((a < 1) && (b > 1))
  695. result = pow(x, a) / (boost::math::tgamma(a, pol) * boost::math::tgamma_delta_ratio(b, a, pol));
  696. else
  697. {
  698. T power = pow(x, a);
  699. T bet = beta_imp(a, b, l, pol);
  700. if (!(boost::math::isnormal)(power) || !(boost::math::isnormal)(bet))
  701. {
  702. result = exp(a * log(x) + boost::math::lgamma(c, pol) - boost::math::lgamma(a, pol) - boost::math::lgamma(b, pol));
  703. }
  704. else
  705. result = power / bet;
  706. }
  707. if(p_derivative)
  708. {
  709. *p_derivative = result * pow(y, b);
  710. BOOST_MATH_ASSERT(*p_derivative >= 0);
  711. }
  712. }
  713. else
  714. {
  715. // Non-normalised, just compute the power:
  716. result = pow(x, a);
  717. }
  718. if(result < tools::min_value<T>())
  719. return s0; // Safeguard: series can't cope with denorms.
  720. ibeta_series_t<T> s(a, b, x, result);
  721. std::uintmax_t max_iter = policies::get_max_series_iterations<Policy>();
  722. result = boost::math::tools::sum_series(s, boost::math::policies::get_epsilon<T, Policy>(), max_iter, s0);
  723. policies::check_series_iterations<T>("boost::math::ibeta<%1%>(%1%, %1%, %1%) in ibeta_series (without lanczos)", max_iter, pol);
  724. return result;
  725. }
  726. //
  727. // Continued fraction for the incomplete beta:
  728. //
  729. template <class T>
  730. struct ibeta_fraction2_t
  731. {
  732. typedef std::pair<T, T> result_type;
  733. ibeta_fraction2_t(T a_, T b_, T x_, T y_) : a(a_), b(b_), x(x_), y(y_), m(0) {}
  734. result_type operator()()
  735. {
  736. T aN = (a + m - 1) * (a + b + m - 1) * m * (b - m) * x * x;
  737. T denom = (a + 2 * m - 1);
  738. aN /= denom * denom;
  739. T bN = static_cast<T>(m);
  740. bN += (m * (b - m) * x) / (a + 2*m - 1);
  741. bN += ((a + m) * (a * y - b * x + 1 + m *(2 - x))) / (a + 2*m + 1);
  742. ++m;
  743. return std::make_pair(aN, bN);
  744. }
  745. private:
  746. T a, b, x, y;
  747. int m;
  748. };
  749. //
  750. // Evaluate the incomplete beta via the continued fraction representation:
  751. //
  752. template <class T, class Policy>
  753. inline T ibeta_fraction2(T a, T b, T x, T y, const Policy& pol, bool normalised, T* p_derivative)
  754. {
  755. typedef typename lanczos::lanczos<T, Policy>::type lanczos_type;
  756. BOOST_MATH_STD_USING
  757. T result = ibeta_power_terms(a, b, x, y, lanczos_type(), normalised, pol);
  758. if(p_derivative)
  759. {
  760. *p_derivative = result;
  761. BOOST_MATH_ASSERT(*p_derivative >= 0);
  762. }
  763. if(result == 0)
  764. return result;
  765. ibeta_fraction2_t<T> f(a, b, x, y);
  766. T fract = boost::math::tools::continued_fraction_b(f, boost::math::policies::get_epsilon<T, Policy>());
  767. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  768. BOOST_MATH_INSTRUMENT_VARIABLE(result);
  769. return result / fract;
  770. }
  771. //
  772. // Computes the difference between ibeta(a,b,x) and ibeta(a+k,b,x):
  773. //
  774. template <class T, class Policy>
  775. T ibeta_a_step(T a, T b, T x, T y, int k, const Policy& pol, bool normalised, T* p_derivative)
  776. {
  777. typedef typename lanczos::lanczos<T, Policy>::type lanczos_type;
  778. BOOST_MATH_INSTRUMENT_VARIABLE(k);
  779. T prefix = ibeta_power_terms(a, b, x, y, lanczos_type(), normalised, pol);
  780. if(p_derivative)
  781. {
  782. *p_derivative = prefix;
  783. BOOST_MATH_ASSERT(*p_derivative >= 0);
  784. }
  785. prefix /= a;
  786. if(prefix == 0)
  787. return prefix;
  788. T sum = 1;
  789. T term = 1;
  790. // series summation from 0 to k-1:
  791. for(int i = 0; i < k-1; ++i)
  792. {
  793. term *= (a+b+i) * x / (a+i+1);
  794. sum += term;
  795. }
  796. prefix *= sum;
  797. return prefix;
  798. }
  799. //
  800. // This function is only needed for the non-regular incomplete beta,
  801. // it computes the delta in:
  802. // beta(a,b,x) = prefix + delta * beta(a+k,b,x)
  803. // it is currently only called for small k.
  804. //
  805. template <class T>
  806. inline T rising_factorial_ratio(T a, T b, int k)
  807. {
  808. // calculate:
  809. // (a)(a+1)(a+2)...(a+k-1)
  810. // _______________________
  811. // (b)(b+1)(b+2)...(b+k-1)
  812. // This is only called with small k, for large k
  813. // it is grossly inefficient, do not use outside it's
  814. // intended purpose!!!
  815. BOOST_MATH_INSTRUMENT_VARIABLE(k);
  816. BOOST_MATH_ASSERT(k > 0);
  817. T result = 1;
  818. for(int i = 0; i < k; ++i)
  819. result *= (a+i) / (b+i);
  820. return result;
  821. }
  822. //
  823. // Routine for a > 15, b < 1
  824. //
  825. // Begin by figuring out how large our table of Pn's should be,
  826. // quoted accuracies are "guesstimates" based on empirical observation.
  827. // Note that the table size should never exceed the size of our
  828. // tables of factorials.
  829. //
  830. template <class T>
  831. struct Pn_size
  832. {
  833. // This is likely to be enough for ~35-50 digit accuracy
  834. // but it's hard to quantify exactly:
  835. static constexpr unsigned value =
  836. ::boost::math::max_factorial<T>::value >= 100 ? 50
  837. : ::boost::math::max_factorial<T>::value >= ::boost::math::max_factorial<double>::value ? 30
  838. : ::boost::math::max_factorial<T>::value >= ::boost::math::max_factorial<float>::value ? 15 : 1;
  839. static_assert(::boost::math::max_factorial<T>::value >= ::boost::math::max_factorial<float>::value, "Type does not provide for 35-50 digits of accuracy.");
  840. };
  841. template <>
  842. struct Pn_size<float>
  843. {
  844. static constexpr unsigned value = 15; // ~8-15 digit accuracy
  845. static_assert(::boost::math::max_factorial<float>::value >= 30, "Type does not provide for 8-15 digits of accuracy.");
  846. };
  847. template <>
  848. struct Pn_size<double>
  849. {
  850. static constexpr unsigned value = 30; // 16-20 digit accuracy
  851. static_assert(::boost::math::max_factorial<double>::value >= 60, "Type does not provide for 16-20 digits of accuracy.");
  852. };
  853. template <>
  854. struct Pn_size<long double>
  855. {
  856. static constexpr unsigned value = 50; // ~35-50 digit accuracy
  857. static_assert(::boost::math::max_factorial<long double>::value >= 100, "Type does not provide for ~35-50 digits of accuracy");
  858. };
  859. template <class T, class Policy>
  860. T beta_small_b_large_a_series(T a, T b, T x, T y, T s0, T mult, const Policy& pol, bool normalised)
  861. {
  862. typedef typename lanczos::lanczos<T, Policy>::type lanczos_type;
  863. BOOST_MATH_STD_USING
  864. //
  865. // This is DiDonato and Morris's BGRAT routine, see Eq's 9 through 9.6.
  866. //
  867. // Some values we'll need later, these are Eq 9.1:
  868. //
  869. T bm1 = b - 1;
  870. T t = a + bm1 / 2;
  871. T lx, u; // LCOV_EXCL_LINE
  872. if(y < 0.35)
  873. lx = boost::math::log1p(-y, pol);
  874. else
  875. lx = log(x);
  876. u = -t * lx;
  877. // and from from 9.2:
  878. T prefix; // LCOV_EXCL_LINE
  879. T h = regularised_gamma_prefix(b, u, pol, lanczos_type());
  880. if(h <= tools::min_value<T>())
  881. return s0;
  882. if(normalised)
  883. {
  884. prefix = h / boost::math::tgamma_delta_ratio(a, b, pol);
  885. prefix /= pow(t, b);
  886. }
  887. else
  888. {
  889. prefix = full_igamma_prefix(b, u, pol) / pow(t, b);
  890. }
  891. prefix *= mult;
  892. //
  893. // now we need the quantity Pn, unfortunately this is computed
  894. // recursively, and requires a full history of all the previous values
  895. // so no choice but to declare a big table and hope it's big enough...
  896. //
  897. T p[ ::boost::math::detail::Pn_size<T>::value ] = { 1 }; // see 9.3.
  898. //
  899. // Now an initial value for J, see 9.6:
  900. //
  901. T j = boost::math::gamma_q(b, u, pol) / h;
  902. //
  903. // Now we can start to pull things together and evaluate the sum in Eq 9:
  904. //
  905. T sum = s0 + prefix * j; // Value at N = 0
  906. // some variables we'll need:
  907. unsigned tnp1 = 1; // 2*N+1
  908. T lx2 = lx / 2;
  909. lx2 *= lx2;
  910. T lxp = 1;
  911. T t4 = 4 * t * t;
  912. T b2n = b;
  913. for(unsigned n = 1; n < sizeof(p)/sizeof(p[0]); ++n)
  914. {
  915. /*
  916. // debugging code, enable this if you want to determine whether
  917. // the table of Pn's is large enough...
  918. //
  919. static int max_count = 2;
  920. if(n > max_count)
  921. {
  922. max_count = n;
  923. std::cerr << "Max iterations in BGRAT was " << n << std::endl;
  924. }
  925. */
  926. //
  927. // begin by evaluating the next Pn from Eq 9.4:
  928. //
  929. tnp1 += 2;
  930. p[n] = 0;
  931. T mbn = b - n;
  932. unsigned tmp1 = 3;
  933. for(unsigned m = 1; m < n; ++m)
  934. {
  935. mbn = m * b - n;
  936. p[n] += mbn * p[n-m] / boost::math::unchecked_factorial<T>(tmp1);
  937. tmp1 += 2;
  938. }
  939. p[n] /= n;
  940. p[n] += bm1 / boost::math::unchecked_factorial<T>(tnp1);
  941. //
  942. // Now we want Jn from Jn-1 using Eq 9.6:
  943. //
  944. j = (b2n * (b2n + 1) * j + (u + b2n + 1) * lxp) / t4;
  945. lxp *= lx2;
  946. b2n += 2;
  947. //
  948. // pull it together with Eq 9:
  949. //
  950. T r = prefix * p[n] * j;
  951. sum += r;
  952. // r is always small:
  953. BOOST_MATH_ASSERT(tools::max_value<T>() * tools::epsilon<T>() > fabs(r));
  954. if(fabs(r / tools::epsilon<T>()) < fabs(sum))
  955. break;
  956. }
  957. return sum;
  958. } // template <class T, class Lanczos>T beta_small_b_large_a_series(T a, T b, T x, T y, T s0, T mult, const Lanczos& l, bool normalised)
  959. //
  960. // For integer arguments we can relate the incomplete beta to the
  961. // complement of the binomial distribution cdf and use this finite sum.
  962. //
  963. template <class T, class Policy>
  964. T binomial_ccdf(T n, T k, T x, T y, const Policy& pol)
  965. {
  966. BOOST_MATH_STD_USING // ADL of std names
  967. T result = pow(x, n);
  968. if(result > tools::min_value<T>())
  969. {
  970. T term = result;
  971. for(unsigned i = itrunc(T(n - 1)); i > k; --i)
  972. {
  973. term *= ((i + 1) * y) / ((n - i) * x);
  974. result += term;
  975. }
  976. }
  977. else
  978. {
  979. // First term underflows so we need to start at the mode of the
  980. // distribution and work outwards:
  981. int start = itrunc(n * x);
  982. if(start <= k + 1)
  983. start = itrunc(k + 2);
  984. result = static_cast<T>(pow(x, T(start)) * pow(y, n - T(start)) * boost::math::binomial_coefficient<T>(itrunc(n), itrunc(start), pol));
  985. if(result == 0)
  986. {
  987. // OK, starting slightly above the mode didn't work,
  988. // we'll have to sum the terms the old fashioned way.
  989. // Very hard to get here, possibly only when exponent
  990. // range is very limited (as with type float):
  991. // LCOV_EXCL_START
  992. for(unsigned i = start - 1; i > k; --i)
  993. {
  994. result += static_cast<T>(pow(x, static_cast<T>(i)) * pow(y, n - i) * boost::math::binomial_coefficient<T>(itrunc(n), itrunc(i), pol));
  995. }
  996. // LCOV_EXCL_STOP
  997. }
  998. else
  999. {
  1000. T term = result;
  1001. T start_term = result;
  1002. for(unsigned i = start - 1; i > k; --i)
  1003. {
  1004. term *= ((i + 1) * y) / ((n - i) * x);
  1005. result += term;
  1006. }
  1007. term = start_term;
  1008. for(unsigned i = start + 1; i <= n; ++i)
  1009. {
  1010. term *= (n - i + 1) * x / (i * y);
  1011. result += term;
  1012. }
  1013. }
  1014. }
  1015. return result;
  1016. }
  1017. //
  1018. // The incomplete beta function implementation:
  1019. // This is just a big bunch of spaghetti code to divide up the
  1020. // input range and select the right implementation method for
  1021. // each domain:
  1022. //
  1023. template <class T, class Policy>
  1024. T ibeta_imp(T a, T b, T x, const Policy& pol, bool inv, bool normalised, T* p_derivative)
  1025. {
  1026. static const char* function = "boost::math::ibeta<%1%>(%1%, %1%, %1%)";
  1027. typedef typename lanczos::lanczos<T, Policy>::type lanczos_type;
  1028. BOOST_MATH_STD_USING // for ADL of std math functions.
  1029. BOOST_MATH_INSTRUMENT_VARIABLE(a);
  1030. BOOST_MATH_INSTRUMENT_VARIABLE(b);
  1031. BOOST_MATH_INSTRUMENT_VARIABLE(x);
  1032. BOOST_MATH_INSTRUMENT_VARIABLE(inv);
  1033. BOOST_MATH_INSTRUMENT_VARIABLE(normalised);
  1034. bool invert = inv;
  1035. T fract;
  1036. T y = 1 - x;
  1037. BOOST_MATH_ASSERT((p_derivative == 0) || normalised);
  1038. if(!(boost::math::isfinite)(a))
  1039. return policies::raise_domain_error<T>(function, "The argument a to the incomplete beta function must be finite (got a=%1%).", a, pol);
  1040. if(!(boost::math::isfinite)(b))
  1041. return policies::raise_domain_error<T>(function, "The argument b to the incomplete beta function must be finite (got b=%1%).", b, pol);
  1042. if (!(0 <= x && x <= 1))
  1043. return policies::raise_domain_error<T>(function, "The argument x to the incomplete beta function must be in [0,1] (got x=%1%).", x, pol);
  1044. if(p_derivative)
  1045. *p_derivative = -1; // value not set.
  1046. if(normalised)
  1047. {
  1048. if(a < 0)
  1049. return policies::raise_domain_error<T>(function, "The argument a to the incomplete beta function must be >= zero (got a=%1%).", a, pol);
  1050. if(b < 0)
  1051. return policies::raise_domain_error<T>(function, "The argument b to the incomplete beta function must be >= zero (got b=%1%).", b, pol);
  1052. // extend to a few very special cases:
  1053. if(a == 0)
  1054. {
  1055. if(b == 0)
  1056. return policies::raise_domain_error<T>(function, "The arguments a and b to the incomplete beta function cannot both be zero, with x=%1%.", x, pol);
  1057. if(b > 0)
  1058. return static_cast<T>(inv ? 0 : 1);
  1059. }
  1060. else if(b == 0)
  1061. {
  1062. if(a > 0)
  1063. return static_cast<T>(inv ? 1 : 0);
  1064. }
  1065. }
  1066. else
  1067. {
  1068. if(a <= 0)
  1069. return policies::raise_domain_error<T>(function, "The argument a to the incomplete beta function must be greater than zero (got a=%1%).", a, pol);
  1070. if(b <= 0)
  1071. return policies::raise_domain_error<T>(function, "The argument b to the incomplete beta function must be greater than zero (got b=%1%).", b, pol);
  1072. }
  1073. if(x == 0)
  1074. {
  1075. if(p_derivative)
  1076. {
  1077. *p_derivative = (a == 1) ? (T)1 : (a < 1) ? T(tools::max_value<T>() / 2) : T(tools::min_value<T>() * 2);
  1078. }
  1079. return (invert ? (normalised ? T(1) : boost::math::beta(a, b, pol)) : T(0));
  1080. }
  1081. if(x == 1)
  1082. {
  1083. if(p_derivative)
  1084. {
  1085. *p_derivative = (b == 1) ? T(1) : (b < 1) ? T(tools::max_value<T>() / 2) : T(tools::min_value<T>() * 2);
  1086. }
  1087. return (invert == 0 ? (normalised ? 1 : boost::math::beta(a, b, pol)) : 0);
  1088. }
  1089. if((a == 0.5f) && (b == 0.5f))
  1090. {
  1091. // We have an arcsine distribution:
  1092. if(p_derivative)
  1093. {
  1094. *p_derivative = 1 / (constants::pi<T>() * sqrt(y * x));
  1095. }
  1096. T p = invert ? asin(sqrt(y)) / constants::half_pi<T>() : asin(sqrt(x)) / constants::half_pi<T>();
  1097. if(!normalised)
  1098. p *= constants::pi<T>();
  1099. return p;
  1100. }
  1101. if(a == 1)
  1102. {
  1103. std::swap(a, b);
  1104. std::swap(x, y);
  1105. invert = !invert;
  1106. }
  1107. if(b == 1)
  1108. {
  1109. //
  1110. // Special case see: http://functions.wolfram.com/GammaBetaErf/BetaRegularized/03/01/01/
  1111. //
  1112. if(a == 1)
  1113. {
  1114. if(p_derivative)
  1115. *p_derivative = 1;
  1116. return invert ? y : x;
  1117. }
  1118. if(p_derivative)
  1119. {
  1120. *p_derivative = a * pow(x, a - 1);
  1121. }
  1122. T p; // LCOV_EXCL_LINE
  1123. if(y < 0.5)
  1124. p = invert ? T(-boost::math::expm1(a * boost::math::log1p(-y, pol), pol)) : T(exp(a * boost::math::log1p(-y, pol)));
  1125. else
  1126. p = invert ? T(-boost::math::powm1(x, a, pol)) : T(pow(x, a));
  1127. if(!normalised)
  1128. p /= a;
  1129. return p;
  1130. }
  1131. if((std::min)(a, b) <= 1)
  1132. {
  1133. if(x > 0.5)
  1134. {
  1135. std::swap(a, b);
  1136. std::swap(x, y);
  1137. invert = !invert;
  1138. BOOST_MATH_INSTRUMENT_VARIABLE(invert);
  1139. }
  1140. if((std::max)(a, b) <= 1)
  1141. {
  1142. // Both a,b < 1:
  1143. if((a >= (std::min)(T(0.2), b)) || (pow(x, a) <= 0.9))
  1144. {
  1145. if(!invert)
  1146. {
  1147. fract = ibeta_series(a, b, x, T(0), lanczos_type(), normalised, p_derivative, y, pol);
  1148. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1149. }
  1150. else
  1151. {
  1152. fract = -(normalised ? 1 : boost::math::beta(a, b, pol));
  1153. invert = false;
  1154. fract = -ibeta_series(a, b, x, fract, lanczos_type(), normalised, p_derivative, y, pol);
  1155. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1156. }
  1157. }
  1158. else
  1159. {
  1160. std::swap(a, b);
  1161. std::swap(x, y);
  1162. invert = !invert;
  1163. if(y >= 0.3)
  1164. {
  1165. if(!invert)
  1166. {
  1167. fract = ibeta_series(a, b, x, T(0), lanczos_type(), normalised, p_derivative, y, pol);
  1168. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1169. }
  1170. else
  1171. {
  1172. fract = -(normalised ? 1 : boost::math::beta(a, b, pol));
  1173. invert = false;
  1174. fract = -ibeta_series(a, b, x, fract, lanczos_type(), normalised, p_derivative, y, pol);
  1175. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1176. }
  1177. }
  1178. else
  1179. {
  1180. // Sidestep on a, and then use the series representation:
  1181. T prefix; // LCOV_EXCL_LINE
  1182. if(!normalised)
  1183. {
  1184. prefix = rising_factorial_ratio(T(a+b), a, 20);
  1185. }
  1186. else
  1187. {
  1188. prefix = 1;
  1189. }
  1190. fract = ibeta_a_step(a, b, x, y, 20, pol, normalised, p_derivative);
  1191. if(!invert)
  1192. {
  1193. fract = beta_small_b_large_a_series(T(a + 20), b, x, y, fract, prefix, pol, normalised);
  1194. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1195. }
  1196. else
  1197. {
  1198. fract -= (normalised ? 1 : boost::math::beta(a, b, pol));
  1199. invert = false;
  1200. fract = -beta_small_b_large_a_series(T(a + 20), b, x, y, fract, prefix, pol, normalised);
  1201. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1202. }
  1203. }
  1204. }
  1205. }
  1206. else
  1207. {
  1208. // One of a, b < 1 only:
  1209. if((b <= 1) || ((x < 0.1) && (pow(b * x, a) <= 0.7)))
  1210. {
  1211. if(!invert)
  1212. {
  1213. fract = ibeta_series(a, b, x, T(0), lanczos_type(), normalised, p_derivative, y, pol);
  1214. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1215. }
  1216. else
  1217. {
  1218. fract = -(normalised ? 1 : boost::math::beta(a, b, pol));
  1219. invert = false;
  1220. fract = -ibeta_series(a, b, x, fract, lanczos_type(), normalised, p_derivative, y, pol);
  1221. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1222. }
  1223. }
  1224. else
  1225. {
  1226. std::swap(a, b);
  1227. std::swap(x, y);
  1228. invert = !invert;
  1229. if(y >= 0.3)
  1230. {
  1231. if(!invert)
  1232. {
  1233. fract = ibeta_series(a, b, x, T(0), lanczos_type(), normalised, p_derivative, y, pol);
  1234. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1235. }
  1236. else
  1237. {
  1238. fract = -(normalised ? 1 : boost::math::beta(a, b, pol));
  1239. invert = false;
  1240. fract = -ibeta_series(a, b, x, fract, lanczos_type(), normalised, p_derivative, y, pol);
  1241. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1242. }
  1243. }
  1244. else if(a >= 15)
  1245. {
  1246. if(!invert)
  1247. {
  1248. fract = beta_small_b_large_a_series(a, b, x, y, T(0), T(1), pol, normalised);
  1249. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1250. }
  1251. else
  1252. {
  1253. fract = -(normalised ? 1 : boost::math::beta(a, b, pol));
  1254. invert = false;
  1255. fract = -beta_small_b_large_a_series(a, b, x, y, fract, T(1), pol, normalised);
  1256. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1257. }
  1258. }
  1259. else
  1260. {
  1261. // Sidestep to improve errors:
  1262. T prefix; // LCOV_EXCL_LINE
  1263. if(!normalised)
  1264. {
  1265. prefix = rising_factorial_ratio(T(a+b), a, 20);
  1266. }
  1267. else
  1268. {
  1269. prefix = 1;
  1270. }
  1271. fract = ibeta_a_step(a, b, x, y, 20, pol, normalised, p_derivative);
  1272. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1273. if(!invert)
  1274. {
  1275. fract = beta_small_b_large_a_series(T(a + 20), b, x, y, fract, prefix, pol, normalised);
  1276. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1277. }
  1278. else
  1279. {
  1280. fract -= (normalised ? 1 : boost::math::beta(a, b, pol));
  1281. invert = false;
  1282. fract = -beta_small_b_large_a_series(T(a + 20), b, x, y, fract, prefix, pol, normalised);
  1283. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1284. }
  1285. }
  1286. }
  1287. }
  1288. }
  1289. else
  1290. {
  1291. // Both a,b >= 1:
  1292. T lambda; // LCOV_EXCL_LINE
  1293. if(a < b)
  1294. {
  1295. lambda = a - (a + b) * x;
  1296. }
  1297. else
  1298. {
  1299. lambda = (a + b) * y - b;
  1300. }
  1301. if(lambda < 0)
  1302. {
  1303. std::swap(a, b);
  1304. std::swap(x, y);
  1305. invert = !invert;
  1306. BOOST_MATH_INSTRUMENT_VARIABLE(invert);
  1307. }
  1308. if(b < 40)
  1309. {
  1310. if((floor(a) == a) && (floor(b) == b) && (a < static_cast<T>((std::numeric_limits<int>::max)() - 100)) && (y != 1))
  1311. {
  1312. // relate to the binomial distribution and use a finite sum:
  1313. T k = a - 1;
  1314. T n = b + k;
  1315. fract = binomial_ccdf(n, k, x, y, pol);
  1316. if(!normalised)
  1317. fract *= boost::math::beta(a, b, pol);
  1318. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1319. }
  1320. else if(b * x <= 0.7)
  1321. {
  1322. if(!invert)
  1323. {
  1324. fract = ibeta_series(a, b, x, T(0), lanczos_type(), normalised, p_derivative, y, pol);
  1325. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1326. }
  1327. else
  1328. {
  1329. fract = -(normalised ? 1 : boost::math::beta(a, b, pol));
  1330. invert = false;
  1331. fract = -ibeta_series(a, b, x, fract, lanczos_type(), normalised, p_derivative, y, pol);
  1332. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1333. }
  1334. }
  1335. else if(a > 15)
  1336. {
  1337. // sidestep so we can use the series representation:
  1338. int n = itrunc(T(floor(b)), pol);
  1339. if(n == b)
  1340. --n;
  1341. T bbar = b - n;
  1342. T prefix; // LCOV_EXCL_LINE
  1343. if(!normalised)
  1344. {
  1345. prefix = rising_factorial_ratio(T(a+bbar), bbar, n);
  1346. }
  1347. else
  1348. {
  1349. prefix = 1;
  1350. }
  1351. fract = ibeta_a_step(bbar, a, y, x, n, pol, normalised, static_cast<T*>(nullptr));
  1352. fract = beta_small_b_large_a_series(a, bbar, x, y, fract, T(1), pol, normalised);
  1353. fract /= prefix;
  1354. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1355. }
  1356. else if(normalised)
  1357. {
  1358. // The formula here for the non-normalised case is tricky to figure
  1359. // out (for me!!), and requires two pochhammer calculations rather
  1360. // than one, so leave it for now and only use this in the normalized case....
  1361. int n = itrunc(T(floor(b)), pol);
  1362. T bbar = b - n;
  1363. if(bbar <= 0)
  1364. {
  1365. --n;
  1366. bbar += 1;
  1367. }
  1368. fract = ibeta_a_step(bbar, a, y, x, n, pol, normalised, static_cast<T*>(nullptr));
  1369. fract += ibeta_a_step(a, bbar, x, y, 20, pol, normalised, static_cast<T*>(nullptr));
  1370. if(invert)
  1371. fract -= 1; // Note this line would need changing if we ever enable this branch in non-normalized case
  1372. fract = beta_small_b_large_a_series(T(a+20), bbar, x, y, fract, T(1), pol, normalised);
  1373. if(invert)
  1374. {
  1375. fract = -fract;
  1376. invert = false;
  1377. }
  1378. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1379. }
  1380. else
  1381. {
  1382. fract = ibeta_fraction2(a, b, x, y, pol, normalised, p_derivative);
  1383. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1384. }
  1385. }
  1386. else
  1387. {
  1388. fract = ibeta_fraction2(a, b, x, y, pol, normalised, p_derivative);
  1389. BOOST_MATH_INSTRUMENT_VARIABLE(fract);
  1390. }
  1391. }
  1392. if(p_derivative)
  1393. {
  1394. if(*p_derivative < 0)
  1395. {
  1396. *p_derivative = ibeta_power_terms(a, b, x, y, lanczos_type(), true, pol);
  1397. }
  1398. T div = y * x;
  1399. if(*p_derivative != 0)
  1400. {
  1401. if((tools::max_value<T>() * div < *p_derivative))
  1402. {
  1403. // overflow, return an arbitrarily large value:
  1404. *p_derivative = tools::max_value<T>() / 2; // LCOV_EXCL_LINE Probably can only get here with denormalized x.
  1405. }
  1406. else
  1407. {
  1408. *p_derivative /= div;
  1409. }
  1410. }
  1411. }
  1412. return invert ? (normalised ? 1 : boost::math::beta(a, b, pol)) - fract : fract;
  1413. } // template <class T, class Lanczos>T ibeta_imp(T a, T b, T x, const Lanczos& l, bool inv, bool normalised)
  1414. template <class T, class Policy>
  1415. inline T ibeta_imp(T a, T b, T x, const Policy& pol, bool inv, bool normalised)
  1416. {
  1417. return ibeta_imp(a, b, x, pol, inv, normalised, static_cast<T*>(nullptr));
  1418. }
  1419. template <class T, class Policy>
  1420. T ibeta_derivative_imp(T a, T b, T x, const Policy& pol)
  1421. {
  1422. static const char* function = "ibeta_derivative<%1%>(%1%,%1%,%1%)";
  1423. //
  1424. // start with the usual error checks:
  1425. //
  1426. if (!(boost::math::isfinite)(a))
  1427. return policies::raise_domain_error<T>(function, "The argument a to the incomplete beta function must be finite (got a=%1%).", a, pol);
  1428. if (!(boost::math::isfinite)(b))
  1429. return policies::raise_domain_error<T>(function, "The argument b to the incomplete beta function must be finite (got b=%1%).", b, pol);
  1430. if (!(0 <= x && x <= 1))
  1431. return policies::raise_domain_error<T>(function, "The argument x to the incomplete beta function must be in [0,1] (got x=%1%).", x, pol);
  1432. if(a <= 0)
  1433. return policies::raise_domain_error<T>(function, "The argument a to the incomplete beta function must be greater than zero (got a=%1%).", a, pol);
  1434. if(b <= 0)
  1435. return policies::raise_domain_error<T>(function, "The argument b to the incomplete beta function must be greater than zero (got b=%1%).", b, pol);
  1436. //
  1437. // Now the corner cases:
  1438. //
  1439. if(x == 0)
  1440. {
  1441. return (a > 1) ? 0 :
  1442. (a == 1) ? 1 / boost::math::beta(a, b, pol) : policies::raise_overflow_error<T>(function, nullptr, pol);
  1443. }
  1444. else if(x == 1)
  1445. {
  1446. return (b > 1) ? 0 :
  1447. (b == 1) ? 1 / boost::math::beta(a, b, pol) : policies::raise_overflow_error<T>(function, nullptr, pol);
  1448. }
  1449. //
  1450. // Now the regular cases:
  1451. //
  1452. typedef typename lanczos::lanczos<T, Policy>::type lanczos_type;
  1453. T y = (1 - x) * x;
  1454. T f1;
  1455. if (!(boost::math::isinf)(1 / y))
  1456. {
  1457. f1 = ibeta_power_terms<T>(a, b, x, 1 - x, lanczos_type(), true, pol, 1 / y, function);
  1458. }
  1459. else
  1460. {
  1461. return (a > 1) ? 0 : (a == 1) ? 1 / boost::math::beta(a, b, pol) : policies::raise_overflow_error<T>(function, nullptr, pol);
  1462. }
  1463. return f1;
  1464. }
  1465. //
  1466. // Some forwarding functions that disambiguate the third argument type:
  1467. //
  1468. template <class RT1, class RT2, class Policy>
  1469. inline typename tools::promote_args<RT1, RT2>::type
  1470. beta(RT1 a, RT2 b, const Policy&, const std::true_type*)
  1471. {
  1472. BOOST_FPU_EXCEPTION_GUARD
  1473. typedef typename tools::promote_args<RT1, RT2>::type result_type;
  1474. typedef typename policies::evaluation<result_type, Policy>::type value_type;
  1475. typedef typename lanczos::lanczos<value_type, Policy>::type evaluation_type;
  1476. typedef typename policies::normalise<
  1477. Policy,
  1478. policies::promote_float<false>,
  1479. policies::promote_double<false>,
  1480. policies::discrete_quantile<>,
  1481. policies::assert_undefined<> >::type forwarding_policy;
  1482. return policies::checked_narrowing_cast<result_type, forwarding_policy>(detail::beta_imp(static_cast<value_type>(a), static_cast<value_type>(b), evaluation_type(), forwarding_policy()), "boost::math::beta<%1%>(%1%,%1%)");
  1483. }
  1484. template <class RT1, class RT2, class RT3>
  1485. inline typename tools::promote_args<RT1, RT2, RT3>::type
  1486. beta(RT1 a, RT2 b, RT3 x, const std::false_type*)
  1487. {
  1488. return boost::math::beta(a, b, x, policies::policy<>());
  1489. }
  1490. } // namespace detail
  1491. //
  1492. // The actual function entry-points now follow, these just figure out
  1493. // which Lanczos approximation to use
  1494. // and forward to the implementation functions:
  1495. //
  1496. template <class RT1, class RT2, class A>
  1497. inline typename tools::promote_args<RT1, RT2, A>::type
  1498. beta(RT1 a, RT2 b, A arg)
  1499. {
  1500. using tag = typename policies::is_policy<A>::type;
  1501. using ReturnType = tools::promote_args_t<RT1, RT2, A>;
  1502. return static_cast<ReturnType>(boost::math::detail::beta(a, b, arg, static_cast<tag*>(nullptr)));
  1503. }
  1504. template <class RT1, class RT2>
  1505. inline typename tools::promote_args<RT1, RT2>::type
  1506. beta(RT1 a, RT2 b)
  1507. {
  1508. return boost::math::beta(a, b, policies::policy<>());
  1509. }
  1510. template <class RT1, class RT2, class RT3, class Policy>
  1511. inline typename tools::promote_args<RT1, RT2, RT3>::type
  1512. beta(RT1 a, RT2 b, RT3 x, const Policy&)
  1513. {
  1514. BOOST_FPU_EXCEPTION_GUARD
  1515. typedef typename tools::promote_args<RT1, RT2, RT3>::type result_type;
  1516. typedef typename policies::evaluation<result_type, Policy>::type value_type;
  1517. typedef typename policies::normalise<
  1518. Policy,
  1519. policies::promote_float<false>,
  1520. policies::promote_double<false>,
  1521. policies::discrete_quantile<>,
  1522. policies::assert_undefined<> >::type forwarding_policy;
  1523. return policies::checked_narrowing_cast<result_type, forwarding_policy>(detail::ibeta_imp(static_cast<value_type>(a), static_cast<value_type>(b), static_cast<value_type>(x), forwarding_policy(), false, false), "boost::math::beta<%1%>(%1%,%1%,%1%)");
  1524. }
  1525. template <class RT1, class RT2, class RT3, class Policy>
  1526. inline typename tools::promote_args<RT1, RT2, RT3>::type
  1527. betac(RT1 a, RT2 b, RT3 x, const Policy&)
  1528. {
  1529. BOOST_FPU_EXCEPTION_GUARD
  1530. typedef typename tools::promote_args<RT1, RT2, RT3>::type result_type;
  1531. typedef typename policies::evaluation<result_type, Policy>::type value_type;
  1532. typedef typename policies::normalise<
  1533. Policy,
  1534. policies::promote_float<false>,
  1535. policies::promote_double<false>,
  1536. policies::discrete_quantile<>,
  1537. policies::assert_undefined<> >::type forwarding_policy;
  1538. return policies::checked_narrowing_cast<result_type, forwarding_policy>(detail::ibeta_imp(static_cast<value_type>(a), static_cast<value_type>(b), static_cast<value_type>(x), forwarding_policy(), true, false), "boost::math::betac<%1%>(%1%,%1%,%1%)");
  1539. }
  1540. template <class RT1, class RT2, class RT3>
  1541. inline typename tools::promote_args<RT1, RT2, RT3>::type
  1542. betac(RT1 a, RT2 b, RT3 x)
  1543. {
  1544. return boost::math::betac(a, b, x, policies::policy<>());
  1545. }
  1546. template <class RT1, class RT2, class RT3, class Policy>
  1547. inline typename tools::promote_args<RT1, RT2, RT3>::type
  1548. ibeta(RT1 a, RT2 b, RT3 x, const Policy&)
  1549. {
  1550. BOOST_FPU_EXCEPTION_GUARD
  1551. typedef typename tools::promote_args<RT1, RT2, RT3>::type result_type;
  1552. typedef typename policies::evaluation<result_type, Policy>::type value_type;
  1553. typedef typename policies::normalise<
  1554. Policy,
  1555. policies::promote_float<false>,
  1556. policies::promote_double<false>,
  1557. policies::discrete_quantile<>,
  1558. policies::assert_undefined<> >::type forwarding_policy;
  1559. return policies::checked_narrowing_cast<result_type, forwarding_policy>(detail::ibeta_imp(static_cast<value_type>(a), static_cast<value_type>(b), static_cast<value_type>(x), forwarding_policy(), false, true), "boost::math::ibeta<%1%>(%1%,%1%,%1%)");
  1560. }
  1561. template <class RT1, class RT2, class RT3>
  1562. inline typename tools::promote_args<RT1, RT2, RT3>::type
  1563. ibeta(RT1 a, RT2 b, RT3 x)
  1564. {
  1565. return boost::math::ibeta(a, b, x, policies::policy<>());
  1566. }
  1567. template <class RT1, class RT2, class RT3, class Policy>
  1568. inline typename tools::promote_args<RT1, RT2, RT3>::type
  1569. ibetac(RT1 a, RT2 b, RT3 x, const Policy&)
  1570. {
  1571. BOOST_FPU_EXCEPTION_GUARD
  1572. typedef typename tools::promote_args<RT1, RT2, RT3>::type result_type;
  1573. typedef typename policies::evaluation<result_type, Policy>::type value_type;
  1574. typedef typename policies::normalise<
  1575. Policy,
  1576. policies::promote_float<false>,
  1577. policies::promote_double<false>,
  1578. policies::discrete_quantile<>,
  1579. policies::assert_undefined<> >::type forwarding_policy;
  1580. return policies::checked_narrowing_cast<result_type, forwarding_policy>(detail::ibeta_imp(static_cast<value_type>(a), static_cast<value_type>(b), static_cast<value_type>(x), forwarding_policy(), true, true), "boost::math::ibetac<%1%>(%1%,%1%,%1%)");
  1581. }
  1582. template <class RT1, class RT2, class RT3>
  1583. inline typename tools::promote_args<RT1, RT2, RT3>::type
  1584. ibetac(RT1 a, RT2 b, RT3 x)
  1585. {
  1586. return boost::math::ibetac(a, b, x, policies::policy<>());
  1587. }
  1588. template <class RT1, class RT2, class RT3, class Policy>
  1589. inline typename tools::promote_args<RT1, RT2, RT3>::type
  1590. ibeta_derivative(RT1 a, RT2 b, RT3 x, const Policy&)
  1591. {
  1592. BOOST_FPU_EXCEPTION_GUARD
  1593. typedef typename tools::promote_args<RT1, RT2, RT3>::type result_type;
  1594. typedef typename policies::evaluation<result_type, Policy>::type value_type;
  1595. typedef typename policies::normalise<
  1596. Policy,
  1597. policies::promote_float<false>,
  1598. policies::promote_double<false>,
  1599. policies::discrete_quantile<>,
  1600. policies::assert_undefined<> >::type forwarding_policy;
  1601. return policies::checked_narrowing_cast<result_type, forwarding_policy>(detail::ibeta_derivative_imp(static_cast<value_type>(a), static_cast<value_type>(b), static_cast<value_type>(x), forwarding_policy()), "boost::math::ibeta_derivative<%1%>(%1%,%1%,%1%)");
  1602. }
  1603. template <class RT1, class RT2, class RT3>
  1604. inline typename tools::promote_args<RT1, RT2, RT3>::type
  1605. ibeta_derivative(RT1 a, RT2 b, RT3 x)
  1606. {
  1607. return boost::math::ibeta_derivative(a, b, x, policies::policy<>());
  1608. }
  1609. } // namespace math
  1610. } // namespace boost
  1611. #include <boost/math/special_functions/detail/ibeta_inverse.hpp>
  1612. #include <boost/math/special_functions/detail/ibeta_inv_ab.hpp>
  1613. #endif // BOOST_MATH_SPECIAL_BETA_HPP