10 #ifndef EIGEN_SPECIAL_FUNCTIONS_H
11 #define EIGEN_SPECIAL_FUNCTIONS_H
46 template <
typename Scalar>
49 THIS_TYPE_IS_NOT_SUPPORTED)
56 template <
typename Scalar>
57 struct lgamma_retval {
61 #if EIGEN_HAS_C99_MATH
63 #if defined(__GLIBC__) && ((__GLIBC__>=2 && __GLIBC_MINOR__ >= 19) || __GLIBC__>2) \
64 && (defined(_DEFAULT_SOURCE) || defined(_BSD_SOURCE) || defined(_SVID_SOURCE))
65 #define EIGEN_HAS_LGAMMA_R
69 #if defined(__GLIBC__) && ((__GLIBC__==2 && __GLIBC_MINOR__ < 19) || __GLIBC__<2) \
70 && (defined(_BSD_SOURCE) || defined(_SVID_SOURCE))
71 #define EIGEN_HAS_LGAMMA_R
75 struct lgamma_impl<float> {
77 static EIGEN_STRONG_INLINE
float run(
float x) {
78 #if !defined(EIGEN_GPU_COMPILE_PHASE) && defined (EIGEN_HAS_LGAMMA_R) && !defined(__APPLE__)
80 return ::lgammaf_r(
x, &dummy);
81 #elif defined(SYCL_DEVICE_ONLY)
82 return cl::sycl::lgamma(
x);
90 struct lgamma_impl<double> {
92 static EIGEN_STRONG_INLINE
double run(
double x) {
93 #if !defined(EIGEN_GPU_COMPILE_PHASE) && defined(EIGEN_HAS_LGAMMA_R) && !defined(__APPLE__)
95 return ::lgamma_r(
x, &dummy);
96 #elif defined(SYCL_DEVICE_ONLY)
97 return cl::sycl::lgamma(
x);
104 #undef EIGEN_HAS_LGAMMA_R
111 template <
typename Scalar>
112 struct digamma_retval {
129 template <
typename Scalar>
130 struct digamma_impl_maybe_poly {
132 THIS_TYPE_IS_NOT_SUPPORTED)
141 struct digamma_impl_maybe_poly<float> {
143 static EIGEN_STRONG_INLINE
float run(
const float s) {
145 -4.16666666666666666667E-3f,
146 3.96825396825396825397E-3f,
147 -8.33333333333333333333E-3f,
148 8.33333333333333333333E-2f
154 return z * internal::ppolevl<float, 3>::run(z, A);
160 struct digamma_impl_maybe_poly<double> {
162 static EIGEN_STRONG_INLINE
double run(
const double s) {
164 8.33333333333333333333E-2,
165 -2.10927960927960927961E-2,
166 7.57575757575757575758E-3,
167 -4.16666666666666666667E-3,
168 3.96825396825396825397E-3,
169 -8.33333333333333333333E-3,
170 8.33333333333333333333E-2
176 return z * internal::ppolevl<double, 6>::run(z, A);
182 template <
typename Scalar>
183 struct digamma_impl {
185 static Scalar run(Scalar
x) {
243 Scalar
p,
q, nz, s,
w,
y;
244 bool negative =
false;
246 const Scalar nan = NumTraits<Scalar>::quiet_NaN();
247 const Scalar m_pi = Scalar(EIGEN_PI);
249 const Scalar zero = Scalar(0);
250 const Scalar one = Scalar(1);
251 const Scalar half = Scalar(0.5);
281 while (s < Scalar(10)) {
286 y = digamma_impl_maybe_poly<Scalar>::run(s);
290 return (negative) ?
y - nz :
y;
307 template <
typename T>
309 constexpr
float kErfInvOneMinusHalfULP = 3.832506856900711f;
310 const T clamp =
pcmp_le(pset1<T>(kErfInvOneMinusHalfULP),
pabs(
x));
312 const T alpha_1 = pset1<T>(1.128379143519084f);
313 const T alpha_3 = pset1<T>(0.18520832239976145f);
314 const T alpha_5 = pset1<T>(0.050955695062380861f);
315 const T alpha_7 = pset1<T>(0.0034082910107109506f);
316 const T alpha_9 = pset1<T>(0.00022905065861350646f);
319 const T beta_0 = pset1<T>(1.0f);
320 const T beta_2 = pset1<T>(0.49746925110067538f);
321 const T beta_4 = pset1<T>(0.11098505178285362f);
322 const T beta_6 = pset1<T>(0.014070470171167667f);
323 const T beta_8 = pset1<T>(0.0010179625278914885f);
324 const T beta_10 = pset1<T>(0.000023547966471313185f);
325 const T beta_12 = pset1<T>(-1.1791602954361697e-7f);
331 T p =
pmadd(x2, alpha_9, alpha_7);
338 T q =
pmadd(x2, beta_12, beta_10);
349 template <
typename T>
352 static EIGEN_STRONG_INLINE
T run(
const T&
x) {
357 template <
typename Scalar>
362 #if EIGEN_HAS_C99_MATH
364 struct erf_impl<float> {
366 static EIGEN_STRONG_INLINE
float run(
float x) {
367 #if defined(SYCL_DEVICE_ONLY)
368 return cl::sycl::erf(
x);
376 struct erf_impl<double> {
378 static EIGEN_STRONG_INLINE
double run(
double x) {
379 #if defined(SYCL_DEVICE_ONLY)
380 return cl::sycl::erf(
x);
392 template <
typename Scalar>
395 THIS_TYPE_IS_NOT_SUPPORTED)
402 template <
typename Scalar>
407 #if EIGEN_HAS_C99_MATH
409 struct erfc_impl<float> {
411 static EIGEN_STRONG_INLINE
float run(
const float x) {
412 #if defined(SYCL_DEVICE_ONLY)
413 return cl::sycl::erfc(
x);
421 struct erfc_impl<double> {
423 static EIGEN_STRONG_INLINE
double run(
const double x) {
424 #if defined(SYCL_DEVICE_ONLY)
425 return cl::sycl::erfc(
x);
493 const T& should_flipsign,
const T&
x) {
494 typedef typename unpacket_traits<T>::type Scalar;
495 const T sign_mask = pset1<T>(Scalar(-0.0));
496 T sign_bit = pand<T>(should_flipsign, sign_mask);
497 return pxor<T>(sign_bit,
x);
502 const double& should_flipsign,
const double&
x) {
503 return should_flipsign == 0 ?
x : -
x;
508 const float& should_flipsign,
const float&
x) {
509 return should_flipsign == 0 ?
x : -
x;
516 template <
typename T,
typename ScalarType>
518 const ScalarType
p0[] = {
519 ScalarType(-5.99633501014107895267e1),
520 ScalarType(9.80010754185999661536e1),
521 ScalarType(-5.66762857469070293439e1),
522 ScalarType(1.39312609387279679503e1),
523 ScalarType(-1.23916583867381258016e0)
525 const ScalarType q0[] = {
527 ScalarType(1.95448858338141759834e0),
528 ScalarType(4.67627912898881538453e0),
529 ScalarType(8.63602421390890590575e1),
530 ScalarType(-2.25462687854119370527e2),
531 ScalarType(2.00260212380060660359e2),
532 ScalarType(-8.20372256168333339912e1),
533 ScalarType(1.59056225126211695515e1),
534 ScalarType(-1.18331621121330003142e0)
536 const T sqrt2pi = pset1<T>(ScalarType(2.50662827463100050242e0));
537 const T half = pset1<T>(ScalarType(0.5));
538 T c, c2, ndtri_gt_exp_neg_two;
544 internal::ppolevl<T, 4>::run(c2,
p0),
545 internal::ppolevl<T, 8>::run(c2, q0))),
c);
546 return pmul(ndtri_gt_exp_neg_two, sqrt2pi);
549 template <
typename T,
typename ScalarType>
551 const T&
b,
const T& should_flipsign) {
555 const ScalarType
p1[] = {
556 ScalarType(4.05544892305962419923e0),
557 ScalarType(3.15251094599893866154e1),
558 ScalarType(5.71628192246421288162e1),
559 ScalarType(4.40805073893200834700e1),
560 ScalarType(1.46849561928858024014e1),
561 ScalarType(2.18663306850790267539e0),
562 ScalarType(-1.40256079171354495875e-1),
563 ScalarType(-3.50424626827848203418e-2),
564 ScalarType(-8.57456785154685413611e-4)
566 const ScalarType q1[] = {
568 ScalarType(1.57799883256466749731e1),
569 ScalarType(4.53907635128879210584e1),
570 ScalarType(4.13172038254672030440e1),
571 ScalarType(1.50425385692907503408e1),
572 ScalarType(2.50464946208309415979e0),
573 ScalarType(-1.42182922854787788574e-1),
574 ScalarType(-3.80806407691578277194e-2),
575 ScalarType(-9.33259480895457427372e-4)
580 const ScalarType p2[] = {
581 ScalarType(3.23774891776946035970e0),
582 ScalarType(6.91522889068984211695e0),
583 ScalarType(3.93881025292474443415e0),
584 ScalarType(1.33303460815807542389e0),
585 ScalarType(2.01485389549179081538e-1),
586 ScalarType(1.23716634817820021358e-2),
587 ScalarType(3.01581553508235416007e-4),
588 ScalarType(2.65806974686737550832e-6),
589 ScalarType(6.23974539184983293730e-9)
591 const ScalarType q2[] = {
593 ScalarType(6.02427039364742014255e0),
594 ScalarType(3.67983563856160859403e0),
595 ScalarType(1.37702099489081330271e0),
596 ScalarType(2.16236993594496635890e-1),
597 ScalarType(1.34204006088543189037e-2),
598 ScalarType(3.28014464682127739104e-4),
599 ScalarType(2.89247864745380683936e-6),
600 ScalarType(6.79019408009981274425e-9)
602 const T eight = pset1<T>(ScalarType(8.0));
603 const T neg_two = pset1<T>(ScalarType(-2));
612 pdiv(internal::ppolevl<T, 8>::run(z,
p1),
613 internal::ppolevl<T, 8>::run(z, q1)),
614 pdiv(internal::ppolevl<T, 8>::run(z, p2),
615 internal::ppolevl<T, 8>::run(z, q2))));
619 template <
typename T,
typename ScalarType>
625 const T zero = pset1<T>(ScalarType(0));
626 const T one = pset1<T>(ScalarType(1));
628 const T exp_neg_two = pset1<T>(ScalarType(0.13533528323661269189));
636 generic_ndtri_gt_exp_neg_two<T, ScalarType>(
b),
637 generic_ndtri_lt_exp_neg_two<T, ScalarType>(
b, should_flipsign));
644 template <
typename Scalar>
645 struct ndtri_retval {
649 #if !EIGEN_HAS_C99_MATH
651 template <
typename Scalar>
654 THIS_TYPE_IS_NOT_SUPPORTED)
663 template <
typename Scalar>
666 static EIGEN_STRONG_INLINE Scalar run(
const Scalar
x) {
667 return generic_ndtri<Scalar, Scalar>(
x);
678 template <
typename Scalar>
679 struct igammac_retval {
684 template <
typename Scalar>
685 struct cephes_helper {
687 static EIGEN_STRONG_INLINE Scalar machep() {
eigen_assert(
false &&
"machep not supported for this type");
return 0.0; }
689 static EIGEN_STRONG_INLINE Scalar big() {
eigen_assert(
false &&
"big not supported for this type");
return 0.0; }
691 static EIGEN_STRONG_INLINE Scalar biginv() {
eigen_assert(
false &&
"biginv not supported for this type");
return 0.0; }
695 struct cephes_helper<float> {
697 static EIGEN_STRONG_INLINE
float machep() {
698 return NumTraits<float>::epsilon() / 2;
701 static EIGEN_STRONG_INLINE
float big() {
703 return 1.0f / (NumTraits<float>::epsilon() / 2);
706 static EIGEN_STRONG_INLINE
float biginv() {
713 struct cephes_helper<double> {
715 static EIGEN_STRONG_INLINE
double machep() {
716 return NumTraits<double>::epsilon() / 2;
719 static EIGEN_STRONG_INLINE
double big() {
720 return 1.0 / NumTraits<double>::epsilon();
723 static EIGEN_STRONG_INLINE
double biginv() {
725 return NumTraits<double>::epsilon();
731 template <
typename Scalar>
744 template <
typename Scalar, IgammaComputationMode mode>
753 if (internal::is_same<Scalar, float>::value) {
755 }
else if (internal::is_same<Scalar, double>::value) {
762 template <
typename Scalar, IgammaComputationMode mode>
763 struct igammac_cf_impl {
774 static Scalar run(Scalar a, Scalar
x) {
775 const Scalar zero = 0;
776 const Scalar one = 1;
777 const Scalar two = 2;
778 const Scalar machep = cephes_helper<Scalar>::machep();
779 const Scalar big = cephes_helper<Scalar>::big();
780 const Scalar biginv = cephes_helper<Scalar>::biginv();
786 Scalar ax = main_igamma_term<Scalar>(a,
x);
797 Scalar z =
x +
y + one;
801 Scalar pkm1 =
x + one;
803 Scalar ans = pkm1 / qkm1;
805 Scalar dpkm2_da = zero;
806 Scalar dqkm2_da = zero;
807 Scalar dpkm1_da = zero;
808 Scalar dqkm1_da = -
x;
809 Scalar dans_da = (dpkm1_da - ans * dqkm1_da) / qkm1;
811 for (
int i = 0; i < igamma_num_iterations<Scalar, mode>();
i++) {
817 Scalar pk = pkm1 * z - pkm2 * yc;
818 Scalar qk = qkm1 * z - qkm2 * yc;
820 Scalar dpk_da = dpkm1_da * z - pkm1 - dpkm2_da * yc + pkm2 *
c;
821 Scalar dqk_da = dqkm1_da * z - qkm1 - dqkm2_da * yc + qkm2 *
c;
824 Scalar ans_prev = ans;
827 Scalar dans_da_prev = dans_da;
828 dans_da = (dpk_da - ans * dqk_da) / qk;
835 if (
numext::abs(dans_da - dans_da_prev) <= machep) {
865 Scalar dlogax_da =
numext::log(
x) - digamma_impl<Scalar>::run(a);
866 Scalar dax_da = ax * dlogax_da;
872 return ans * dax_da + dans_da * ax;
875 return -(dans_da + ans * dlogax_da) *
x;
880 template <
typename Scalar, IgammaComputationMode mode>
881 struct igamma_series_impl {
891 static Scalar run(Scalar a, Scalar
x) {
892 const Scalar zero = 0;
893 const Scalar one = 1;
894 const Scalar machep = cephes_helper<Scalar>::machep();
896 Scalar ax = main_igamma_term<Scalar>(a,
x);
914 Scalar dans_da = zero;
916 for (
int i = 0; i < igamma_num_iterations<Scalar, mode>();
i++) {
919 Scalar dterm_da = -
x / (r * r);
920 dc_da = term * dc_da + dterm_da *
c;
926 if (c <= machep * ans) {
936 Scalar dlogax_da =
numext::log(
x) - digamma_impl<Scalar>::run(a + one);
937 Scalar dax_da = ax * dlogax_da;
943 return ans * dax_da + dans_da * ax;
946 return -(dans_da + ans * dlogax_da) *
x / a;
951 #if !EIGEN_HAS_C99_MATH
953 template <
typename Scalar>
954 struct igammac_impl {
956 THIS_TYPE_IS_NOT_SUPPORTED)
965 template <
typename Scalar>
966 struct igammac_impl {
968 static Scalar run(Scalar a, Scalar
x) {
1023 const Scalar zero = 0;
1024 const Scalar one = 1;
1025 const Scalar nan = NumTraits<Scalar>::quiet_NaN();
1027 if ((
x < zero) || (a <= zero)) {
1036 if ((
x < one) || (
x < a)) {
1037 return (one - igamma_series_impl<Scalar, VALUE>::run(a,
x));
1040 return igammac_cf_impl<Scalar, VALUE>::run(a,
x);
1050 #if !EIGEN_HAS_C99_MATH
1052 template <
typename Scalar, IgammaComputationMode mode>
1053 struct igamma_generic_impl {
1055 THIS_TYPE_IS_NOT_SUPPORTED)
1064 template <
typename Scalar, IgammaComputationMode mode>
1065 struct igamma_generic_impl {
1067 static Scalar run(Scalar a, Scalar
x) {
1076 const Scalar zero = 0;
1077 const Scalar one = 1;
1078 const Scalar nan = NumTraits<Scalar>::quiet_NaN();
1080 if (
x == zero)
return zero;
1082 if ((
x < zero) || (a <= zero)) {
1090 if ((
x > one) && (
x > a)) {
1091 Scalar ret = igammac_cf_impl<Scalar, mode>::run(a,
x);
1092 if (mode ==
VALUE) {
1099 return igamma_series_impl<Scalar, mode>::run(a,
x);
1105 template <
typename Scalar>
1106 struct igamma_retval {
1107 typedef Scalar type;
1110 template <
typename Scalar>
1111 struct igamma_impl : igamma_generic_impl<Scalar, VALUE> {
1181 template <
typename Scalar>
1182 struct igamma_der_a_retval : igamma_retval<Scalar> {};
1184 template <
typename Scalar>
1185 struct igamma_der_a_impl : igamma_generic_impl<Scalar, DERIVATIVE> {
1202 template <
typename Scalar>
1203 struct gamma_sample_der_alpha_retval : igamma_retval<Scalar> {};
1205 template <
typename Scalar>
1206 struct gamma_sample_der_alpha_impl
1207 : igamma_generic_impl<Scalar, SAMPLE_DERIVATIVE> {
1251 template <
typename Scalar>
1252 struct zeta_retval {
1253 typedef Scalar type;
1256 template <
typename Scalar>
1257 struct zeta_impl_series {
1259 THIS_TYPE_IS_NOT_SUPPORTED)
1267 struct zeta_impl_series<float> {
1269 static EIGEN_STRONG_INLINE
bool run(
float& a,
float&
b,
float& s,
const float x,
const float machep) {
1287 struct zeta_impl_series<double> {
1289 static EIGEN_STRONG_INLINE
bool run(
double& a,
double&
b,
double& s,
const double x,
const double machep) {
1291 while( (i < 9) || (a <= 9.0) )
1306 template <
typename Scalar>
1309 static Scalar run(Scalar
x, Scalar
q) {
1372 Scalar
p, r,
a,
b, k, s, t,
w;
1374 const Scalar
A[] = {
1380 Scalar(-1.8924375803183791606e9),
1381 Scalar(7.47242496e10),
1382 Scalar(-2.950130727918164224e12),
1383 Scalar(1.1646782814350067249e14),
1384 Scalar(-4.5979787224074726105e15),
1385 Scalar(1.8152105401943546773e17),
1386 Scalar(-7.1661652561756670113e18)
1389 const Scalar maxnum = NumTraits<Scalar>::infinity();
1390 const Scalar zero = Scalar(0.0), half = Scalar(0.5), one = Scalar(1.0);
1391 const Scalar machep = cephes_helper<Scalar>::machep();
1392 const Scalar nan = NumTraits<Scalar>::quiet_NaN();
1428 if (zeta_impl_series<Scalar>::run(a,
b, s,
x, machep)) {
1445 for( i=0;
i<12;
i++ )
1468 template <
typename Scalar>
1469 struct polygamma_retval {
1470 typedef Scalar type;
1473 #if !EIGEN_HAS_C99_MATH
1475 template <
typename Scalar>
1476 struct polygamma_impl {
1478 THIS_TYPE_IS_NOT_SUPPORTED)
1487 template <
typename Scalar>
1488 struct polygamma_impl {
1490 static Scalar run(Scalar n, Scalar
x) {
1491 Scalar zero = 0.0, one = 1.0;
1492 Scalar nplus =
n + one;
1493 const Scalar nan = NumTraits<Scalar>::quiet_NaN();
1500 else if (n == zero) {
1501 return digamma_impl<Scalar>::run(
x);
1505 Scalar factorial =
numext::exp(lgamma_impl<Scalar>::run(nplus));
1506 return numext::pow(-one, nplus) * factorial * zeta_impl<Scalar>::run(nplus,
x);
1517 template <
typename Scalar>
1518 struct betainc_retval {
1519 typedef Scalar type;
1522 #if !EIGEN_HAS_C99_MATH
1524 template <
typename Scalar>
1525 struct betainc_impl {
1527 THIS_TYPE_IS_NOT_SUPPORTED)
1536 template <
typename Scalar>
1537 struct betainc_impl {
1539 THIS_TYPE_IS_NOT_SUPPORTED)
1541 EIGEN_DEVICE_FUNC static EIGEN_STRONG_INLINE Scalar run(Scalar, Scalar, Scalar) {
1617 template <
typename Scalar>
1618 struct incbeta_cfe {
1620 internal::is_same<Scalar, double>::value),
1621 THIS_TYPE_IS_NOT_SUPPORTED)
1623 EIGEN_DEVICE_FUNC static EIGEN_STRONG_INLINE Scalar run(Scalar a, Scalar
b, Scalar
x,
bool small_branch) {
1624 const Scalar big = cephes_helper<Scalar>::big();
1625 const Scalar machep = cephes_helper<Scalar>::machep();
1626 const Scalar biginv = cephes_helper<Scalar>::biginv();
1628 const Scalar zero = 0;
1629 const Scalar one = 1;
1630 const Scalar two = 2;
1632 Scalar xk, pk, pkm1, pkm2, qk, qkm1, qkm2;
1633 Scalar k1, k2, k3, k4, k5, k6, k7, k8, k26update;
1637 const int num_iters = (internal::is_same<Scalar, float>::value) ? 100 : 300;
1638 const Scalar thresh =
1639 (internal::is_same<Scalar, float>::value) ? machep : Scalar(3) * machep;
1640 Scalar r = (internal::is_same<Scalar, float>::value) ? zero : one;
1673 xk = -(
x * k1 * k2) / (k3 * k4);
1674 pk = pkm1 + pkm2 * xk;
1675 qk = qkm1 + qkm2 * xk;
1681 xk = (
x * k5 * k6) / (k7 * k8);
1682 pk = pkm1 + pkm2 * xk;
1683 qk = qkm1 + qkm2 * xk;
1718 }
while (++n < num_iters);
1725 template <
typename Scalar>
1726 struct betainc_helper {};
1729 struct betainc_helper<float> {
1733 float ans,
a,
b, t,
x, onemx;
1734 bool reversed_a_b =
false;
1739 if (xx > (aa / (aa + bb))) {
1740 reversed_a_b =
true;
1755 t = betainc_helper<float>::incbps(a,
b,
x);
1756 if (reversed_a_b) t = 1.0f - t;
1761 ans =
x * (
a +
b - 2.0f) / (a - 1.0f);
1763 ans = incbeta_cfe<float>::run(a,
b,
x,
true );
1766 ans = incbeta_cfe<float>::run(a,
b,
x,
false );
1771 lgamma_impl<float>::run(a) - lgamma_impl<float>::run(
b);
1775 if (reversed_a_b) t = 1.0f - t;
1780 static EIGEN_STRONG_INLINE
float incbps(
float a,
float b,
float x) {
1782 const float machep = cephes_helper<float>::machep();
1785 y -= lgamma_impl<float>::run(a) + lgamma_impl<float>::run(
b);
1786 y += lgamma_impl<float>::run(a +
b);
1806 struct betainc_impl<float> {
1808 static float run(
float a,
float b,
float x) {
1809 const float nan = NumTraits<float>::quiet_NaN();
1812 if (a <= 0.0f)
return nan;
1813 if (
b <= 0.0f)
return nan;
1814 if ((
x <= 0.0f) || (
x >= 1.0f)) {
1815 if (
x == 0.0f)
return 0.0f;
1816 if (
x == 1.0f)
return 1.0f;
1823 ans = betainc_helper<float>::incbsa(a + 1.0f,
b,
x);
1825 lgamma_impl<float>::run(a +
b) - lgamma_impl<float>::run(a + 1.0f) -
1826 lgamma_impl<float>::run(
b);
1829 return betainc_helper<float>::incbsa(a,
b,
x);
1835 struct betainc_helper<double> {
1837 static EIGEN_STRONG_INLINE
double incbps(
double a,
double b,
double x) {
1838 const double machep = cephes_helper<double>::machep();
1840 double s, t, u,
v,
n, t1, z, ai;
1851 u = (
n -
b) *
x / n;
1868 t = lgamma_impl<double>::run(a +
b) - lgamma_impl<double>::run(a) -
1875 struct betainc_impl<double> {
1877 static double run(
double aa,
double bb,
double xx) {
1878 const double nan = NumTraits<double>::quiet_NaN();
1879 const double machep = cephes_helper<double>::machep();
1882 double a,
b, t,
x, xc,
w,
y;
1883 bool reversed_a_b =
false;
1885 if (aa <= 0.0 || bb <= 0.0) {
1889 if ((xx <= 0.0) || (xx >= 1.0)) {
1890 if (xx == 0.0)
return (0.0);
1891 if (xx == 1.0)
return (1.0);
1896 if ((bb * xx) <= 1.0 && xx <= 0.95) {
1897 return betainc_helper<double>::incbps(aa, bb, xx);
1903 if (xx > (aa / (aa + bb))) {
1904 reversed_a_b =
true;
1916 if (reversed_a_b && (
b *
x) <= 1.0 &&
x <= 0.95) {
1917 t = betainc_helper<double>::incbps(a,
b,
x);
1927 y =
x * (
a +
b - 2.0) - (a - 1.0);
1929 w = incbeta_cfe<double>::run(a,
b,
x,
true );
1931 w = incbeta_cfe<double>::run(a,
b,
x,
false ) / xc;
1952 y += t + lgamma_impl<double>::run(a +
b) - lgamma_impl<double>::run(a) -
1953 lgamma_impl<double>::run(
b);
1977 template <
typename Scalar>
1983 template <
typename Scalar>
1989 template <
typename Scalar>
1995 template <
typename Scalar>
2001 template <
typename Scalar>
2003 erf(
const Scalar&
x) {
2007 template <
typename Scalar>
2009 erfc(
const Scalar&
x) {
2013 template <
typename Scalar>
2019 template <
typename Scalar>
2021 igamma(
const Scalar& a,
const Scalar&
x) {
2025 template <
typename Scalar>
2031 template <
typename Scalar>
2037 template <
typename Scalar>
2039 igammac(
const Scalar& a,
const Scalar&
x) {
2043 template <
typename Scalar>
2045 betainc(const Scalar& a, const Scalar&
b, const Scalar&
x) {
Array< int, Dynamic, 1 > v
SparseMatrix< double > A(n, n)
#define EIGEN_ALWAYS_INLINE
#define EIGEN_DEVICE_FUNC
#define EIGEN_MATHFUNC_IMPL(func, scalar)
#define EIGEN_STATIC_ASSERT(X, MSG)
static Scalar main_igamma_term(Scalar a, Scalar x)
Packet pselect(const Packet &mask, const Packet &a, const Packet &b)
bool pmul(const bool &a, const bool &b)
T generic_fast_erf_float(const T &x)
EIGEN_DECLARE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet plog(const Packet &a)
Packet pabs(const Packet &a)
Packet pcmp_lt(const Packet &a, const Packet &b)
double flipsign< double >(const double &should_flipsign, const double &x)
T flipsign(const T &should_flipsign, const T &x)
Packet pcmp_le(const Packet &a, const Packet &b)
T generic_ndtri_gt_exp_neg_two(const T &b)
EIGEN_DECLARE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet psqrt(const Packet &a)
float flipsign< float >(const float &should_flipsign, const float &x)
Packet psub(const Packet &a, const Packet &b)
Packet pcmp_eq(const Packet &a, const Packet &b)
EIGEN_ALWAYS_INLINE T generic_ndtri(const T &a)
T generic_ndtri_lt_exp_neg_two(const T &b, const T &should_flipsign)
bool psign(const bool &a)
Packet pdiv(const Packet &a, const Packet &b)
Packet preciprocal(const Packet &a)
Packet pmadd(const Packet &a, const Packet &b, const Packet &c)
int igamma_num_iterations()
EIGEN_ALWAYS_INLINE T tan(const T &x)
EIGEN_ALWAYS_INLINE bool() isinf(const Eigen::bfloat16 &h)
EIGEN_MATHFUNC_RETVAL(abs2, Scalar) abs2(const Scalar &x)
internal::pow_impl< ScalarX, ScalarY >::result_type pow(const ScalarX &x, const ScalarY &y)
bool equal_strict(const double &x, const double &y)
EIGEN_ALWAYS_INLINE T exp(const T &x)
Scalar() floor(const Scalar &x)
EIGEN_ALWAYS_INLINE bool() isnan(const Eigen::bfloat16 &h)
EIGEN_ALWAYS_INLINE std::enable_if_t< NumTraits< T >::IsSigned||NumTraits< T >::IsComplex, typename NumTraits< T >::Real > abs(const T &x)
EIGEN_ALWAYS_INLINE T log(const T &x)
: TensorContractionSycl.h, provides various tensor contraction kernel for SYCL backend
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_erfc_op< typename Derived::Scalar >, const Derived > erfc(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseBinaryOp< Eigen::internal::scalar_igammac_op< typename Derived::Scalar >, const Derived, const ExponentDerived > igammac(const Eigen::ArrayBase< Derived > &a, const Eigen::ArrayBase< ExponentDerived > &x)
const TensorCwiseTernaryOp< internal::scalar_betainc_op< typename XDerived::Scalar >, const ADerived, const BDerived, const XDerived > betainc(const ADerived &a, const BDerived &b, const XDerived &x)
const Eigen::CwiseBinaryOp< Eigen::internal::scalar_igamma_op< typename Derived::Scalar >, const Derived, const ExponentDerived > igamma(const Eigen::ArrayBase< Derived > &a, const Eigen::ArrayBase< ExponentDerived > &x)
const Eigen::CwiseBinaryOp< Eigen::internal::scalar_polygamma_op< typename DerivedX::Scalar >, const DerivedN, const DerivedX > polygamma(const Eigen::ArrayBase< DerivedN > &n, const Eigen::ArrayBase< DerivedX > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_ndtri_op< typename Derived::Scalar >, const Derived > ndtri(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_lgamma_op< typename Derived::Scalar >, const Derived > lgamma(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_erf_op< typename Derived::Scalar >, const Derived > erf(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseBinaryOp< Eigen::internal::scalar_igamma_der_a_op< typename Derived::Scalar >, const Derived, const ExponentDerived > igamma_der_a(const Eigen::ArrayBase< Derived > &a, const Eigen::ArrayBase< ExponentDerived > &x)
const Eigen::CwiseBinaryOp< Eigen::internal::scalar_zeta_op< typename DerivedX::Scalar >, const DerivedX, const DerivedQ > zeta(const Eigen::ArrayBase< DerivedX > &x, const Eigen::ArrayBase< DerivedQ > &q)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_digamma_op< typename Derived::Scalar >, const Derived > digamma(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseBinaryOp< Eigen::internal::scalar_gamma_sample_der_alpha_op< typename AlphaDerived::Scalar >, const AlphaDerived, const SampleDerived > gamma_sample_der_alpha(const Eigen::ArrayBase< AlphaDerived > &alpha, const Eigen::ArrayBase< SampleDerived > &sample)