special_functions.cpp 22 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2016 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. //
  6. // This Source Code Form is subject to the terms of the Mozilla
  7. // Public License v. 2.0. If a copy of the MPL was not distributed
  8. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  9. #include <limits.h>
  10. #include "main.h"
  11. #include "../Eigen/SpecialFunctions"
  12. // Hack to allow "implicit" conversions from double to Scalar via comma-initialization.
  13. template<typename Derived>
  14. Eigen::CommaInitializer<Derived> operator<<(Eigen::DenseBase<Derived>& dense, double v) {
  15. return (dense << static_cast<typename Derived::Scalar>(v));
  16. }
  17. template<typename XprType>
  18. Eigen::CommaInitializer<XprType>& operator,(Eigen::CommaInitializer<XprType>& ci, double v) {
  19. return (ci, static_cast<typename XprType::Scalar>(v));
  20. }
  21. template<typename X, typename Y>
  22. void verify_component_wise(const X& x, const Y& y)
  23. {
  24. for(Index i=0; i<x.size(); ++i)
  25. {
  26. if((numext::isfinite)(y(i)))
  27. VERIFY_IS_APPROX( x(i), y(i) );
  28. else if((numext::isnan)(y(i)))
  29. VERIFY((numext::isnan)(x(i)));
  30. else
  31. VERIFY_IS_EQUAL( x(i), y(i) );
  32. }
  33. }
  34. template<typename ArrayType> void array_special_functions()
  35. {
  36. using std::abs;
  37. using std::sqrt;
  38. typedef typename ArrayType::Scalar Scalar;
  39. typedef typename NumTraits<Scalar>::Real RealScalar;
  40. Scalar plusinf = std::numeric_limits<Scalar>::infinity();
  41. Scalar nan = std::numeric_limits<Scalar>::quiet_NaN();
  42. Index rows = internal::random<Index>(1,30);
  43. Index cols = 1;
  44. // API
  45. {
  46. ArrayType m1 = ArrayType::Random(rows,cols);
  47. #if EIGEN_HAS_C99_MATH
  48. VERIFY_IS_APPROX(m1.lgamma(), lgamma(m1));
  49. VERIFY_IS_APPROX(m1.digamma(), digamma(m1));
  50. VERIFY_IS_APPROX(m1.erf(), erf(m1));
  51. VERIFY_IS_APPROX(m1.erfc(), erfc(m1));
  52. #endif // EIGEN_HAS_C99_MATH
  53. }
  54. #if EIGEN_HAS_C99_MATH
  55. // check special functions (comparing against numpy implementation)
  56. if (!NumTraits<Scalar>::IsComplex)
  57. {
  58. {
  59. ArrayType m1 = ArrayType::Random(rows,cols);
  60. ArrayType m2 = ArrayType::Random(rows,cols);
  61. // Test various propreties of igamma & igammac. These are normalized
  62. // gamma integrals where
  63. // igammac(a, x) = Gamma(a, x) / Gamma(a)
  64. // igamma(a, x) = gamma(a, x) / Gamma(a)
  65. // where Gamma and gamma are considered the standard unnormalized
  66. // upper and lower incomplete gamma functions, respectively.
  67. ArrayType a = m1.abs() + Scalar(2);
  68. ArrayType x = m2.abs() + Scalar(2);
  69. ArrayType zero = ArrayType::Zero(rows, cols);
  70. ArrayType one = ArrayType::Constant(rows, cols, Scalar(1.0));
  71. ArrayType a_m1 = a - one;
  72. ArrayType Gamma_a_x = Eigen::igammac(a, x) * a.lgamma().exp();
  73. ArrayType Gamma_a_m1_x = Eigen::igammac(a_m1, x) * a_m1.lgamma().exp();
  74. ArrayType gamma_a_x = Eigen::igamma(a, x) * a.lgamma().exp();
  75. ArrayType gamma_a_m1_x = Eigen::igamma(a_m1, x) * a_m1.lgamma().exp();
  76. // Gamma(a, 0) == Gamma(a)
  77. VERIFY_IS_APPROX(Eigen::igammac(a, zero), one);
  78. // Gamma(a, x) + gamma(a, x) == Gamma(a)
  79. VERIFY_IS_APPROX(Gamma_a_x + gamma_a_x, a.lgamma().exp());
  80. // Gamma(a, x) == (a - 1) * Gamma(a-1, x) + x^(a-1) * exp(-x)
  81. VERIFY_IS_APPROX(Gamma_a_x, (a - Scalar(1)) * Gamma_a_m1_x + x.pow(a-Scalar(1)) * (-x).exp());
  82. // gamma(a, x) == (a - 1) * gamma(a-1, x) - x^(a-1) * exp(-x)
  83. VERIFY_IS_APPROX(gamma_a_x, (a - Scalar(1)) * gamma_a_m1_x - x.pow(a-Scalar(1)) * (-x).exp());
  84. }
  85. {
  86. // Verify for large a and x that values are between 0 and 1.
  87. ArrayType m1 = ArrayType::Random(rows,cols);
  88. ArrayType m2 = ArrayType::Random(rows,cols);
  89. int max_exponent = std::numeric_limits<Scalar>::max_exponent10;
  90. ArrayType a = m1.abs() * Scalar(pow(10., max_exponent - 1));
  91. ArrayType x = m2.abs() * Scalar(pow(10., max_exponent - 1));
  92. for (int i = 0; i < a.size(); ++i) {
  93. Scalar igam = numext::igamma(a(i), x(i));
  94. VERIFY(0 <= igam);
  95. VERIFY(igam <= 1);
  96. }
  97. }
  98. {
  99. // Check exact values of igamma and igammac against a third party calculation.
  100. Scalar a_s[] = {Scalar(0), Scalar(1), Scalar(1.5), Scalar(4), Scalar(0.0001), Scalar(1000.5)};
  101. Scalar x_s[] = {Scalar(0), Scalar(1), Scalar(1.5), Scalar(4), Scalar(0.0001), Scalar(1000.5)};
  102. // location i*6+j corresponds to a_s[i], x_s[j].
  103. Scalar igamma_s[][6] = {
  104. {Scalar(0.0), nan, nan, nan, nan, nan},
  105. {Scalar(0.0), Scalar(0.6321205588285578), Scalar(0.7768698398515702),
  106. Scalar(0.9816843611112658), Scalar(9.999500016666262e-05),
  107. Scalar(1.0)},
  108. {Scalar(0.0), Scalar(0.4275932955291202), Scalar(0.608374823728911),
  109. Scalar(0.9539882943107686), Scalar(7.522076445089201e-07),
  110. Scalar(1.0)},
  111. {Scalar(0.0), Scalar(0.01898815687615381),
  112. Scalar(0.06564245437845008), Scalar(0.5665298796332909),
  113. Scalar(4.166333347221828e-18), Scalar(1.0)},
  114. {Scalar(0.0), Scalar(0.9999780593618628), Scalar(0.9999899967080838),
  115. Scalar(0.9999996219837988), Scalar(0.9991370418689945), Scalar(1.0)},
  116. {Scalar(0.0), Scalar(0.0), Scalar(0.0), Scalar(0.0), Scalar(0.0),
  117. Scalar(0.5042041932513908)}};
  118. Scalar igammac_s[][6] = {
  119. {nan, nan, nan, nan, nan, nan},
  120. {Scalar(1.0), Scalar(0.36787944117144233),
  121. Scalar(0.22313016014842982), Scalar(0.018315638888734182),
  122. Scalar(0.9999000049998333), Scalar(0.0)},
  123. {Scalar(1.0), Scalar(0.5724067044708798), Scalar(0.3916251762710878),
  124. Scalar(0.04601170568923136), Scalar(0.9999992477923555),
  125. Scalar(0.0)},
  126. {Scalar(1.0), Scalar(0.9810118431238462), Scalar(0.9343575456215499),
  127. Scalar(0.4334701203667089), Scalar(1.0), Scalar(0.0)},
  128. {Scalar(1.0), Scalar(2.1940638138146658e-05),
  129. Scalar(1.0003291916285e-05), Scalar(3.7801620118431334e-07),
  130. Scalar(0.0008629581310054535), Scalar(0.0)},
  131. {Scalar(1.0), Scalar(1.0), Scalar(1.0), Scalar(1.0), Scalar(1.0),
  132. Scalar(0.49579580674813944)}};
  133. for (int i = 0; i < 6; ++i) {
  134. for (int j = 0; j < 6; ++j) {
  135. if ((std::isnan)(igamma_s[i][j])) {
  136. VERIFY((std::isnan)(numext::igamma(a_s[i], x_s[j])));
  137. } else {
  138. VERIFY_IS_APPROX(numext::igamma(a_s[i], x_s[j]), igamma_s[i][j]);
  139. }
  140. if ((std::isnan)(igammac_s[i][j])) {
  141. VERIFY((std::isnan)(numext::igammac(a_s[i], x_s[j])));
  142. } else {
  143. VERIFY_IS_APPROX(numext::igammac(a_s[i], x_s[j]), igammac_s[i][j]);
  144. }
  145. }
  146. }
  147. }
  148. }
  149. #endif // EIGEN_HAS_C99_MATH
  150. // Check the ndtri function against scipy.special.ndtri
  151. {
  152. ArrayType x(7), res(7), ref(7);
  153. x << 0.5, 0.2, 0.8, 0.9, 0.1, 0.99, 0.01;
  154. ref << 0., -0.8416212335729142, 0.8416212335729142, 1.2815515655446004, -1.2815515655446004, 2.3263478740408408, -2.3263478740408408;
  155. CALL_SUBTEST( verify_component_wise(ref, ref); );
  156. CALL_SUBTEST( res = x.ndtri(); verify_component_wise(res, ref); );
  157. CALL_SUBTEST( res = ndtri(x); verify_component_wise(res, ref); );
  158. // ndtri(normal_cdf(x)) ~= x
  159. CALL_SUBTEST(
  160. ArrayType m1 = ArrayType::Random(32);
  161. using std::sqrt;
  162. ArrayType cdf_val = (m1 / Scalar(sqrt(2.))).erf();
  163. cdf_val = (cdf_val + Scalar(1)) / Scalar(2);
  164. verify_component_wise(cdf_val.ndtri(), m1););
  165. }
  166. // Check the zeta function against scipy.special.zeta
  167. {
  168. ArrayType x(10), q(10), res(10), ref(10);
  169. x << 1.5, 4, 10.5, 10000.5, 3, 1, 0.9, 2, 3, 4;
  170. q << 2, 1.5, 3, 1.0001, -2.5, 1.2345, 1.2345, -1, -2, -3;
  171. ref << 1.61237534869, 0.234848505667, 1.03086757337e-5, 0.367879440865, 0.054102025820864097, plusinf, nan, plusinf, nan, plusinf;
  172. CALL_SUBTEST( verify_component_wise(ref, ref); );
  173. CALL_SUBTEST( res = x.zeta(q); verify_component_wise(res, ref); );
  174. CALL_SUBTEST( res = zeta(x,q); verify_component_wise(res, ref); );
  175. }
  176. // digamma
  177. {
  178. ArrayType x(9), res(9), ref(9);
  179. x << 1, 1.5, 4, -10.5, 10000.5, 0, -1, -2, -3;
  180. ref << -0.5772156649015329, 0.03648997397857645, 1.2561176684318, 2.398239129535781, 9.210340372392849, nan, nan, nan, nan;
  181. CALL_SUBTEST( verify_component_wise(ref, ref); );
  182. CALL_SUBTEST( res = x.digamma(); verify_component_wise(res, ref); );
  183. CALL_SUBTEST( res = digamma(x); verify_component_wise(res, ref); );
  184. }
  185. #if EIGEN_HAS_C99_MATH
  186. {
  187. ArrayType n(16), x(16), res(16), ref(16);
  188. n << 1, 1, 1, 1.5, 17, 31, 28, 8, 42, 147, 170, -1, 0, 1, 2, 3;
  189. x << 2, 3, 25.5, 1.5, 4.7, 11.8, 17.7, 30.2, 15.8, 54.1, 64, -1, -2, -3, -4, -5;
  190. ref << 0.644934066848, 0.394934066848, 0.0399946696496, nan, 293.334565435, 0.445487887616, -2.47810300902e-07, -8.29668781082e-09, -0.434562276666, 0.567742190178, -0.0108615497927, nan, nan, plusinf, nan, plusinf;
  191. CALL_SUBTEST( verify_component_wise(ref, ref); );
  192. if(sizeof(RealScalar)>=8) { // double
  193. // Reason for commented line: http://eigen.tuxfamily.org/bz/show_bug.cgi?id=1232
  194. // CALL_SUBTEST( res = x.polygamma(n); verify_component_wise(res, ref); );
  195. CALL_SUBTEST( res = polygamma(n,x); verify_component_wise(res, ref); );
  196. }
  197. else {
  198. // CALL_SUBTEST( res = x.polygamma(n); verify_component_wise(res.head(8), ref.head(8)); );
  199. CALL_SUBTEST( res = polygamma(n,x); verify_component_wise(res.head(8), ref.head(8)); );
  200. }
  201. }
  202. #endif
  203. #if EIGEN_HAS_C99_MATH
  204. {
  205. // Inputs and ground truth generated with scipy via:
  206. // a = np.logspace(-3, 3, 5) - 1e-3
  207. // b = np.logspace(-3, 3, 5) - 1e-3
  208. // x = np.linspace(-0.1, 1.1, 5)
  209. // (full_a, full_b, full_x) = np.vectorize(lambda a, b, x: (a, b, x))(*np.ix_(a, b, x))
  210. // full_a = full_a.flatten().tolist() # same for full_b, full_x
  211. // v = scipy.special.betainc(full_a, full_b, full_x).flatten().tolist()
  212. //
  213. // Note in Eigen, we call betainc with arguments in the order (x, a, b).
  214. ArrayType a(125);
  215. ArrayType b(125);
  216. ArrayType x(125);
  217. ArrayType v(125);
  218. ArrayType res(125);
  219. a << 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,
  220. 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,
  221. 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
  222. 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
  223. 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
  224. 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
  225. 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
  226. 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
  227. 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
  228. 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
  229. 0.03062277660168379, 0.999, 0.999, 0.999, 0.999, 0.999, 0.999, 0.999,
  230. 0.999, 0.999, 0.999, 0.999, 0.999, 0.999, 0.999, 0.999, 0.999, 0.999,
  231. 0.999, 0.999, 0.999, 0.999, 0.999, 0.999, 0.999, 0.999,
  232. 31.62177660168379, 31.62177660168379, 31.62177660168379,
  233. 31.62177660168379, 31.62177660168379, 31.62177660168379,
  234. 31.62177660168379, 31.62177660168379, 31.62177660168379,
  235. 31.62177660168379, 31.62177660168379, 31.62177660168379,
  236. 31.62177660168379, 31.62177660168379, 31.62177660168379,
  237. 31.62177660168379, 31.62177660168379, 31.62177660168379,
  238. 31.62177660168379, 31.62177660168379, 31.62177660168379,
  239. 31.62177660168379, 31.62177660168379, 31.62177660168379,
  240. 31.62177660168379, 999.999, 999.999, 999.999, 999.999, 999.999, 999.999,
  241. 999.999, 999.999, 999.999, 999.999, 999.999, 999.999, 999.999, 999.999,
  242. 999.999, 999.999, 999.999, 999.999, 999.999, 999.999, 999.999, 999.999,
  243. 999.999, 999.999, 999.999;
  244. b << 0.0, 0.0, 0.0, 0.0, 0.0, 0.03062277660168379, 0.03062277660168379,
  245. 0.03062277660168379, 0.03062277660168379, 0.03062277660168379, 0.999,
  246. 0.999, 0.999, 0.999, 0.999, 31.62177660168379, 31.62177660168379,
  247. 31.62177660168379, 31.62177660168379, 31.62177660168379, 999.999,
  248. 999.999, 999.999, 999.999, 999.999, 0.0, 0.0, 0.0, 0.0, 0.0,
  249. 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
  250. 0.03062277660168379, 0.03062277660168379, 0.999, 0.999, 0.999, 0.999,
  251. 0.999, 31.62177660168379, 31.62177660168379, 31.62177660168379,
  252. 31.62177660168379, 31.62177660168379, 999.999, 999.999, 999.999,
  253. 999.999, 999.999, 0.0, 0.0, 0.0, 0.0, 0.0, 0.03062277660168379,
  254. 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
  255. 0.03062277660168379, 0.999, 0.999, 0.999, 0.999, 0.999,
  256. 31.62177660168379, 31.62177660168379, 31.62177660168379,
  257. 31.62177660168379, 31.62177660168379, 999.999, 999.999, 999.999,
  258. 999.999, 999.999, 0.0, 0.0, 0.0, 0.0, 0.0, 0.03062277660168379,
  259. 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
  260. 0.03062277660168379, 0.999, 0.999, 0.999, 0.999, 0.999,
  261. 31.62177660168379, 31.62177660168379, 31.62177660168379,
  262. 31.62177660168379, 31.62177660168379, 999.999, 999.999, 999.999,
  263. 999.999, 999.999, 0.0, 0.0, 0.0, 0.0, 0.0, 0.03062277660168379,
  264. 0.03062277660168379, 0.03062277660168379, 0.03062277660168379,
  265. 0.03062277660168379, 0.999, 0.999, 0.999, 0.999, 0.999,
  266. 31.62177660168379, 31.62177660168379, 31.62177660168379,
  267. 31.62177660168379, 31.62177660168379, 999.999, 999.999, 999.999,
  268. 999.999, 999.999;
  269. x << -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5,
  270. 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2,
  271. 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1,
  272. 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1,
  273. -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8,
  274. 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5,
  275. 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2,
  276. 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1,
  277. 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5, 0.8, 1.1, -0.1, 0.2, 0.5,
  278. 0.8, 1.1;
  279. v << nan, nan, nan, nan, nan, nan, nan, nan, nan, nan, nan, nan, nan, nan,
  280. nan, nan, nan, nan, nan, nan, nan, nan, nan, nan, nan, nan, nan, nan,
  281. nan, nan, nan, 0.47972119876364683, 0.5, 0.5202788012363533, nan, nan,
  282. 0.9518683957740043, 0.9789663010413743, 0.9931729188073435, nan, nan,
  283. 0.999995949033062, 0.9999999999993698, 0.9999999999999999, nan, nan,
  284. 0.9999999999999999, 0.9999999999999999, 0.9999999999999999, nan, nan,
  285. nan, nan, nan, nan, nan, 0.006827081192655869, 0.0210336989586256,
  286. 0.04813160422599567, nan, nan, 0.20014344256217678, 0.5000000000000001,
  287. 0.7998565574378232, nan, nan, 0.9991401428435834, 0.999999999698403,
  288. 0.9999999999999999, nan, nan, 0.9999999999999999, 0.9999999999999999,
  289. 0.9999999999999999, nan, nan, nan, nan, nan, nan, nan,
  290. 1.0646600232370887e-25, 6.301722877826246e-13, 4.050966937974938e-06,
  291. nan, nan, 7.864342668429763e-23, 3.015969667594166e-10,
  292. 0.0008598571564165444, nan, nan, 6.031987710123844e-08,
  293. 0.5000000000000007, 0.9999999396801229, nan, nan, 0.9999999999999999,
  294. 0.9999999999999999, 0.9999999999999999, nan, nan, nan, nan, nan, nan,
  295. nan, 0.0, 7.029920380986636e-306, 2.2450728208591345e-101, nan, nan,
  296. 0.0, 9.275871147869727e-302, 1.2232913026152827e-97, nan, nan, 0.0,
  297. 3.0891393081932924e-252, 2.9303043666183996e-60, nan, nan,
  298. 2.248913486879199e-196, 0.5000000000004947, 0.9999999999999999, nan;
  299. CALL_SUBTEST(res = betainc(a, b, x);
  300. verify_component_wise(res, v););
  301. }
  302. // Test various properties of betainc
  303. {
  304. ArrayType m1 = ArrayType::Random(32);
  305. ArrayType m2 = ArrayType::Random(32);
  306. ArrayType m3 = ArrayType::Random(32);
  307. ArrayType one = ArrayType::Constant(32, Scalar(1.0));
  308. const Scalar eps = std::numeric_limits<Scalar>::epsilon();
  309. ArrayType a = (m1 * Scalar(4)).exp();
  310. ArrayType b = (m2 * Scalar(4)).exp();
  311. ArrayType x = m3.abs();
  312. // betainc(a, 1, x) == x**a
  313. CALL_SUBTEST(
  314. ArrayType test = betainc(a, one, x);
  315. ArrayType expected = x.pow(a);
  316. verify_component_wise(test, expected););
  317. // betainc(1, b, x) == 1 - (1 - x)**b
  318. CALL_SUBTEST(
  319. ArrayType test = betainc(one, b, x);
  320. ArrayType expected = one - (one - x).pow(b);
  321. verify_component_wise(test, expected););
  322. // betainc(a, b, x) == 1 - betainc(b, a, 1-x)
  323. CALL_SUBTEST(
  324. ArrayType test = betainc(a, b, x) + betainc(b, a, one - x);
  325. ArrayType expected = one;
  326. verify_component_wise(test, expected););
  327. // betainc(a+1, b, x) = betainc(a, b, x) - x**a * (1 - x)**b / (a * beta(a, b))
  328. CALL_SUBTEST(
  329. ArrayType num = x.pow(a) * (one - x).pow(b);
  330. ArrayType denom = a * (a.lgamma() + b.lgamma() - (a + b).lgamma()).exp();
  331. // Add eps to rhs and lhs so that component-wise test doesn't result in
  332. // nans when both outputs are zeros.
  333. ArrayType expected = betainc(a, b, x) - num / denom + eps;
  334. ArrayType test = betainc(a + one, b, x) + eps;
  335. if (sizeof(Scalar) >= 8) { // double
  336. verify_component_wise(test, expected);
  337. } else {
  338. // Reason for limited test: http://eigen.tuxfamily.org/bz/show_bug.cgi?id=1232
  339. verify_component_wise(test.head(8), expected.head(8));
  340. });
  341. // betainc(a, b+1, x) = betainc(a, b, x) + x**a * (1 - x)**b / (b * beta(a, b))
  342. CALL_SUBTEST(
  343. // Add eps to rhs and lhs so that component-wise test doesn't result in
  344. // nans when both outputs are zeros.
  345. ArrayType num = x.pow(a) * (one - x).pow(b);
  346. ArrayType denom = b * (a.lgamma() + b.lgamma() - (a + b).lgamma()).exp();
  347. ArrayType expected = betainc(a, b, x) + num / denom + eps;
  348. ArrayType test = betainc(a, b + one, x) + eps;
  349. verify_component_wise(test, expected););
  350. }
  351. #endif // EIGEN_HAS_C99_MATH
  352. /* Code to generate the data for the following two test cases.
  353. N = 5
  354. np.random.seed(3)
  355. a = np.logspace(-2, 3, 6)
  356. a = np.ravel(np.tile(np.reshape(a, [-1, 1]), [1, N]))
  357. x = np.random.gamma(a, 1.0)
  358. x = np.maximum(x, np.finfo(np.float32).tiny)
  359. def igamma(a, x):
  360. return mpmath.gammainc(a, 0, x, regularized=True)
  361. def igamma_der_a(a, x):
  362. res = mpmath.diff(lambda a_prime: igamma(a_prime, x), a)
  363. return np.float64(res)
  364. def gamma_sample_der_alpha(a, x):
  365. igamma_x = igamma(a, x)
  366. def igammainv_of_igamma(a_prime):
  367. return mpmath.findroot(lambda x_prime: igamma(a_prime, x_prime) -
  368. igamma_x, x, solver='newton')
  369. return np.float64(mpmath.diff(igammainv_of_igamma, a))
  370. v_igamma_der_a = np.vectorize(igamma_der_a)(a, x)
  371. v_gamma_sample_der_alpha = np.vectorize(gamma_sample_der_alpha)(a, x)
  372. */
  373. #if EIGEN_HAS_C99_MATH
  374. // Test igamma_der_a
  375. {
  376. ArrayType a(30);
  377. ArrayType x(30);
  378. ArrayType res(30);
  379. ArrayType v(30);
  380. a << 0.01, 0.01, 0.01, 0.01, 0.01, 0.1, 0.1, 0.1, 0.1, 0.1, 1.0, 1.0, 1.0,
  381. 1.0, 1.0, 10.0, 10.0, 10.0, 10.0, 10.0, 100.0, 100.0, 100.0, 100.0,
  382. 100.0, 1000.0, 1000.0, 1000.0, 1000.0, 1000.0;
  383. x << 1.25668890405e-26, 1.17549435082e-38, 1.20938905072e-05,
  384. 1.17549435082e-38, 1.17549435082e-38, 5.66572070696e-16,
  385. 0.0132865061065, 0.0200034203853, 6.29263709118e-17, 1.37160367764e-06,
  386. 0.333412038288, 1.18135687766, 0.580629033777, 0.170631439426,
  387. 0.786686768458, 7.63873279537, 13.1944344379, 11.896042354,
  388. 10.5830172417, 10.5020942233, 92.8918587747, 95.003720371,
  389. 86.3715926467, 96.0330217672, 82.6389930677, 968.702906754,
  390. 969.463546828, 1001.79726022, 955.047416547, 1044.27458568;
  391. v << -32.7256441441, -36.4394150514, -9.66467612263, -36.4394150514,
  392. -36.4394150514, -1.0891900302, -2.66351229645, -2.48666868596,
  393. -0.929700494428, -3.56327722764, -0.455320135314, -0.391437214323,
  394. -0.491352055991, -0.350454834292, -0.471773162921, -0.104084440522,
  395. -0.0723646747909, -0.0992828975532, -0.121638215446, -0.122619605294,
  396. -0.0317670267286, -0.0359974812869, -0.0154359225363, -0.0375775365921,
  397. -0.00794899153653, -0.00777303219211, -0.00796085782042,
  398. -0.0125850719397, -0.00455500206958, -0.00476436993148;
  399. CALL_SUBTEST(res = igamma_der_a(a, x); verify_component_wise(res, v););
  400. }
  401. // Test gamma_sample_der_alpha
  402. {
  403. ArrayType alpha(30);
  404. ArrayType sample(30);
  405. ArrayType res(30);
  406. ArrayType v(30);
  407. alpha << 0.01, 0.01, 0.01, 0.01, 0.01, 0.1, 0.1, 0.1, 0.1, 0.1, 1.0, 1.0,
  408. 1.0, 1.0, 1.0, 10.0, 10.0, 10.0, 10.0, 10.0, 100.0, 100.0, 100.0, 100.0,
  409. 100.0, 1000.0, 1000.0, 1000.0, 1000.0, 1000.0;
  410. sample << 1.25668890405e-26, 1.17549435082e-38, 1.20938905072e-05,
  411. 1.17549435082e-38, 1.17549435082e-38, 5.66572070696e-16,
  412. 0.0132865061065, 0.0200034203853, 6.29263709118e-17, 1.37160367764e-06,
  413. 0.333412038288, 1.18135687766, 0.580629033777, 0.170631439426,
  414. 0.786686768458, 7.63873279537, 13.1944344379, 11.896042354,
  415. 10.5830172417, 10.5020942233, 92.8918587747, 95.003720371,
  416. 86.3715926467, 96.0330217672, 82.6389930677, 968.702906754,
  417. 969.463546828, 1001.79726022, 955.047416547, 1044.27458568;
  418. v << 7.42424742367e-23, 1.02004297287e-34, 0.0130155240738,
  419. 1.02004297287e-34, 1.02004297287e-34, 1.96505168277e-13, 0.525575786243,
  420. 0.713903991771, 2.32077561808e-14, 0.000179348049886, 0.635500453302,
  421. 1.27561284917, 0.878125852156, 0.41565819538, 1.03606488534,
  422. 0.885964824887, 1.16424049334, 1.10764479598, 1.04590810812,
  423. 1.04193666963, 0.965193152414, 0.976217589464, 0.93008035061,
  424. 0.98153216096, 0.909196397698, 0.98434963993, 0.984738050206,
  425. 1.00106492525, 0.97734200649, 1.02198794179;
  426. CALL_SUBTEST(res = gamma_sample_der_alpha(alpha, sample);
  427. verify_component_wise(res, v););
  428. }
  429. #endif // EIGEN_HAS_C99_MATH
  430. }
  431. EIGEN_DECLARE_TEST(special_functions)
  432. {
  433. CALL_SUBTEST_1(array_special_functions<ArrayXf>());
  434. CALL_SUBTEST_2(array_special_functions<ArrayXd>());
  435. // TODO(cantonios): half/bfloat16 don't have enough precision to reproduce results above.
  436. // CALL_SUBTEST_3(array_special_functions<ArrayX<Eigen::half>>());
  437. // CALL_SUBTEST_4(array_special_functions<ArrayX<Eigen::bfloat16>>());
  438. }