10 #ifndef EIGEN_SPECIAL_FUNCTIONS_H 11 #define EIGEN_SPECIAL_FUNCTIONS_H 80 template <
typename Scalar,
int N>
90 template <
typename Scalar>
104 template <
typename Scalar>
109 THIS_TYPE_IS_NOT_SUPPORTED);
114 template <
typename Scalar>
119 #if EIGEN_HAS_C99_MATH 124 #if !defined(__CUDA_ARCH__) && (defined(_BSD_SOURCE) || defined(_SVID_SOURCE)) && !defined(__APPLE__) 126 return ::lgammaf_r(x, &signgam);
137 #if !defined(__CUDA_ARCH__) && (defined(_BSD_SOURCE) || defined(_SVID_SOURCE)) && !defined(__APPLE__) 139 return ::lgamma_r(x, &signgam);
151 template <
typename Scalar>
169 template <
typename Scalar>
174 THIS_TYPE_IS_NOT_SUPPORTED);
185 -4.16666666666666666667E-3
f,
186 3.96825396825396825397E-3
f,
187 -8.33333333333333333333E-3
f,
188 8.33333333333333333333E-2
f 204 8.33333333333333333333E-2,
205 -2.10927960927960927961E-2,
206 7.57575757575757575758E-3,
207 -4.16666666666666666667E-3,
208 3.96825396825396825397E-3,
209 -8.33333333333333333333E-3,
210 8.33333333333333333333E-2
222 template <
typename Scalar>
225 static Scalar
run(Scalar x) {
283 Scalar p,
q, nz,
s,
w,
y;
284 bool negative =
false;
287 const Scalar m_pi = Scalar(
EIGEN_PI);
289 const Scalar zero = Scalar(0);
290 const Scalar one = Scalar(1);
291 const Scalar
half = Scalar(0.5);
321 while (s < Scalar(10)) {
330 return (negative) ? y - nz :
y;
338 template <
typename Scalar>
343 THIS_TYPE_IS_NOT_SUPPORTED);
348 template <
typename Scalar>
353 #if EIGEN_HAS_C99_MATH 365 #endif // EIGEN_HAS_C99_MATH 371 template <
typename Scalar>
376 THIS_TYPE_IS_NOT_SUPPORTED);
381 template <
typename Scalar>
386 #if EIGEN_HAS_C99_MATH 398 #endif // EIGEN_HAS_C99_MATH 404 template <
typename Scalar>
410 template <
typename Scalar>
455 #if !EIGEN_HAS_C99_MATH 457 template <
typename Scalar>
460 static Scalar
run(Scalar a, Scalar x) {
462 THIS_TYPE_IS_NOT_SUPPORTED);
471 template <
typename Scalar>
474 static Scalar
run(Scalar a, Scalar x) {
529 const Scalar zero = 0;
530 const Scalar one = 1;
533 if ((x < zero) || (a <= zero)) {
538 if ((x < one) || (x < a)) {
563 EIGEN_DEVICE_FUNC
static Scalar Impl(Scalar a, Scalar x) {
564 const Scalar zero = 0;
565 const Scalar one = 1;
566 const Scalar two = 2;
573 Scalar ans, ax, c, yc, r, t,
y,
z;
574 Scalar pk, pkm1, pkm2, qk, qkm1, qkm2;
576 if (x == inf)
return zero;
600 pk = pkm1 * z - pkm2 * yc;
601 qk = qkm1 * z - qkm2 * yc;
628 #endif // EIGEN_HAS_C99_MATH 634 template <
typename Scalar>
639 #if !EIGEN_HAS_C99_MATH 641 template <
typename Scalar>
646 THIS_TYPE_IS_NOT_SUPPORTED);
653 template <
typename Scalar>
656 static Scalar
run(Scalar a, Scalar x) {
717 const Scalar zero = 0;
718 const Scalar one = 1;
721 if (x == zero)
return zero;
723 if ((x < zero) || (a <= zero)) {
727 if ((x > one) && (x > a)) {
752 EIGEN_DEVICE_FUNC
static Scalar Impl(Scalar a, Scalar x) {
753 const Scalar zero = 0;
754 const Scalar one = 1;
758 Scalar ans, ax, c, r;
777 if (c/ans <= machep) {
782 return (ans * ax / a);
786 #endif // EIGEN_HAS_C99_MATH 792 template <
typename Scalar>
797 template <
typename Scalar>
802 THIS_TYPE_IS_NOT_SUPPORTED);
832 while( (i < 9) || (a <= 9.0) )
847 template <
typename Scalar>
850 static Scalar
run(Scalar x, Scalar
q) {
913 Scalar p, r, a,
b, k,
s, t,
w;
921 Scalar(-1.8924375803183791606e9),
922 Scalar(7.47242496e10),
923 Scalar(-2.950130727918164224e12),
924 Scalar(1.1646782814350067249e14),
925 Scalar(-4.5979787224074726105e15),
926 Scalar(1.8152105401943546773e17),
927 Scalar(-7.1661652561756670113e18)
931 const Scalar zero = 0.0,
half = 0.5, one = 1.0;
973 for( i=0; i<12; i++ )
996 template <
typename Scalar>
1001 #if !EIGEN_HAS_C99_MATH 1003 template <
typename Scalar>
1008 THIS_TYPE_IS_NOT_SUPPORTED);
1015 template <
typename Scalar>
1018 static Scalar
run(Scalar n, Scalar x) {
1019 Scalar zero = 0.0, one = 1.0;
1020 Scalar nplus = n + one;
1028 else if (n == zero) {
1039 #endif // EIGEN_HAS_C99_MATH 1045 template <
typename Scalar>
1050 #if !EIGEN_HAS_C99_MATH 1052 template <
typename Scalar>
1057 THIS_TYPE_IS_NOT_SUPPORTED);
1064 template <
typename Scalar>
1138 THIS_TYPE_IS_NOT_SUPPORTED);
1146 template <
typename Scalar>
1147 struct incbeta_cfe {
1152 THIS_TYPE_IS_NOT_SUPPORTED);
1157 const Scalar zero = 0;
1158 const Scalar one = 1;
1159 const Scalar two = 2;
1161 Scalar xk, pk, pkm1, pkm2, qk, qkm1, qkm2;
1162 Scalar k1, k2, k3, k4, k5, k6, k7, k8, k26update;
1167 const Scalar thresh =
1202 xk = -(x * k1 * k2) / (k3 * k4);
1203 pk = pkm1 + pkm2 * xk;
1204 qk = qkm1 + qkm2 * xk;
1210 xk = (x * k5 * k6) / (k7 * k8);
1211 pk = pkm1 + pkm2 * xk;
1212 qk = qkm1 + qkm2 * xk;
1247 }
while (++n < num_iters);
1254 template <
typename Scalar>
1255 struct betainc_helper {};
1258 struct betainc_helper<float> {
1262 float ans, a,
b, t, x, onemx;
1263 bool reversed_a_b =
false;
1268 if (xx > (aa / (aa + bb))) {
1269 reversed_a_b =
true;
1284 t = betainc_helper<float>::incbps(a, b, x);
1285 if (reversed_a_b) t = 1.0f - t;
1290 ans = x * (a + b - 2.0f) / (a - 1.0
f);
1304 if (reversed_a_b) t = 1.0f - t;
1337 static float run(
float a,
float b,
float x) {
1341 if (a <= 0.0
f)
return nan;
1342 if (b <= 0.0
f)
return nan;
1343 if ((x <= 0.0
f) || (x >= 1.0
f)) {
1344 if (x == 0.0
f)
return 0.0f;
1345 if (x == 1.0
f)
return 1.0f;
1352 ans = betainc_helper<float>::incbsa(a + 1.0
f, b, x);
1358 return betainc_helper<float>::incbsa(a, b, x);
1364 struct betainc_helper<double> {
1369 double s, t, u, v, n, t1,
z,
ai;
1380 u = (n -
b) * x / n;
1406 static double run(
double aa,
double bb,
double xx) {
1411 double a,
b, t, x, xc,
w,
y;
1412 bool reversed_a_b =
false;
1414 if (aa <= 0.0 || bb <= 0.0) {
1418 if ((xx <= 0.0) || (xx >= 1.0)) {
1419 if (xx == 0.0)
return (0.0);
1420 if (xx == 1.0)
return (1.0);
1425 if ((bb * xx) <= 1.0 && xx <= 0.95) {
1426 return betainc_helper<double>::incbps(aa, bb, xx);
1432 if (xx > (aa / (aa + bb))) {
1433 reversed_a_b =
true;
1445 if (reversed_a_b && (b * x) <= 1.0 && x <= 0.95) {
1446 t = betainc_helper<double>::incbps(a, b, x);
1456 y = x * (a + b - 2.0) - (a - 1.0);
1500 #endif // EIGEN_HAS_C99_MATH 1506 template <
typename Scalar>
1508 lgamma(const Scalar& x) {
1512 template <
typename Scalar>
1518 template <
typename Scalar>
1520 zeta(const Scalar& x, const Scalar&
q) {
1524 template <
typename Scalar>
1530 template <
typename Scalar>
1532 erf(const Scalar& x) {
1536 template <
typename Scalar>
1538 erfc(const Scalar& x) {
1542 template <
typename Scalar>
1544 igamma(
const Scalar& a,
const Scalar&
x) {
1548 template <
typename Scalar>
1550 igammac(
const Scalar& a,
const Scalar&
x) {
1554 template <
typename Scalar>
1556 betainc(const Scalar& a, const Scalar&
b, const Scalar& x) {
1565 #endif // EIGEN_SPECIAL_FUNCTIONS_H static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar biginv()
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar big()
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC half pow(const half &a, const half &b)
#define EIGEN_STRONG_INLINE
const mpreal ai(const mpreal &x, mp_rnd_t r=mpreal::get_default_rnd())
EIGEN_DEVICE_FUNC const ExpReturnType exp() const
EIGEN_DEVICE_FUNC const ErfReturnType erf() const
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar run(const Scalar)
static int f(const TensorMap< Tensor< int, 3 > > &tensor)
EIGEN_DEVICE_FUNC const Scalar & x
EIGEN_DEVICE_FUNC const LogReturnType log() const
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
#define EIGEN_STATIC_ASSERT(CONDITION, MSG)
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const AbsReturnType abs() const
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 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)
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE float machep()
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar run(Scalar a, Scalar x)
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar run(const Scalar x, const Scalar coef[])
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar run(Scalar a, Scalar b, Scalar x)
static EIGEN_DEVICE_FUNC Scalar run(Scalar x, Scalar q)
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE double biginv()
EIGEN_DEVICE_FUNC const Log1pReturnType log1p() const
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar run(const Scalar)
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar run(const Scalar)
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar run(Scalar n, Scalar x)
EIGEN_DEVICE_FUNC const Scalar & q
EIGEN_DEVICE_FUNC const TanReturnType tan() const
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE float biginv()
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar run(const Scalar, const Scalar coef[])
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE double machep()
EIGEN_DEVICE_FUNC const ErfcReturnType erfc() const
TFSIMD_FORCE_INLINE const tfScalar & z() const
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE float run(const float s)
static EIGEN_DEVICE_FUNC Scalar run(Scalar x)
TFSIMD_FORCE_INLINE const tfScalar & w() const
#define EIGEN_MATHFUNC_RETVAL(func, scalar)
EIGEN_DEVICE_FUNC const FloorReturnType floor() const
static EIGEN_DEVICE_FUNC Scalar run(Scalar a, Scalar x)
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE bool run(double &a, double &b, double &s, const double x, const double machep)
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)
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE float big()
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE bool run(float &a, float &b, float &s, const float x, const float machep)
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)
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE double run(const double s)
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE 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)
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE double big()
void run(Expr &expr, Dev &dev)
EIGEN_DEVICE_FUNC const Scalar & b
#define EIGEN_MATHFUNC_IMPL(func, scalar)
EIGEN_DEVICE_FUNC const DigammaReturnType digamma() const
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar run(const Scalar)
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar run(const Scalar)
EIGEN_DEVICE_FUNC const LgammaReturnType lgamma() const
static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar machep()