non_central_chi_squared.hpp 40 KB

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