beta.hpp 63 KB

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