non_central_chi_squared.hpp 42 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999
  1. // boost\math\distributions\non_central_chi_squared.hpp
  2. // Copyright John Maddock 2008.
  3. // Copyright Matt Borland 2024.
  4. // Use, modification and distribution are subject to the
  5. // Boost Software License, Version 1.0.
  6. // (See accompanying file LICENSE_1_0.txt
  7. // or copy at http://www.boost.org/LICENSE_1_0.txt)
  8. #ifndef BOOST_MATH_SPECIAL_NON_CENTRAL_CHI_SQUARE_HPP
  9. #define BOOST_MATH_SPECIAL_NON_CENTRAL_CHI_SQUARE_HPP
  10. #include <boost/math/tools/config.hpp>
  11. #include <boost/math/tools/tuple.hpp>
  12. #include <boost/math/tools/cstdint.hpp>
  13. #include <boost/math/tools/numeric_limits.hpp>
  14. #include <boost/math/distributions/fwd.hpp>
  15. #include <boost/math/special_functions/gamma.hpp> // for incomplete gamma. gamma_q
  16. #include <boost/math/special_functions/bessel.hpp> // for cyl_bessel_i
  17. #include <boost/math/special_functions/round.hpp> // for llround
  18. #include <boost/math/distributions/complement.hpp> // complements
  19. #include <boost/math/distributions/chi_squared.hpp> // central distribution
  20. #include <boost/math/distributions/detail/common_error_handling.hpp> // error checks
  21. #include <boost/math/special_functions/fpclassify.hpp> // isnan.
  22. #include <boost/math/tools/roots.hpp> // for root finding.
  23. #include <boost/math/distributions/detail/generic_mode.hpp>
  24. #include <boost/math/distributions/detail/generic_quantile.hpp>
  25. #include <boost/math/policies/policy.hpp>
  26. namespace boost
  27. {
  28. namespace math
  29. {
  30. template <class RealType, class Policy>
  31. class non_central_chi_squared_distribution;
  32. namespace detail{
  33. template <class T, class Policy>
  34. BOOST_MATH_GPU_ENABLED T non_central_chi_square_q(T x, T f, T theta, const Policy& pol, T init_sum = 0)
  35. {
  36. //
  37. // Computes the complement of the Non-Central Chi-Square
  38. // Distribution CDF by summing a weighted sum of complements
  39. // of the central-distributions. The weighting factor is
  40. // a Poisson Distribution.
  41. //
  42. // This is an application of the technique described in:
  43. //
  44. // Computing discrete mixtures of continuous
  45. // distributions: noncentral chisquare, noncentral t
  46. // and the distribution of the square of the sample
  47. // multiple correlation coefficient.
  48. // D. Benton, K. Krishnamoorthy.
  49. // Computational Statistics & Data Analysis 43 (2003) 249 - 267
  50. //
  51. BOOST_MATH_STD_USING
  52. // Special case:
  53. if(x == 0)
  54. return 1;
  55. //
  56. // Initialize the variables we'll be using:
  57. //
  58. T lambda = theta / 2;
  59. T del = f / 2;
  60. T y = x / 2;
  61. boost::math::uintmax_t max_iter = policies::get_max_series_iterations<Policy>();
  62. T errtol = boost::math::policies::get_epsilon<T, Policy>();
  63. T sum = init_sum;
  64. //
  65. // k is the starting location for iteration, we'll
  66. // move both forwards and backwards from this point.
  67. // k is chosen as the peek of the Poisson weights, which
  68. // will occur *before* the largest term.
  69. //
  70. long long k = llround(lambda, pol);
  71. // Forwards and backwards Poisson weights:
  72. T poisf = boost::math::gamma_p_derivative(static_cast<T>(1 + k), lambda, pol);
  73. T poisb = poisf * k / lambda;
  74. // Initial forwards central chi squared term:
  75. T gamf = boost::math::gamma_q(del + k, y, pol);
  76. // Forwards and backwards recursion terms on the central chi squared:
  77. T xtermf = boost::math::gamma_p_derivative(del + 1 + k, y, pol);
  78. T xtermb = xtermf * (del + k) / y;
  79. // Initial backwards central chi squared term:
  80. T gamb = gamf - xtermb;
  81. //
  82. // Forwards iteration first, this is the
  83. // stable direction for the gamma function
  84. // recurrences:
  85. //
  86. long long i;
  87. for(i = k; static_cast<boost::math::uintmax_t>(i-k) < max_iter; ++i)
  88. {
  89. T term = poisf * gamf;
  90. sum += term;
  91. poisf *= lambda / (i + 1);
  92. gamf += xtermf;
  93. xtermf *= y / (del + i + 1);
  94. if(((sum == 0) || (fabs(term / sum) < errtol)) && (term >= poisf * gamf))
  95. break;
  96. }
  97. //Error check:
  98. if(static_cast<boost::math::uintmax_t>(i-k) >= max_iter)
  99. return policies::raise_evaluation_error("cdf(non_central_chi_squared_distribution<%1%>, %1%)", "Series did not converge, closest value was %1%", sum, pol); // LCOV_EXCL_LINE
  100. //
  101. // Now backwards iteration: the gamma
  102. // function recurrences are unstable in this
  103. // direction, we rely on the terms diminishing in size
  104. // faster than we introduce cancellation errors.
  105. // For this reason it's very important that we start
  106. // *before* the largest term so that backwards iteration
  107. // is strictly converging.
  108. //
  109. for(i = k - 1; i >= 0; --i)
  110. {
  111. T term = poisb * gamb;
  112. sum += term;
  113. poisb *= i / lambda;
  114. xtermb *= (del + i) / y;
  115. gamb -= xtermb;
  116. if((sum == 0) || (fabs(term / sum) < errtol))
  117. break;
  118. }
  119. return sum;
  120. }
  121. template <class T, class Policy>
  122. BOOST_MATH_GPU_ENABLED T non_central_chi_square_p_ding(T x, T f, T theta, const Policy& pol, T init_sum = 0)
  123. {
  124. //
  125. // This is an implementation of:
  126. //
  127. // Algorithm AS 275:
  128. // Computing the Non-Central #2 Distribution Function
  129. // Cherng G. Ding
  130. // Applied Statistics, Vol. 41, No. 2. (1992), pp. 478-482.
  131. //
  132. // This uses a stable forward iteration to sum the
  133. // CDF, unfortunately this can not be used for large
  134. // values of the non-centrality parameter because:
  135. // * The first term may underflow to zero.
  136. // * We may need an extra-ordinary number of terms
  137. // before we reach the first *significant* term.
  138. //
  139. BOOST_MATH_STD_USING
  140. // Special case:
  141. if(x == 0)
  142. return 0;
  143. T tk = boost::math::gamma_p_derivative(f/2 + 1, x/2, pol);
  144. T lambda = theta / 2;
  145. T vk = exp(-lambda);
  146. T uk = vk;
  147. T sum = init_sum + tk * vk;
  148. if(sum == 0)
  149. return sum;
  150. boost::math::uintmax_t max_iter = policies::get_max_series_iterations<Policy>();
  151. T errtol = boost::math::policies::get_epsilon<T, Policy>();
  152. int i;
  153. T lterm(0), term(0);
  154. for(i = 1; static_cast<boost::math::uintmax_t>(i) < max_iter; ++i)
  155. {
  156. tk = tk * x / (f + 2 * i);
  157. uk = uk * lambda / i;
  158. vk = vk + uk;
  159. lterm = term;
  160. term = vk * tk;
  161. sum += term;
  162. if((fabs(term / sum) < errtol) && (term <= lterm))
  163. break;
  164. }
  165. //Error check:
  166. if(static_cast<boost::math::uintmax_t>(i) >= max_iter)
  167. return policies::raise_evaluation_error("cdf(non_central_chi_squared_distribution<%1%>, %1%)", "Series did not converge, closest value was %1%", sum, pol); // LCOV_EXCL_LINE
  168. return sum;
  169. }
  170. template <class T, class Policy>
  171. BOOST_MATH_GPU_ENABLED T non_central_chi_square_p(T y, T n, T lambda, const Policy& pol, T init_sum)
  172. {
  173. //
  174. // This is taken more or less directly from:
  175. //
  176. // Computing discrete mixtures of continuous
  177. // distributions: noncentral chisquare, noncentral t
  178. // and the distribution of the square of the sample
  179. // multiple correlation coefficient.
  180. // D. Benton, K. Krishnamoorthy.
  181. // Computational Statistics & Data Analysis 43 (2003) 249 - 267
  182. //
  183. // We're summing a Poisson weighting term multiplied by
  184. // a central chi squared distribution.
  185. //
  186. BOOST_MATH_STD_USING
  187. // Special case:
  188. if(y == 0)
  189. return 0;
  190. boost::math::uintmax_t max_iter = policies::get_max_series_iterations<Policy>();
  191. T errtol = boost::math::policies::get_epsilon<T, Policy>();
  192. T errorf(0), errorb(0);
  193. T x = y / 2;
  194. T del = lambda / 2;
  195. //
  196. // Starting location for the iteration, we'll iterate
  197. // both forwards and backwards from this point. The
  198. // location chosen is the maximum of the Poisson weight
  199. // function, which ocurrs *after* the largest term in the
  200. // sum.
  201. //
  202. long long k = llround(del, pol);
  203. T a = n / 2 + k;
  204. // Central chi squared term for forward iteration:
  205. T gamkf = boost::math::gamma_p(a, x, pol);
  206. if(lambda == 0)
  207. return gamkf;
  208. // Central chi squared term for backward iteration:
  209. T gamkb = gamkf;
  210. // Forwards Poisson weight:
  211. T poiskf = gamma_p_derivative(static_cast<T>(k+1), del, pol);
  212. // Backwards Poisson weight:
  213. T poiskb = poiskf;
  214. // Forwards gamma function recursion term:
  215. T xtermf = boost::math::gamma_p_derivative(a, x, pol);
  216. // Backwards gamma function recursion term:
  217. T xtermb = xtermf * x / a;
  218. T sum = init_sum + poiskf * gamkf;
  219. if(sum == 0)
  220. return sum;
  221. int i = 1;
  222. //
  223. // Backwards recursion first, this is the stable
  224. // direction for gamma function recurrences:
  225. //
  226. while(i <= k)
  227. {
  228. xtermb *= (a - i + 1) / x;
  229. gamkb += xtermb;
  230. poiskb = poiskb * (k - i + 1) / del;
  231. errorf = errorb;
  232. errorb = gamkb * poiskb;
  233. sum += errorb;
  234. if((fabs(errorb / sum) < errtol) && (errorb <= errorf))
  235. break;
  236. ++i;
  237. }
  238. i = 1;
  239. //
  240. // Now forwards recursion, the gamma function
  241. // recurrence relation is unstable in this direction,
  242. // so we rely on the magnitude of successive terms
  243. // decreasing faster than we introduce cancellation error.
  244. // For this reason it's vital that k is chosen to be *after*
  245. // the largest term, so that successive forward iterations
  246. // are strictly (and rapidly) converging.
  247. //
  248. do
  249. {
  250. xtermf = xtermf * x / (a + i - 1);
  251. gamkf = gamkf - xtermf;
  252. poiskf = poiskf * del / (k + i);
  253. errorf = poiskf * gamkf;
  254. sum += errorf;
  255. ++i;
  256. }while((fabs(errorf / sum) > errtol) && (static_cast<boost::math::uintmax_t>(i) < max_iter));
  257. //Error check:
  258. if(static_cast<boost::math::uintmax_t>(i) >= max_iter)
  259. return policies::raise_evaluation_error("cdf(non_central_chi_squared_distribution<%1%>, %1%)", "Series did not converge, closest value was %1%", sum, pol); // LCOV_EXCL_LINE
  260. return sum;
  261. }
  262. template <class T, class Policy>
  263. BOOST_MATH_GPU_ENABLED T non_central_chi_square_pdf(T x, T n, T lambda, const Policy& pol)
  264. {
  265. //
  266. // As above but for the PDF:
  267. //
  268. BOOST_MATH_STD_USING
  269. boost::math::uintmax_t max_iter = policies::get_max_series_iterations<Policy>();
  270. T errtol = boost::math::policies::get_epsilon<T, Policy>();
  271. T x2 = x / 2;
  272. T n2 = n / 2;
  273. T l2 = lambda / 2;
  274. T sum = 0;
  275. long long k = lltrunc(l2);
  276. T pois = gamma_p_derivative(static_cast<T>(k + 1), l2, pol) * gamma_p_derivative(static_cast<T>(n2 + k), x2);
  277. if(pois == 0)
  278. return 0;
  279. T poisb = pois;
  280. for(long long i = k; ; ++i)
  281. {
  282. sum += pois;
  283. if(pois / sum < errtol)
  284. break;
  285. if(static_cast<boost::math::uintmax_t>(i - k) >= max_iter)
  286. return policies::raise_evaluation_error("pdf(non_central_chi_squared_distribution<%1%>, %1%)", "Series did not converge, closest value was %1%", sum, pol); // LCOV_EXCL_LINE
  287. pois *= l2 * x2 / ((i + 1) * (n2 + i));
  288. }
  289. for(long long i = k - 1; i >= 0; --i)
  290. {
  291. poisb *= (i + 1) * (n2 + i) / (l2 * x2);
  292. sum += poisb;
  293. if(poisb / sum < errtol)
  294. break;
  295. }
  296. return sum / 2;
  297. }
  298. template <class RealType, class Policy>
  299. BOOST_MATH_GPU_ENABLED inline RealType non_central_chi_squared_cdf(RealType x, RealType k, RealType l, bool invert, const Policy&)
  300. {
  301. typedef typename policies::evaluation<RealType, Policy>::type value_type;
  302. typedef typename policies::normalise<
  303. Policy,
  304. policies::promote_float<false>,
  305. policies::promote_double<false>,
  306. policies::discrete_quantile<>,
  307. policies::assert_undefined<> >::type forwarding_policy;
  308. BOOST_MATH_STD_USING
  309. value_type result;
  310. if(l == 0)
  311. return invert == false ? cdf(boost::math::chi_squared_distribution<RealType, Policy>(k), x) : cdf(complement(boost::math::chi_squared_distribution<RealType, Policy>(k), x));
  312. else if(x > k + l)
  313. {
  314. // Complement is the smaller of the two:
  315. result = detail::non_central_chi_square_q(
  316. static_cast<value_type>(x),
  317. static_cast<value_type>(k),
  318. static_cast<value_type>(l),
  319. forwarding_policy(),
  320. static_cast<value_type>(invert ? 0 : -1));
  321. invert = !invert;
  322. }
  323. else if(l < 200)
  324. {
  325. // For small values of the non-centrality parameter
  326. // we can use Ding's method:
  327. result = detail::non_central_chi_square_p_ding(
  328. static_cast<value_type>(x),
  329. static_cast<value_type>(k),
  330. static_cast<value_type>(l),
  331. forwarding_policy(),
  332. static_cast<value_type>(invert ? -1 : 0));
  333. }
  334. else
  335. {
  336. // For largers values of the non-centrality
  337. // parameter Ding's method will consume an
  338. // extra-ordinary number of terms, and worse
  339. // may return zero when the result is in fact
  340. // finite, use Krishnamoorthy's method instead:
  341. result = detail::non_central_chi_square_p(
  342. static_cast<value_type>(x),
  343. static_cast<value_type>(k),
  344. static_cast<value_type>(l),
  345. forwarding_policy(),
  346. static_cast<value_type>(invert ? -1 : 0));
  347. }
  348. if(invert)
  349. result = -result;
  350. return policies::checked_narrowing_cast<RealType, forwarding_policy>(
  351. result,
  352. "boost::math::non_central_chi_squared_cdf<%1%>(%1%, %1%, %1%)");
  353. }
  354. template <class T, class Policy>
  355. struct nccs_quantile_functor
  356. {
  357. BOOST_MATH_GPU_ENABLED nccs_quantile_functor(const non_central_chi_squared_distribution<T,Policy>& d, T t, bool c)
  358. : dist(d), target(t), comp(c) {}
  359. BOOST_MATH_GPU_ENABLED T operator()(const T& x)
  360. {
  361. return comp ?
  362. target - cdf(complement(dist, x))
  363. : cdf(dist, x) - target;
  364. }
  365. private:
  366. non_central_chi_squared_distribution<T,Policy> dist;
  367. T target;
  368. bool comp;
  369. };
  370. template <class RealType, class Policy>
  371. BOOST_MATH_GPU_ENABLED RealType nccs_quantile(const non_central_chi_squared_distribution<RealType, Policy>& dist, const RealType& p, bool comp)
  372. {
  373. BOOST_MATH_STD_USING
  374. constexpr auto function = "quantile(non_central_chi_squared_distribution<%1%>, %1%)";
  375. typedef typename policies::evaluation<RealType, Policy>::type value_type;
  376. typedef typename policies::normalise<
  377. Policy,
  378. policies::promote_float<false>,
  379. policies::promote_double<false>,
  380. policies::discrete_quantile<>,
  381. policies::assert_undefined<> >::type forwarding_policy;
  382. value_type k = dist.degrees_of_freedom();
  383. value_type l = dist.non_centrality();
  384. value_type r;
  385. if(!detail::check_df(
  386. function,
  387. k, &r, Policy())
  388. ||
  389. !detail::check_non_centrality(
  390. function,
  391. l,
  392. &r,
  393. Policy())
  394. ||
  395. !detail::check_probability(
  396. function,
  397. static_cast<value_type>(p),
  398. &r,
  399. Policy()))
  400. return static_cast<RealType>(r);
  401. //
  402. // Special cases get short-circuited first:
  403. //
  404. if(p == 0)
  405. return comp ? policies::raise_overflow_error<RealType>(function, 0, Policy()) : 0;
  406. if(p == 1)
  407. return comp ? 0 : policies::raise_overflow_error<RealType>(function, 0, Policy());
  408. //
  409. // This is Pearson's approximation to the quantile, see
  410. // Pearson, E. S. (1959) "Note on an approximation to the distribution of
  411. // noncentral chi squared", Biometrika 46: 364.
  412. // See also:
  413. // "A comparison of approximations to percentiles of the noncentral chi2-distribution",
  414. // Hardeo Sahai and Mario Miguel Ojeda, Revista de Matematica: Teoria y Aplicaciones 2003 10(1-2) : 57-76.
  415. // Note that the latter reference refers to an approximation of the CDF, when they really mean the quantile.
  416. //
  417. value_type b = -(l * l) / (k + 3 * l);
  418. value_type c = (k + 3 * l) / (k + 2 * l);
  419. value_type ff = (k + 2 * l) / (c * c);
  420. value_type guess;
  421. if(comp)
  422. {
  423. guess = b + c * quantile(complement(chi_squared_distribution<value_type, forwarding_policy>(ff), p));
  424. }
  425. else
  426. {
  427. guess = b + c * quantile(chi_squared_distribution<value_type, forwarding_policy>(ff), p);
  428. }
  429. //
  430. // Sometimes guess goes very small or negative, in that case we have
  431. // to do something else for the initial guess, this approximation
  432. // was provided in a private communication from Thomas Luu, PhD candidate,
  433. // University College London. It's an asymptotic expansion for the
  434. // quantile which usually gets us within an order of magnitude of the
  435. // correct answer.
  436. // Fast and accurate parallel computation of quantile functions for random number generation,
  437. // Thomas LuuDoctorial Thesis 2016
  438. // http://discovery.ucl.ac.uk/1482128/
  439. //
  440. if(guess < 0.005)
  441. {
  442. value_type pp = comp ? 1 - p : p;
  443. //guess = pow(pow(value_type(2), (k / 2 - 1)) * exp(l / 2) * pp * k, 2 / k);
  444. guess = pow(pow(value_type(2), (k / 2 - 1)) * exp(l / 2) * pp * k * boost::math::tgamma(k / 2, forwarding_policy()), (2 / k));
  445. if(guess == 0)
  446. guess = tools::min_value<value_type>();
  447. }
  448. value_type result = detail::generic_quantile(
  449. non_central_chi_squared_distribution<value_type, forwarding_policy>(k, l),
  450. p,
  451. guess,
  452. comp,
  453. function);
  454. return policies::checked_narrowing_cast<RealType, forwarding_policy>(
  455. result,
  456. function);
  457. }
  458. template <class RealType, class Policy>
  459. BOOST_MATH_GPU_ENABLED RealType nccs_pdf(const non_central_chi_squared_distribution<RealType, Policy>& dist, const RealType& x)
  460. {
  461. BOOST_MATH_STD_USING
  462. constexpr auto function = "pdf(non_central_chi_squared_distribution<%1%>, %1%)";
  463. typedef typename policies::evaluation<RealType, Policy>::type value_type;
  464. typedef typename policies::normalise<
  465. Policy,
  466. policies::promote_float<false>,
  467. policies::promote_double<false>,
  468. policies::discrete_quantile<>,
  469. policies::assert_undefined<> >::type forwarding_policy;
  470. value_type k = dist.degrees_of_freedom();
  471. value_type l = dist.non_centrality();
  472. value_type r;
  473. if(!detail::check_df(
  474. function,
  475. k, &r, Policy())
  476. ||
  477. !detail::check_non_centrality(
  478. function,
  479. l,
  480. &r,
  481. Policy())
  482. ||
  483. !detail::check_positive_x(
  484. function,
  485. (value_type)x,
  486. &r,
  487. Policy()))
  488. return static_cast<RealType>(r);
  489. if(l == 0)
  490. return pdf(boost::math::chi_squared_distribution<RealType, forwarding_policy>(dist.degrees_of_freedom()), x);
  491. // Special case:
  492. if(x == 0)
  493. return 0;
  494. if(l > 50)
  495. {
  496. r = non_central_chi_square_pdf(static_cast<value_type>(x), k, l, forwarding_policy());
  497. }
  498. else
  499. {
  500. r = log(x / l) * (k / 4 - 0.5f) - (x + l) / 2;
  501. if(fabs(r) >= tools::log_max_value<RealType>() / 4)
  502. {
  503. r = non_central_chi_square_pdf(static_cast<value_type>(x), k, l, forwarding_policy());
  504. }
  505. else
  506. {
  507. r = exp(r);
  508. r = 0.5f * r
  509. * boost::math::cyl_bessel_i(k/2 - 1, sqrt(l * x), forwarding_policy());
  510. }
  511. }
  512. return policies::checked_narrowing_cast<RealType, forwarding_policy>(
  513. r,
  514. function);
  515. }
  516. template <class RealType, class Policy>
  517. struct degrees_of_freedom_finder
  518. {
  519. BOOST_MATH_GPU_ENABLED degrees_of_freedom_finder(
  520. RealType lam_, RealType x_, RealType p_, bool c)
  521. : lam(lam_), x(x_), p(p_), comp(c) {}
  522. BOOST_MATH_GPU_ENABLED RealType operator()(const RealType& v)
  523. {
  524. non_central_chi_squared_distribution<RealType, Policy> d(v, lam);
  525. return comp ?
  526. RealType(p - cdf(complement(d, x)))
  527. : RealType(cdf(d, x) - p);
  528. }
  529. private:
  530. RealType lam;
  531. RealType x;
  532. RealType p;
  533. bool comp;
  534. };
  535. template <class RealType, class Policy>
  536. BOOST_MATH_GPU_ENABLED inline RealType find_degrees_of_freedom(
  537. RealType lam, RealType x, RealType p, RealType q, const Policy& pol)
  538. {
  539. constexpr auto function = "non_central_chi_squared<%1%>::find_degrees_of_freedom";
  540. if((p == 0) || (q == 0))
  541. {
  542. //
  543. // Can't a thing if one of p and q is zero:
  544. //
  545. return policies::raise_evaluation_error<RealType>(function, "Can't find degrees of freedom when the probability is 0 or 1, only possible answer is %1%", // LCOV_EXCL_LINE
  546. RealType(boost::math::numeric_limits<RealType>::quiet_NaN()), Policy()); // LCOV_EXCL_LINE
  547. }
  548. degrees_of_freedom_finder<RealType, Policy> f(lam, x, p < q ? p : q, p < q ? false : true);
  549. tools::eps_tolerance<RealType> tol(policies::digits<RealType, Policy>());
  550. boost::math::uintmax_t max_iter = policies::get_max_root_iterations<Policy>();
  551. //
  552. // Pick an initial guess that we know will give us a probability
  553. // right around 0.5.
  554. //
  555. RealType guess = x - lam;
  556. if(guess < 1)
  557. guess = 1;
  558. boost::math::pair<RealType, RealType> ir = tools::bracket_and_solve_root(
  559. f, guess, RealType(2), false, tol, max_iter, pol);
  560. RealType result = ir.first + (ir.second - ir.first) / 2;
  561. if(max_iter >= policies::get_max_root_iterations<Policy>())
  562. {
  563. return policies::raise_evaluation_error<RealType>(function, "Unable to locate solution in a reasonable time:" // LCOV_EXCL_LINE
  564. " or there is no answer to problem. Current best guess is %1%", result, Policy()); // LCOV_EXCL_LINE
  565. }
  566. return result;
  567. }
  568. template <class RealType, class Policy>
  569. struct non_centrality_finder
  570. {
  571. BOOST_MATH_GPU_ENABLED non_centrality_finder(
  572. RealType v_, RealType x_, RealType p_, bool c)
  573. : v(v_), x(x_), p(p_), comp(c) {}
  574. BOOST_MATH_GPU_ENABLED RealType operator()(const RealType& lam)
  575. {
  576. non_central_chi_squared_distribution<RealType, Policy> d(v, lam);
  577. return comp ?
  578. RealType(p - cdf(complement(d, x)))
  579. : RealType(cdf(d, x) - p);
  580. }
  581. private:
  582. RealType v;
  583. RealType x;
  584. RealType p;
  585. bool comp;
  586. };
  587. template <class RealType, class Policy>
  588. BOOST_MATH_GPU_ENABLED inline RealType find_non_centrality(
  589. RealType v, RealType x, RealType p, RealType q, const Policy& pol)
  590. {
  591. constexpr auto function = "non_central_chi_squared<%1%>::find_non_centrality";
  592. if((p == 0) || (q == 0))
  593. {
  594. //
  595. // Can't do a thing if one of p and q is zero:
  596. //
  597. return policies::raise_evaluation_error<RealType>(function, "Can't find non centrality parameter when the probability is 0 or 1, only possible answer is %1%", // LCOV_EXCL_LINE
  598. RealType(boost::math::numeric_limits<RealType>::quiet_NaN()), Policy()); // LCOV_EXCL_LINE
  599. }
  600. non_centrality_finder<RealType, Policy> f(v, x, p < q ? p : q, p < q ? false : true);
  601. tools::eps_tolerance<RealType> tol(policies::digits<RealType, Policy>());
  602. boost::math::uintmax_t max_iter = policies::get_max_root_iterations<Policy>();
  603. //
  604. // Pick an initial guess that we know will give us a probability
  605. // right around 0.5.
  606. //
  607. RealType guess = x - v;
  608. if(guess < 1)
  609. guess = 1;
  610. boost::math::pair<RealType, RealType> ir = tools::bracket_and_solve_root(
  611. f, guess, RealType(2), false, tol, max_iter, pol);
  612. RealType result = ir.first + (ir.second - ir.first) / 2;
  613. if(max_iter >= policies::get_max_root_iterations<Policy>())
  614. {
  615. return policies::raise_evaluation_error<RealType>(function, "Unable to locate solution in a reasonable time:" // LCOV_EXCL_LINE
  616. " or there is no answer to problem. Current best guess is %1%", result, Policy()); // LCOV_EXCL_LINE
  617. }
  618. return result;
  619. }
  620. }
  621. template <class RealType = double, class Policy = policies::policy<> >
  622. class non_central_chi_squared_distribution
  623. {
  624. public:
  625. typedef RealType value_type;
  626. typedef Policy policy_type;
  627. BOOST_MATH_GPU_ENABLED non_central_chi_squared_distribution(RealType df_, RealType lambda) : df(df_), ncp(lambda)
  628. {
  629. constexpr auto function = "boost::math::non_central_chi_squared_distribution<%1%>::non_central_chi_squared_distribution(%1%,%1%)";
  630. RealType r;
  631. detail::check_df(
  632. function,
  633. df, &r, Policy());
  634. detail::check_non_centrality(
  635. function,
  636. ncp,
  637. &r,
  638. Policy());
  639. } // non_central_chi_squared_distribution constructor.
  640. BOOST_MATH_GPU_ENABLED RealType degrees_of_freedom() const
  641. { // Private data getter function.
  642. return df;
  643. }
  644. BOOST_MATH_GPU_ENABLED RealType non_centrality() const
  645. { // Private data getter function.
  646. return ncp;
  647. }
  648. BOOST_MATH_GPU_ENABLED static RealType find_degrees_of_freedom(RealType lam, RealType x, RealType p)
  649. {
  650. constexpr auto function = "non_central_chi_squared<%1%>::find_degrees_of_freedom";
  651. typedef typename policies::evaluation<RealType, Policy>::type eval_type;
  652. typedef typename policies::normalise<
  653. Policy,
  654. policies::promote_float<false>,
  655. policies::promote_double<false>,
  656. policies::discrete_quantile<>,
  657. policies::assert_undefined<> >::type forwarding_policy;
  658. eval_type result = detail::find_degrees_of_freedom(
  659. static_cast<eval_type>(lam),
  660. static_cast<eval_type>(x),
  661. static_cast<eval_type>(p),
  662. static_cast<eval_type>(1-p),
  663. forwarding_policy());
  664. return policies::checked_narrowing_cast<RealType, forwarding_policy>(
  665. result,
  666. function);
  667. }
  668. template <class A, class B, class C>
  669. BOOST_MATH_GPU_ENABLED static RealType find_degrees_of_freedom(const complemented3_type<A,B,C>& c)
  670. {
  671. constexpr auto function = "non_central_chi_squared<%1%>::find_degrees_of_freedom";
  672. typedef typename policies::evaluation<RealType, Policy>::type eval_type;
  673. typedef typename policies::normalise<
  674. Policy,
  675. policies::promote_float<false>,
  676. policies::promote_double<false>,
  677. policies::discrete_quantile<>,
  678. policies::assert_undefined<> >::type forwarding_policy;
  679. eval_type result = detail::find_degrees_of_freedom(
  680. static_cast<eval_type>(c.dist),
  681. static_cast<eval_type>(c.param1),
  682. static_cast<eval_type>(1-c.param2),
  683. static_cast<eval_type>(c.param2),
  684. forwarding_policy());
  685. return policies::checked_narrowing_cast<RealType, forwarding_policy>(
  686. result,
  687. function);
  688. }
  689. BOOST_MATH_GPU_ENABLED static RealType find_non_centrality(RealType v, RealType x, RealType p)
  690. {
  691. constexpr auto function = "non_central_chi_squared<%1%>::find_non_centrality";
  692. typedef typename policies::evaluation<RealType, Policy>::type eval_type;
  693. typedef typename policies::normalise<
  694. Policy,
  695. policies::promote_float<false>,
  696. policies::promote_double<false>,
  697. policies::discrete_quantile<>,
  698. policies::assert_undefined<> >::type forwarding_policy;
  699. eval_type result = detail::find_non_centrality(
  700. static_cast<eval_type>(v),
  701. static_cast<eval_type>(x),
  702. static_cast<eval_type>(p),
  703. static_cast<eval_type>(1-p),
  704. forwarding_policy());
  705. return policies::checked_narrowing_cast<RealType, forwarding_policy>(
  706. result,
  707. function);
  708. }
  709. template <class A, class B, class C>
  710. BOOST_MATH_GPU_ENABLED static RealType find_non_centrality(const complemented3_type<A,B,C>& c)
  711. {
  712. constexpr auto function = "non_central_chi_squared<%1%>::find_non_centrality";
  713. typedef typename policies::evaluation<RealType, Policy>::type eval_type;
  714. typedef typename policies::normalise<
  715. Policy,
  716. policies::promote_float<false>,
  717. policies::promote_double<false>,
  718. policies::discrete_quantile<>,
  719. policies::assert_undefined<> >::type forwarding_policy;
  720. eval_type result = detail::find_non_centrality(
  721. static_cast<eval_type>(c.dist),
  722. static_cast<eval_type>(c.param1),
  723. static_cast<eval_type>(1-c.param2),
  724. static_cast<eval_type>(c.param2),
  725. forwarding_policy());
  726. return policies::checked_narrowing_cast<RealType, forwarding_policy>(
  727. result,
  728. function);
  729. }
  730. private:
  731. // Data member, initialized by constructor.
  732. RealType df; // degrees of freedom.
  733. RealType ncp; // non-centrality parameter
  734. }; // template <class RealType, class Policy> class non_central_chi_squared_distribution
  735. typedef non_central_chi_squared_distribution<double> non_central_chi_squared; // Reserved name of type double.
  736. #ifdef __cpp_deduction_guides
  737. template <class RealType>
  738. non_central_chi_squared_distribution(RealType,RealType)->non_central_chi_squared_distribution<typename boost::math::tools::promote_args<RealType>::type>;
  739. #endif
  740. // Non-member functions to give properties of the distribution.
  741. template <class RealType, class Policy>
  742. BOOST_MATH_GPU_ENABLED inline const boost::math::pair<RealType, RealType> range(const non_central_chi_squared_distribution<RealType, Policy>& /* dist */)
  743. { // Range of permissible values for random variable k.
  744. using boost::math::tools::max_value;
  745. return boost::math::pair<RealType, RealType>(static_cast<RealType>(0), max_value<RealType>()); // Max integer?
  746. }
  747. template <class RealType, class Policy>
  748. BOOST_MATH_GPU_ENABLED inline const boost::math::pair<RealType, RealType> support(const non_central_chi_squared_distribution<RealType, Policy>& /* dist */)
  749. { // Range of supported values for random variable k.
  750. // This is range where cdf rises from 0 to 1, and outside it, the pdf is zero.
  751. using boost::math::tools::max_value;
  752. return boost::math::pair<RealType, RealType>(static_cast<RealType>(0), max_value<RealType>());
  753. }
  754. template <class RealType, class Policy>
  755. BOOST_MATH_GPU_ENABLED inline RealType mean(const non_central_chi_squared_distribution<RealType, Policy>& dist)
  756. { // Mean of poisson distribution = lambda.
  757. constexpr auto function = "boost::math::non_central_chi_squared_distribution<%1%>::mean()";
  758. RealType k = dist.degrees_of_freedom();
  759. RealType l = dist.non_centrality();
  760. RealType r;
  761. if(!detail::check_df(
  762. function,
  763. k, &r, Policy())
  764. ||
  765. !detail::check_non_centrality(
  766. function,
  767. l,
  768. &r,
  769. Policy()))
  770. return static_cast<RealType>(r);
  771. return k + l;
  772. } // mean
  773. template <class RealType, class Policy>
  774. BOOST_MATH_GPU_ENABLED inline RealType mode(const non_central_chi_squared_distribution<RealType, Policy>& dist)
  775. { // mode.
  776. constexpr auto function = "mode(non_central_chi_squared_distribution<%1%> const&)";
  777. RealType k = dist.degrees_of_freedom();
  778. RealType l = dist.non_centrality();
  779. RealType r;
  780. if(!detail::check_df(
  781. function,
  782. k, &r, Policy())
  783. ||
  784. !detail::check_non_centrality(
  785. function,
  786. l,
  787. &r,
  788. Policy()))
  789. return static_cast<RealType>(r);
  790. bool asymptotic_mode = k < l/4;
  791. RealType starting_point = asymptotic_mode ? k + l - RealType(3) : RealType(1) + k;
  792. return detail::generic_find_mode(dist, starting_point, function);
  793. }
  794. template <class RealType, class Policy>
  795. BOOST_MATH_GPU_ENABLED inline RealType variance(const non_central_chi_squared_distribution<RealType, Policy>& dist)
  796. { // variance.
  797. constexpr auto function = "boost::math::non_central_chi_squared_distribution<%1%>::variance()";
  798. RealType k = dist.degrees_of_freedom();
  799. RealType l = dist.non_centrality();
  800. RealType r;
  801. if(!detail::check_df(
  802. function,
  803. k, &r, Policy())
  804. ||
  805. !detail::check_non_centrality(
  806. function,
  807. l,
  808. &r,
  809. Policy()))
  810. return static_cast<RealType>(r);
  811. return 2 * (2 * l + k);
  812. }
  813. // RealType standard_deviation(const non_central_chi_squared_distribution<RealType, Policy>& dist)
  814. // standard_deviation provided by derived accessors.
  815. template <class RealType, class Policy>
  816. BOOST_MATH_GPU_ENABLED inline RealType skewness(const non_central_chi_squared_distribution<RealType, Policy>& dist)
  817. { // skewness = sqrt(l).
  818. constexpr auto function = "boost::math::non_central_chi_squared_distribution<%1%>::skewness()";
  819. RealType k = dist.degrees_of_freedom();
  820. RealType l = dist.non_centrality();
  821. RealType r;
  822. if(!detail::check_df(
  823. function,
  824. k, &r, Policy())
  825. ||
  826. !detail::check_non_centrality(
  827. function,
  828. l,
  829. &r,
  830. Policy()))
  831. return static_cast<RealType>(r);
  832. BOOST_MATH_STD_USING
  833. return pow(2 / (k + 2 * l), RealType(3)/2) * (k + 3 * l);
  834. }
  835. template <class RealType, class Policy>
  836. BOOST_MATH_GPU_ENABLED inline RealType kurtosis_excess(const non_central_chi_squared_distribution<RealType, Policy>& dist)
  837. {
  838. constexpr auto function = "boost::math::non_central_chi_squared_distribution<%1%>::kurtosis_excess()";
  839. RealType k = dist.degrees_of_freedom();
  840. RealType l = dist.non_centrality();
  841. RealType r;
  842. if(!detail::check_df(
  843. function,
  844. k, &r, Policy())
  845. ||
  846. !detail::check_non_centrality(
  847. function,
  848. l,
  849. &r,
  850. Policy()))
  851. return static_cast<RealType>(r);
  852. return 12 * (k + 4 * l) / ((k + 2 * l) * (k + 2 * l));
  853. } // kurtosis_excess
  854. template <class RealType, class Policy>
  855. BOOST_MATH_GPU_ENABLED inline RealType kurtosis(const non_central_chi_squared_distribution<RealType, Policy>& dist)
  856. {
  857. return kurtosis_excess(dist) + 3;
  858. }
  859. template <class RealType, class Policy>
  860. BOOST_MATH_GPU_ENABLED inline RealType pdf(const non_central_chi_squared_distribution<RealType, Policy>& dist, const RealType& x)
  861. { // Probability Density/Mass Function.
  862. return detail::nccs_pdf(dist, x);
  863. } // pdf
  864. template <class RealType, class Policy>
  865. BOOST_MATH_GPU_ENABLED RealType cdf(const non_central_chi_squared_distribution<RealType, Policy>& dist, const RealType& x)
  866. {
  867. constexpr auto function = "boost::math::non_central_chi_squared_distribution<%1%>::cdf(%1%)";
  868. RealType k = dist.degrees_of_freedom();
  869. RealType l = dist.non_centrality();
  870. RealType r;
  871. if(!detail::check_df(
  872. function,
  873. k, &r, Policy())
  874. ||
  875. !detail::check_non_centrality(
  876. function,
  877. l,
  878. &r,
  879. Policy())
  880. ||
  881. !detail::check_positive_x(
  882. function,
  883. x,
  884. &r,
  885. Policy()))
  886. return static_cast<RealType>(r);
  887. return detail::non_central_chi_squared_cdf(x, k, l, false, Policy());
  888. } // cdf
  889. template <class RealType, class Policy>
  890. BOOST_MATH_GPU_ENABLED RealType cdf(const complemented2_type<non_central_chi_squared_distribution<RealType, Policy>, RealType>& c)
  891. { // Complemented Cumulative Distribution Function
  892. constexpr auto function = "boost::math::non_central_chi_squared_distribution<%1%>::cdf(%1%)";
  893. non_central_chi_squared_distribution<RealType, Policy> const& dist = c.dist;
  894. RealType x = c.param;
  895. RealType k = dist.degrees_of_freedom();
  896. RealType l = dist.non_centrality();
  897. RealType r;
  898. if(!detail::check_df(
  899. function,
  900. k, &r, Policy())
  901. ||
  902. !detail::check_non_centrality(
  903. function,
  904. l,
  905. &r,
  906. Policy())
  907. ||
  908. !detail::check_positive_x(
  909. function,
  910. x,
  911. &r,
  912. Policy()))
  913. return static_cast<RealType>(r);
  914. return detail::non_central_chi_squared_cdf(x, k, l, true, Policy());
  915. } // ccdf
  916. template <class RealType, class Policy>
  917. BOOST_MATH_GPU_ENABLED inline RealType quantile(const non_central_chi_squared_distribution<RealType, Policy>& dist, const RealType& p)
  918. { // Quantile (or Percent Point) function.
  919. return detail::nccs_quantile(dist, p, false);
  920. } // quantile
  921. template <class RealType, class Policy>
  922. BOOST_MATH_GPU_ENABLED inline RealType quantile(const complemented2_type<non_central_chi_squared_distribution<RealType, Policy>, RealType>& c)
  923. { // Quantile (or Percent Point) function.
  924. return detail::nccs_quantile(c.dist, c.param, true);
  925. } // quantile complement.
  926. } // namespace math
  927. } // namespace boost
  928. // This include must be at the end, *after* the accessors
  929. // for this distribution have been defined, in order to
  930. // keep compilers that support two-phase lookup happy.
  931. #include <boost/math/distributions/detail/derived_accessors.hpp>
  932. #endif // BOOST_MATH_SPECIAL_NON_CENTRAL_CHI_SQUARE_HPP