Eigen-unsupported  3.4.1 (git rev 28ded8800c26864e537852658428ab44c8399e87)
 
Loading...
Searching...
No Matches
SpecialFunctionsImpl.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2015 Eugene Brevdo <ebrevdo@gmail.com>
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
10#ifndef EIGEN_SPECIAL_FUNCTIONS_H
11#define EIGEN_SPECIAL_FUNCTIONS_H
12
13namespace Eigen {
14namespace internal {
15
16// Parts of this code are based on the Cephes Math Library.
17//
18// Cephes Math Library Release 2.8: June, 2000
19// Copyright 1984, 1987, 1992, 2000 by Stephen L. Moshier
20//
21// Permission has been kindly provided by the original author
22// to incorporate the Cephes software into the Eigen codebase:
23//
24// From: Stephen Moshier
25// To: Eugene Brevdo
26// Subject: Re: Permission to wrap several cephes functions in Eigen
27//
28// Hello Eugene,
29//
30// Thank you for writing.
31//
32// If your licensing is similar to BSD, the formal way that has been
33// handled is simply to add a statement to the effect that you are incorporating
34// the Cephes software by permission of the author.
35//
36// Good luck with your project,
37// Steve
38
39
40/****************************************************************************
41 * Implementation of lgamma, requires C++11/C99 *
42 ****************************************************************************/
43
44template <typename Scalar>
45struct lgamma_impl {
46 EIGEN_DEVICE_FUNC
47 static EIGEN_STRONG_INLINE Scalar run(const Scalar) {
48 EIGEN_STATIC_ASSERT((internal::is_same<Scalar, Scalar>::value == false),
49 THIS_TYPE_IS_NOT_SUPPORTED);
50 return Scalar(0);
51 }
52};
53
54template <typename Scalar>
55struct lgamma_retval {
56 typedef Scalar type;
57};
58
59#if EIGEN_HAS_C99_MATH
60// Since glibc 2.19
61#if defined(__GLIBC__) && ((__GLIBC__>=2 && __GLIBC_MINOR__ >= 19) || __GLIBC__>2) \
62 && (defined(_DEFAULT_SOURCE) || defined(_BSD_SOURCE) || defined(_SVID_SOURCE))
63#define EIGEN_HAS_LGAMMA_R
64#endif
65
66// Glibc versions before 2.19
67#if defined(__GLIBC__) && ((__GLIBC__==2 && __GLIBC_MINOR__ < 19) || __GLIBC__<2) \
68 && (defined(_BSD_SOURCE) || defined(_SVID_SOURCE))
69#define EIGEN_HAS_LGAMMA_R
70#endif
71
72template <>
73struct lgamma_impl<float> {
74 EIGEN_DEVICE_FUNC
75 static EIGEN_STRONG_INLINE float run(float x) {
76#if !defined(EIGEN_GPU_COMPILE_PHASE) && defined (EIGEN_HAS_LGAMMA_R) && !defined(__APPLE__)
77 int dummy;
78 return ::lgammaf_r(x, &dummy);
79#elif defined(SYCL_DEVICE_ONLY)
80 return cl::sycl::lgamma(x);
81#else
82 return ::lgammaf(x);
83#endif
84 }
85};
86
87template <>
88struct lgamma_impl<double> {
89 EIGEN_DEVICE_FUNC
90 static EIGEN_STRONG_INLINE double run(double x) {
91#if !defined(EIGEN_GPU_COMPILE_PHASE) && defined(EIGEN_HAS_LGAMMA_R) && !defined(__APPLE__)
92 int dummy;
93 return ::lgamma_r(x, &dummy);
94#elif defined(SYCL_DEVICE_ONLY)
95 return cl::sycl::lgamma(x);
96#else
97 return ::lgamma(x);
98#endif
99 }
100};
101
102#undef EIGEN_HAS_LGAMMA_R
103#endif
104
105/****************************************************************************
106 * Implementation of digamma (psi), based on Cephes *
107 ****************************************************************************/
108
109template <typename Scalar>
110struct digamma_retval {
111 typedef Scalar type;
112};
113
114/*
115 *
116 * Polynomial evaluation helper for the Psi (digamma) function.
117 *
118 * digamma_impl_maybe_poly::run(s) evaluates the asymptotic Psi expansion for
119 * input Scalar s, assuming s is above 10.0.
120 *
121 * If s is above a certain threshold for the given Scalar type, zero
122 * is returned. Otherwise the polynomial is evaluated with enough
123 * coefficients for results matching Scalar machine precision.
124 *
125 *
126 */
127template <typename Scalar>
128struct digamma_impl_maybe_poly {
129 EIGEN_DEVICE_FUNC
130 static EIGEN_STRONG_INLINE Scalar run(const Scalar) {
131 EIGEN_STATIC_ASSERT((internal::is_same<Scalar, Scalar>::value == false),
132 THIS_TYPE_IS_NOT_SUPPORTED);
133 return Scalar(0);
134 }
135};
136
137
138template <>
139struct digamma_impl_maybe_poly<float> {
140 EIGEN_DEVICE_FUNC
141 static EIGEN_STRONG_INLINE float run(const float s) {
142 const float A[] = {
143 -4.16666666666666666667E-3f,
144 3.96825396825396825397E-3f,
145 -8.33333333333333333333E-3f,
146 8.33333333333333333333E-2f
147 };
148
149 float z;
150 if (s < 1.0e8f) {
151 z = 1.0f / (s * s);
152 return z * internal::ppolevl<float, 3>::run(z, A);
153 } else return 0.0f;
154 }
155};
156
157template <>
158struct digamma_impl_maybe_poly<double> {
159 EIGEN_DEVICE_FUNC
160 static EIGEN_STRONG_INLINE double run(const double s) {
161 const double A[] = {
162 8.33333333333333333333E-2,
163 -2.10927960927960927961E-2,
164 7.57575757575757575758E-3,
165 -4.16666666666666666667E-3,
166 3.96825396825396825397E-3,
167 -8.33333333333333333333E-3,
168 8.33333333333333333333E-2
169 };
170
171 double z;
172 if (s < 1.0e17) {
173 z = 1.0 / (s * s);
174 return z * internal::ppolevl<double, 6>::run(z, A);
175 }
176 else return 0.0;
177 }
178};
179
180template <typename Scalar>
181struct digamma_impl {
182 EIGEN_DEVICE_FUNC
183 static Scalar run(Scalar x) {
184 /*
185 *
186 * Psi (digamma) function (modified for Eigen)
187 *
188 *
189 * SYNOPSIS:
190 *
191 * double x, y, psi();
192 *
193 * y = psi( x );
194 *
195 *
196 * DESCRIPTION:
197 *
198 * d -
199 * psi(x) = -- ln | (x)
200 * dx
201 *
202 * is the logarithmic derivative of the gamma function.
203 * For integer x,
204 * n-1
205 * -
206 * psi(n) = -EUL + > 1/k.
207 * -
208 * k=1
209 *
210 * If x is negative, it is transformed to a positive argument by the
211 * reflection formula psi(1-x) = psi(x) + pi cot(pi x).
212 * For general positive x, the argument is made greater than 10
213 * using the recurrence psi(x+1) = psi(x) + 1/x.
214 * Then the following asymptotic expansion is applied:
215 *
216 * inf. B
217 * - 2k
218 * psi(x) = log(x) - 1/2x - > -------
219 * - 2k
220 * k=1 2k x
221 *
222 * where the B2k are Bernoulli numbers.
223 *
224 * ACCURACY (float):
225 * Relative error (except absolute when |psi| < 1):
226 * arithmetic domain # trials peak rms
227 * IEEE 0,30 30000 1.3e-15 1.4e-16
228 * IEEE -30,0 40000 1.5e-15 2.2e-16
229 *
230 * ACCURACY (double):
231 * Absolute error, relative when |psi| > 1 :
232 * arithmetic domain # trials peak rms
233 * IEEE -33,0 30000 8.2e-7 1.2e-7
234 * IEEE 0,33 100000 7.3e-7 7.7e-8
235 *
236 * ERROR MESSAGES:
237 * message condition value returned
238 * psi singularity x integer <=0 INFINITY
239 */
240
241 Scalar p, q, nz, s, w, y;
242 bool negative = false;
243
244 const Scalar nan = NumTraits<Scalar>::quiet_NaN();
245 const Scalar m_pi = Scalar(EIGEN_PI);
246
247 const Scalar zero = Scalar(0);
248 const Scalar one = Scalar(1);
249 const Scalar half = Scalar(0.5);
250 nz = zero;
251
252 if (x <= zero) {
253 negative = true;
254 q = x;
255 p = numext::floor(q);
256 if (p == q) {
257 return nan;
258 }
259 /* Remove the zeros of tan(m_pi x)
260 * by subtracting the nearest integer from x
261 */
262 nz = q - p;
263 if (nz != half) {
264 if (nz > half) {
265 p += one;
266 nz = q - p;
267 }
268 nz = m_pi / numext::tan(m_pi * nz);
269 }
270 else {
271 nz = zero;
272 }
273 x = one - x;
274 }
275
276 /* use the recurrence psi(x+1) = psi(x) + 1/x. */
277 s = x;
278 w = zero;
279 while (s < Scalar(10)) {
280 w += one / s;
281 s += one;
282 }
283
284 y = digamma_impl_maybe_poly<Scalar>::run(s);
285
286 y = numext::log(s) - (half / s) - y - w;
287
288 return (negative) ? y - nz : y;
289 }
290};
291
292/****************************************************************************
293 * Implementation of erf, requires C++11/C99 *
294 ****************************************************************************/
295
303template <typename T>
304EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE T generic_fast_erf_float(const T& x) {
305 const float kErfInvOneMinusHalfULP = 3.832506856900711f;
306 const T clamp = pcmp_le(pset1<T>(kErfInvOneMinusHalfULP), pabs(x));
307 // The monomial coefficients of the numerator polynomial (odd).
308 const T alpha_1 = pset1<T>(-1.60960333262415e-02f);
309 const T alpha_3 = pset1<T>(-2.95459980854025e-03f);
310 const T alpha_5 = pset1<T>(-7.34990630326855e-04f);
311 const T alpha_7 = pset1<T>(-5.69250639462346e-05f);
312 const T alpha_9 = pset1<T>(-2.10102402082508e-06f);
313 const T alpha_11 = pset1<T>(2.77068142495902e-08f);
314 const T alpha_13 = pset1<T>(-2.72614225801306e-10f);
315
316 // The monomial coefficients of the denominator polynomial (even).
317 const T beta_0 = pset1<T>(-1.42647390514189e-02f);
318 const T beta_2 = pset1<T>(-7.37332916720468e-03f);
319 const T beta_4 = pset1<T>(-1.68282697438203e-03f);
320 const T beta_6 = pset1<T>(-2.13374055278905e-04f);
321 const T beta_8 = pset1<T>(-1.45660718464996e-05f);
322
323 // Since the polynomials are odd/even, we need x^2.
324 const T x2 = pmul(x, x);
325
326 // Evaluate the numerator polynomial p.
327 T p = pmadd(x2, alpha_13, alpha_11);
328 p = pmadd(x2, p, alpha_9);
329 p = pmadd(x2, p, alpha_7);
330 p = pmadd(x2, p, alpha_5);
331 p = pmadd(x2, p, alpha_3);
332 p = pmadd(x2, p, alpha_1);
333 p = pmul(x, p);
334
335 // Evaluate the denominator polynomial p.
336 T q = pmadd(x2, beta_8, beta_6);
337 q = pmadd(x2, q, beta_4);
338 q = pmadd(x2, q, beta_2);
339 q = pmadd(x2, q, beta_0);
340
341 // Divide the numerator by the denominator.
342 const T sign = pselect(pcmp_le(x, pset1<T>(0.0f)), pset1<T>(-1.0f), pset1<T>(1.0f));
343 return pselect(clamp, sign, pdiv(p, q));
344}
345
346template <typename T>
347struct erf_impl {
348 EIGEN_DEVICE_FUNC
349 static EIGEN_STRONG_INLINE T run(const T& x) {
350 return generic_fast_erf_float(x);
351 }
352};
353
354template <typename Scalar>
355struct erf_retval {
356 typedef Scalar type;
357};
358
359#if EIGEN_HAS_C99_MATH
360template <>
361struct erf_impl<float> {
362 EIGEN_DEVICE_FUNC
363 static EIGEN_STRONG_INLINE float run(float x) {
364#if defined(SYCL_DEVICE_ONLY)
365 return cl::sycl::erf(x);
366#else
367 return generic_fast_erf_float(x);
368#endif
369 }
370};
371
372template <>
373struct erf_impl<double> {
374 EIGEN_DEVICE_FUNC
375 static EIGEN_STRONG_INLINE double run(double x) {
376#if defined(SYCL_DEVICE_ONLY)
377 return cl::sycl::erf(x);
378#else
379 return ::erf(x);
380#endif
381 }
382};
383#endif // EIGEN_HAS_C99_MATH
384
385/***************************************************************************
386* Implementation of erfc, requires C++11/C99 *
387****************************************************************************/
388
389template <typename Scalar>
390struct erfc_impl {
391 EIGEN_DEVICE_FUNC
392 static EIGEN_STRONG_INLINE Scalar run(const Scalar) {
393 EIGEN_STATIC_ASSERT((internal::is_same<Scalar, Scalar>::value == false),
394 THIS_TYPE_IS_NOT_SUPPORTED);
395 return Scalar(0);
396 }
397};
398
399template <typename Scalar>
400struct erfc_retval {
401 typedef Scalar type;
402};
403
404#if EIGEN_HAS_C99_MATH
405template <>
406struct erfc_impl<float> {
407 EIGEN_DEVICE_FUNC
408 static EIGEN_STRONG_INLINE float run(const float x) {
409#if defined(SYCL_DEVICE_ONLY)
410 return cl::sycl::erfc(x);
411#else
412 return ::erfcf(x);
413#endif
414 }
415};
416
417template <>
418struct erfc_impl<double> {
419 EIGEN_DEVICE_FUNC
420 static EIGEN_STRONG_INLINE double run(const double x) {
421#if defined(SYCL_DEVICE_ONLY)
422 return cl::sycl::erfc(x);
423#else
424 return ::erfc(x);
425#endif
426 }
427};
428#endif // EIGEN_HAS_C99_MATH
429
430
431/***************************************************************************
432* Implementation of ndtri. *
433****************************************************************************/
434
435/* Inverse of Normal distribution function (modified for Eigen).
436 *
437 *
438 * SYNOPSIS:
439 *
440 * double x, y, ndtri();
441 *
442 * x = ndtri( y );
443 *
444 *
445 *
446 * DESCRIPTION:
447 *
448 * Returns the argument, x, for which the area under the
449 * Gaussian probability density function (integrated from
450 * minus infinity to x) is equal to y.
451 *
452 *
453 * For small arguments 0 < y < exp(-2), the program computes
454 * z = sqrt( -2.0 * log(y) ); then the approximation is
455 * x = z - log(z)/z - (1/z) P(1/z) / Q(1/z).
456 * There are two rational functions P/Q, one for 0 < y < exp(-32)
457 * and the other for y up to exp(-2). For larger arguments,
458 * w = y - 0.5, and x/sqrt(2pi) = w + w**3 R(w**2)/S(w**2)).
459 *
460 *
461 * ACCURACY:
462 *
463 * Relative error:
464 * arithmetic domain # trials peak rms
465 * DEC 0.125, 1 5500 9.5e-17 2.1e-17
466 * DEC 6e-39, 0.135 3500 5.7e-17 1.3e-17
467 * IEEE 0.125, 1 20000 7.2e-16 1.3e-16
468 * IEEE 3e-308, 0.135 50000 4.6e-16 9.8e-17
469 *
470 *
471 * ERROR MESSAGES:
472 *
473 * message condition value returned
474 * ndtri domain x == 0 -INF
475 * ndtri domain x == 1 INF
476 * ndtri domain x < 0, x > 1 NAN
477 */
478 /*
479 Cephes Math Library Release 2.2: June, 1992
480 Copyright 1985, 1987, 1992 by Stephen L. Moshier
481 Direct inquiries to 30 Frost Street, Cambridge, MA 02140
482 */
483
484
485// TODO: Add a cheaper approximation for float.
486
487
488template<typename T>
489EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE T flipsign(
490 const T& should_flipsign, const T& x) {
491 typedef typename unpacket_traits<T>::type Scalar;
492 const T sign_mask = pset1<T>(Scalar(-0.0));
493 T sign_bit = pand<T>(should_flipsign, sign_mask);
494 return pxor<T>(sign_bit, x);
495}
496
497template<>
498EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE double flipsign<double>(
499 const double& should_flipsign, const double& x) {
500 return should_flipsign == 0 ? x : -x;
501}
502
503template<>
504EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE float flipsign<float>(
505 const float& should_flipsign, const float& x) {
506 return should_flipsign == 0 ? x : -x;
507}
508
509// We split this computation in to two so that in the scalar path
510// only one branch is evaluated (due to our template specialization of pselect
511// being an if statement.)
512
513template <typename T, typename ScalarType>
514EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE T generic_ndtri_gt_exp_neg_two(const T& b) {
515 const ScalarType p0[] = {
516 ScalarType(-5.99633501014107895267e1),
517 ScalarType(9.80010754185999661536e1),
518 ScalarType(-5.66762857469070293439e1),
519 ScalarType(1.39312609387279679503e1),
520 ScalarType(-1.23916583867381258016e0)
521 };
522 const ScalarType q0[] = {
523 ScalarType(1.0),
524 ScalarType(1.95448858338141759834e0),
525 ScalarType(4.67627912898881538453e0),
526 ScalarType(8.63602421390890590575e1),
527 ScalarType(-2.25462687854119370527e2),
528 ScalarType(2.00260212380060660359e2),
529 ScalarType(-8.20372256168333339912e1),
530 ScalarType(1.59056225126211695515e1),
531 ScalarType(-1.18331621121330003142e0)
532 };
533 const T sqrt2pi = pset1<T>(ScalarType(2.50662827463100050242e0));
534 const T half = pset1<T>(ScalarType(0.5));
535 T c, c2, ndtri_gt_exp_neg_two;
536
537 c = psub(b, half);
538 c2 = pmul(c, c);
539 ndtri_gt_exp_neg_two = pmadd(c, pmul(
540 c2, pdiv(
541 internal::ppolevl<T, 4>::run(c2, p0),
542 internal::ppolevl<T, 8>::run(c2, q0))), c);
543 return pmul(ndtri_gt_exp_neg_two, sqrt2pi);
544}
545
546template <typename T, typename ScalarType>
547EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE T generic_ndtri_lt_exp_neg_two(
548 const T& b, const T& should_flipsign) {
549 /* Approximation for interval z = sqrt(-2 log a ) between 2 and 8
550 * i.e., a between exp(-2) = .135 and exp(-32) = 1.27e-14.
551 */
552 const ScalarType p1[] = {
553 ScalarType(4.05544892305962419923e0),
554 ScalarType(3.15251094599893866154e1),
555 ScalarType(5.71628192246421288162e1),
556 ScalarType(4.40805073893200834700e1),
557 ScalarType(1.46849561928858024014e1),
558 ScalarType(2.18663306850790267539e0),
559 ScalarType(-1.40256079171354495875e-1),
560 ScalarType(-3.50424626827848203418e-2),
561 ScalarType(-8.57456785154685413611e-4)
562 };
563 const ScalarType q1[] = {
564 ScalarType(1.0),
565 ScalarType(1.57799883256466749731e1),
566 ScalarType(4.53907635128879210584e1),
567 ScalarType(4.13172038254672030440e1),
568 ScalarType(1.50425385692907503408e1),
569 ScalarType(2.50464946208309415979e0),
570 ScalarType(-1.42182922854787788574e-1),
571 ScalarType(-3.80806407691578277194e-2),
572 ScalarType(-9.33259480895457427372e-4)
573 };
574 /* Approximation for interval z = sqrt(-2 log a ) between 8 and 64
575 * i.e., a between exp(-32) = 1.27e-14 and exp(-2048) = 3.67e-890.
576 */
577 const ScalarType p2[] = {
578 ScalarType(3.23774891776946035970e0),
579 ScalarType(6.91522889068984211695e0),
580 ScalarType(3.93881025292474443415e0),
581 ScalarType(1.33303460815807542389e0),
582 ScalarType(2.01485389549179081538e-1),
583 ScalarType(1.23716634817820021358e-2),
584 ScalarType(3.01581553508235416007e-4),
585 ScalarType(2.65806974686737550832e-6),
586 ScalarType(6.23974539184983293730e-9)
587 };
588 const ScalarType q2[] = {
589 ScalarType(1.0),
590 ScalarType(6.02427039364742014255e0),
591 ScalarType(3.67983563856160859403e0),
592 ScalarType(1.37702099489081330271e0),
593 ScalarType(2.16236993594496635890e-1),
594 ScalarType(1.34204006088543189037e-2),
595 ScalarType(3.28014464682127739104e-4),
596 ScalarType(2.89247864745380683936e-6),
597 ScalarType(6.79019408009981274425e-9)
598 };
599 const T eight = pset1<T>(ScalarType(8.0));
600 const T one = pset1<T>(ScalarType(1));
601 const T neg_two = pset1<T>(ScalarType(-2));
602 T x, x0, x1, z;
603
604 x = psqrt(pmul(neg_two, plog(b)));
605 x0 = psub(x, pdiv(plog(x), x));
606 z = pdiv(one, x);
607 x1 = pmul(
608 z, pselect(
609 pcmp_lt(x, eight),
610 pdiv(internal::ppolevl<T, 8>::run(z, p1),
611 internal::ppolevl<T, 8>::run(z, q1)),
612 pdiv(internal::ppolevl<T, 8>::run(z, p2),
613 internal::ppolevl<T, 8>::run(z, q2))));
614 return flipsign(should_flipsign, psub(x0, x1));
615}
616
617template <typename T, typename ScalarType>
618EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE
619T generic_ndtri(const T& a) {
620 const T maxnum = pset1<T>(NumTraits<ScalarType>::infinity());
621 const T neg_maxnum = pset1<T>(-NumTraits<ScalarType>::infinity());
622
623 const T zero = pset1<T>(ScalarType(0));
624 const T one = pset1<T>(ScalarType(1));
625 // exp(-2)
626 const T exp_neg_two = pset1<T>(ScalarType(0.13533528323661269189));
627 T b, ndtri, should_flipsign;
628
629 should_flipsign = pcmp_le(a, psub(one, exp_neg_two));
630 b = pselect(should_flipsign, a, psub(one, a));
631
632 ndtri = pselect(
633 pcmp_lt(exp_neg_two, b),
634 generic_ndtri_gt_exp_neg_two<T, ScalarType>(b),
635 generic_ndtri_lt_exp_neg_two<T, ScalarType>(b, should_flipsign));
636
637 return pselect(
638 pcmp_eq(a, zero), neg_maxnum,
639 pselect(pcmp_eq(one, a), maxnum, ndtri));
640}
641
642template <typename Scalar>
643struct ndtri_retval {
644 typedef Scalar type;
645};
646
647#if !EIGEN_HAS_C99_MATH
648
649template <typename Scalar>
650struct ndtri_impl {
651 EIGEN_DEVICE_FUNC
652 static EIGEN_STRONG_INLINE Scalar run(const Scalar) {
653 EIGEN_STATIC_ASSERT((internal::is_same<Scalar, Scalar>::value == false),
654 THIS_TYPE_IS_NOT_SUPPORTED);
655 return Scalar(0);
656 }
657};
658
659# else
660
661template <typename Scalar>
662struct ndtri_impl {
663 EIGEN_DEVICE_FUNC
664 static EIGEN_STRONG_INLINE Scalar run(const Scalar x) {
665 return generic_ndtri<Scalar, Scalar>(x);
666 }
667};
668
669#endif // EIGEN_HAS_C99_MATH
670
671
672/**************************************************************************************************************
673 * Implementation of igammac (complemented incomplete gamma integral), based on Cephes but requires C++11/C99 *
674 **************************************************************************************************************/
675
676template <typename Scalar>
677struct igammac_retval {
678 typedef Scalar type;
679};
680
681// NOTE: cephes_helper is also used to implement zeta
682template <typename Scalar>
683struct cephes_helper {
684 EIGEN_DEVICE_FUNC
685 static EIGEN_STRONG_INLINE Scalar machep() { assert(false && "machep not supported for this type"); return 0.0; }
686 EIGEN_DEVICE_FUNC
687 static EIGEN_STRONG_INLINE Scalar big() { assert(false && "big not supported for this type"); return 0.0; }
688 EIGEN_DEVICE_FUNC
689 static EIGEN_STRONG_INLINE Scalar biginv() { assert(false && "biginv not supported for this type"); return 0.0; }
690};
691
692template <>
693struct cephes_helper<float> {
694 EIGEN_DEVICE_FUNC
695 static EIGEN_STRONG_INLINE float machep() {
696 return NumTraits<float>::epsilon() / 2; // 1.0 - machep == 1.0
697 }
698 EIGEN_DEVICE_FUNC
699 static EIGEN_STRONG_INLINE float big() {
700 // use epsneg (1.0 - epsneg == 1.0)
701 return 1.0f / (NumTraits<float>::epsilon() / 2);
702 }
703 EIGEN_DEVICE_FUNC
704 static EIGEN_STRONG_INLINE float biginv() {
705 // epsneg
706 return machep();
707 }
708};
709
710template <>
711struct cephes_helper<double> {
712 EIGEN_DEVICE_FUNC
713 static EIGEN_STRONG_INLINE double machep() {
714 return NumTraits<double>::epsilon() / 2; // 1.0 - machep == 1.0
715 }
716 EIGEN_DEVICE_FUNC
717 static EIGEN_STRONG_INLINE double big() {
718 return 1.0 / NumTraits<double>::epsilon();
719 }
720 EIGEN_DEVICE_FUNC
721 static EIGEN_STRONG_INLINE double biginv() {
722 // inverse of eps
723 return NumTraits<double>::epsilon();
724 }
725};
726
727enum IgammaComputationMode { VALUE, DERIVATIVE, SAMPLE_DERIVATIVE };
728
729template <typename Scalar>
730EIGEN_DEVICE_FUNC
731static EIGEN_STRONG_INLINE Scalar main_igamma_term(Scalar a, Scalar x) {
732 /* Compute x**a * exp(-x) / gamma(a) */
733 Scalar logax = a * numext::log(x) - x - lgamma_impl<Scalar>::run(a);
734 if (logax < -numext::log(NumTraits<Scalar>::highest()) ||
735 // Assuming x and a aren't Nan.
736 (numext::isnan)(logax)) {
737 return Scalar(0);
738 }
739 return numext::exp(logax);
740}
741
742template <typename Scalar, IgammaComputationMode mode>
743EIGEN_DEVICE_FUNC
744int igamma_num_iterations() {
745 /* Returns the maximum number of internal iterations for igamma computation.
746 */
747 if (mode == VALUE) {
748 return 2000;
749 }
750
751 if (internal::is_same<Scalar, float>::value) {
752 return 200;
753 } else if (internal::is_same<Scalar, double>::value) {
754 return 500;
755 } else {
756 return 2000;
757 }
758}
759
760template <typename Scalar, IgammaComputationMode mode>
761struct igammac_cf_impl {
762 /* Computes igamc(a, x) or derivative (depending on the mode)
763 * using the continued fraction expansion of the complementary
764 * incomplete Gamma function.
765 *
766 * Preconditions:
767 * a > 0
768 * x >= 1
769 * x >= a
770 */
771 EIGEN_DEVICE_FUNC
772 static Scalar run(Scalar a, Scalar x) {
773 const Scalar zero = 0;
774 const Scalar one = 1;
775 const Scalar two = 2;
776 const Scalar machep = cephes_helper<Scalar>::machep();
777 const Scalar big = cephes_helper<Scalar>::big();
778 const Scalar biginv = cephes_helper<Scalar>::biginv();
779
780 if ((numext::isinf)(x)) {
781 return zero;
782 }
783
784 Scalar ax = main_igamma_term<Scalar>(a, x);
785 // This is independent of mode. If this value is zero,
786 // then the function value is zero. If the function value is zero,
787 // then we are in a neighborhood where the function value evalutes to zero,
788 // so the derivative is zero.
789 if (ax == zero) {
790 return zero;
791 }
792
793 // continued fraction
794 Scalar y = one - a;
795 Scalar z = x + y + one;
796 Scalar c = zero;
797 Scalar pkm2 = one;
798 Scalar qkm2 = x;
799 Scalar pkm1 = x + one;
800 Scalar qkm1 = z * x;
801 Scalar ans = pkm1 / qkm1;
802
803 Scalar dpkm2_da = zero;
804 Scalar dqkm2_da = zero;
805 Scalar dpkm1_da = zero;
806 Scalar dqkm1_da = -x;
807 Scalar dans_da = (dpkm1_da - ans * dqkm1_da) / qkm1;
808
809 for (int i = 0; i < igamma_num_iterations<Scalar, mode>(); i++) {
810 c += one;
811 y += one;
812 z += two;
813
814 Scalar yc = y * c;
815 Scalar pk = pkm1 * z - pkm2 * yc;
816 Scalar qk = qkm1 * z - qkm2 * yc;
817
818 Scalar dpk_da = dpkm1_da * z - pkm1 - dpkm2_da * yc + pkm2 * c;
819 Scalar dqk_da = dqkm1_da * z - qkm1 - dqkm2_da * yc + qkm2 * c;
820
821 if (qk != zero) {
822 Scalar ans_prev = ans;
823 ans = pk / qk;
824
825 Scalar dans_da_prev = dans_da;
826 dans_da = (dpk_da - ans * dqk_da) / qk;
827
828 if (mode == VALUE) {
829 if (numext::abs(ans_prev - ans) <= machep * numext::abs(ans)) {
830 break;
831 }
832 } else {
833 if (numext::abs(dans_da - dans_da_prev) <= machep) {
834 break;
835 }
836 }
837 }
838
839 pkm2 = pkm1;
840 pkm1 = pk;
841 qkm2 = qkm1;
842 qkm1 = qk;
843
844 dpkm2_da = dpkm1_da;
845 dpkm1_da = dpk_da;
846 dqkm2_da = dqkm1_da;
847 dqkm1_da = dqk_da;
848
849 if (numext::abs(pk) > big) {
850 pkm2 *= biginv;
851 pkm1 *= biginv;
852 qkm2 *= biginv;
853 qkm1 *= biginv;
854
855 dpkm2_da *= biginv;
856 dpkm1_da *= biginv;
857 dqkm2_da *= biginv;
858 dqkm1_da *= biginv;
859 }
860 }
861
862 /* Compute x**a * exp(-x) / gamma(a) */
863 Scalar dlogax_da = numext::log(x) - digamma_impl<Scalar>::run(a);
864 Scalar dax_da = ax * dlogax_da;
865
866 switch (mode) {
867 case VALUE:
868 return ans * ax;
869 case DERIVATIVE:
870 return ans * dax_da + dans_da * ax;
871 case SAMPLE_DERIVATIVE:
872 default: // this is needed to suppress clang warning
873 return -(dans_da + ans * dlogax_da) * x;
874 }
875 }
876};
877
878template <typename Scalar, IgammaComputationMode mode>
879struct igamma_series_impl {
880 /* Computes igam(a, x) or its derivative (depending on the mode)
881 * using the series expansion of the incomplete Gamma function.
882 *
883 * Preconditions:
884 * x > 0
885 * a > 0
886 * !(x > 1 && x > a)
887 */
888 EIGEN_DEVICE_FUNC
889 static Scalar run(Scalar a, Scalar x) {
890 const Scalar zero = 0;
891 const Scalar one = 1;
892 const Scalar machep = cephes_helper<Scalar>::machep();
893
894 Scalar ax = main_igamma_term<Scalar>(a, x);
895
896 // This is independent of mode. If this value is zero,
897 // then the function value is zero. If the function value is zero,
898 // then we are in a neighborhood where the function value evalutes to zero,
899 // so the derivative is zero.
900 if (ax == zero) {
901 return zero;
902 }
903
904 ax /= a;
905
906 /* power series */
907 Scalar r = a;
908 Scalar c = one;
909 Scalar ans = one;
910
911 Scalar dc_da = zero;
912 Scalar dans_da = zero;
913
914 for (int i = 0; i < igamma_num_iterations<Scalar, mode>(); i++) {
915 r += one;
916 Scalar term = x / r;
917 Scalar dterm_da = -x / (r * r);
918 dc_da = term * dc_da + dterm_da * c;
919 dans_da += dc_da;
920 c *= term;
921 ans += c;
922
923 if (mode == VALUE) {
924 if (c <= machep * ans) {
925 break;
926 }
927 } else {
928 if (numext::abs(dc_da) <= machep * numext::abs(dans_da)) {
929 break;
930 }
931 }
932 }
933
934 Scalar dlogax_da = numext::log(x) - digamma_impl<Scalar>::run(a + one);
935 Scalar dax_da = ax * dlogax_da;
936
937 switch (mode) {
938 case VALUE:
939 return ans * ax;
940 case DERIVATIVE:
941 return ans * dax_da + dans_da * ax;
942 case SAMPLE_DERIVATIVE:
943 default: // this is needed to suppress clang warning
944 return -(dans_da + ans * dlogax_da) * x / a;
945 }
946 }
947};
948
949#if !EIGEN_HAS_C99_MATH
950
951template <typename Scalar>
952struct igammac_impl {
953 EIGEN_DEVICE_FUNC
954 static Scalar run(Scalar a, Scalar x) {
955 EIGEN_STATIC_ASSERT((internal::is_same<Scalar, Scalar>::value == false),
956 THIS_TYPE_IS_NOT_SUPPORTED);
957 return Scalar(0);
958 }
959};
960
961#else
962
963template <typename Scalar>
964struct igammac_impl {
965 EIGEN_DEVICE_FUNC
966 static Scalar run(Scalar a, Scalar x) {
967 /* igamc()
968 *
969 * Incomplete gamma integral (modified for Eigen)
970 *
971 *
972 *
973 * SYNOPSIS:
974 *
975 * double a, x, y, igamc();
976 *
977 * y = igamc( a, x );
978 *
979 * DESCRIPTION:
980 *
981 * The function is defined by
982 *
983 *
984 * igamc(a,x) = 1 - igam(a,x)
985 *
986 * inf.
987 * -
988 * 1 | | -t a-1
989 * = ----- | e t dt.
990 * - | |
991 * | (a) -
992 * x
993 *
994 *
995 * In this implementation both arguments must be positive.
996 * The integral is evaluated by either a power series or
997 * continued fraction expansion, depending on the relative
998 * values of a and x.
999 *
1000 * ACCURACY (float):
1001 *
1002 * Relative error:
1003 * arithmetic domain # trials peak rms
1004 * IEEE 0,30 30000 7.8e-6 5.9e-7
1005 *
1006 *
1007 * ACCURACY (double):
1008 *
1009 * Tested at random a, x.
1010 * a x Relative error:
1011 * arithmetic domain domain # trials peak rms
1012 * IEEE 0.5,100 0,100 200000 1.9e-14 1.7e-15
1013 * IEEE 0.01,0.5 0,100 200000 1.4e-13 1.6e-15
1014 *
1015 */
1016 /*
1017 Cephes Math Library Release 2.2: June, 1992
1018 Copyright 1985, 1987, 1992 by Stephen L. Moshier
1019 Direct inquiries to 30 Frost Street, Cambridge, MA 02140
1020 */
1021 const Scalar zero = 0;
1022 const Scalar one = 1;
1023 const Scalar nan = NumTraits<Scalar>::quiet_NaN();
1024
1025 if ((x < zero) || (a <= zero)) {
1026 // domain error
1027 return nan;
1028 }
1029
1030 if ((numext::isnan)(a) || (numext::isnan)(x)) { // propagate nans
1031 return nan;
1032 }
1033
1034 if ((x < one) || (x < a)) {
1035 return (one - igamma_series_impl<Scalar, VALUE>::run(a, x));
1036 }
1037
1038 return igammac_cf_impl<Scalar, VALUE>::run(a, x);
1039 }
1040};
1041
1042#endif // EIGEN_HAS_C99_MATH
1043
1044/************************************************************************************************
1045 * Implementation of igamma (incomplete gamma integral), based on Cephes but requires C++11/C99 *
1046 ************************************************************************************************/
1047
1048#if !EIGEN_HAS_C99_MATH
1049
1050template <typename Scalar, IgammaComputationMode mode>
1051struct igamma_generic_impl {
1052 EIGEN_DEVICE_FUNC
1053 static EIGEN_STRONG_INLINE Scalar run(Scalar a, Scalar x) {
1054 EIGEN_STATIC_ASSERT((internal::is_same<Scalar, Scalar>::value == false),
1055 THIS_TYPE_IS_NOT_SUPPORTED);
1056 return Scalar(0);
1057 }
1058};
1059
1060#else
1061
1062template <typename Scalar, IgammaComputationMode mode>
1063struct igamma_generic_impl {
1064 EIGEN_DEVICE_FUNC
1065 static Scalar run(Scalar a, Scalar x) {
1066 /* Depending on the mode, returns
1067 * - VALUE: incomplete Gamma function igamma(a, x)
1068 * - DERIVATIVE: derivative of incomplete Gamma function d/da igamma(a, x)
1069 * - SAMPLE_DERIVATIVE: implicit derivative of a Gamma random variable
1070 * x ~ Gamma(x | a, 1), dx/da = -1 / Gamma(x | a, 1) * d igamma(a, x) / dx
1071 *
1072 * Derivatives are implemented by forward-mode differentiation.
1073 */
1074 const Scalar zero = 0;
1075 const Scalar one = 1;
1076 const Scalar nan = NumTraits<Scalar>::quiet_NaN();
1077
1078 if (x == zero) return zero;
1079
1080 if ((x < zero) || (a <= zero)) { // domain error
1081 return nan;
1082 }
1083
1084 if ((numext::isnan)(a) || (numext::isnan)(x)) { // propagate nans
1085 return nan;
1086 }
1087
1088 if ((x > one) && (x > a)) {
1089 Scalar ret = igammac_cf_impl<Scalar, mode>::run(a, x);
1090 if (mode == VALUE) {
1091 return one - ret;
1092 } else {
1093 return -ret;
1094 }
1095 }
1096
1097 return igamma_series_impl<Scalar, mode>::run(a, x);
1098 }
1099};
1100
1101#endif // EIGEN_HAS_C99_MATH
1102
1103template <typename Scalar>
1104struct igamma_retval {
1105 typedef Scalar type;
1106};
1107
1108template <typename Scalar>
1109struct igamma_impl : igamma_generic_impl<Scalar, VALUE> {
1110 /* igam()
1111 * Incomplete gamma integral.
1112 *
1113 * The CDF of Gamma(a, 1) random variable at the point x.
1114 *
1115 * Accuracy estimation. For each a in [10^-2, 10^-1...10^3] we sample
1116 * 50 Gamma random variables x ~ Gamma(x | a, 1), a total of 300 points.
1117 * The ground truth is computed by mpmath. Mean absolute error:
1118 * float: 1.26713e-05
1119 * double: 2.33606e-12
1120 *
1121 * Cephes documentation below.
1122 *
1123 * SYNOPSIS:
1124 *
1125 * double a, x, y, igam();
1126 *
1127 * y = igam( a, x );
1128 *
1129 * DESCRIPTION:
1130 *
1131 * The function is defined by
1132 *
1133 * x
1134 * -
1135 * 1 | | -t a-1
1136 * igam(a,x) = ----- | e t dt.
1137 * - | |
1138 * | (a) -
1139 * 0
1140 *
1141 *
1142 * In this implementation both arguments must be positive.
1143 * The integral is evaluated by either a power series or
1144 * continued fraction expansion, depending on the relative
1145 * values of a and x.
1146 *
1147 * ACCURACY (double):
1148 *
1149 * Relative error:
1150 * arithmetic domain # trials peak rms
1151 * IEEE 0,30 200000 3.6e-14 2.9e-15
1152 * IEEE 0,100 300000 9.9e-14 1.5e-14
1153 *
1154 *
1155 * ACCURACY (float):
1156 *
1157 * Relative error:
1158 * arithmetic domain # trials peak rms
1159 * IEEE 0,30 20000 7.8e-6 5.9e-7
1160 *
1161 */
1162 /*
1163 Cephes Math Library Release 2.2: June, 1992
1164 Copyright 1985, 1987, 1992 by Stephen L. Moshier
1165 Direct inquiries to 30 Frost Street, Cambridge, MA 02140
1166 */
1167
1168 /* left tail of incomplete gamma function:
1169 *
1170 * inf. k
1171 * a -x - x
1172 * x e > ----------
1173 * - -
1174 * k=0 | (a+k+1)
1175 *
1176 */
1177};
1178
1179template <typename Scalar>
1180struct igamma_der_a_retval : igamma_retval<Scalar> {};
1181
1182template <typename Scalar>
1183struct igamma_der_a_impl : igamma_generic_impl<Scalar, DERIVATIVE> {
1184 /* Derivative of the incomplete Gamma function with respect to a.
1185 *
1186 * Computes d/da igamma(a, x) by forward differentiation of the igamma code.
1187 *
1188 * Accuracy estimation. For each a in [10^-2, 10^-1...10^3] we sample
1189 * 50 Gamma random variables x ~ Gamma(x | a, 1), a total of 300 points.
1190 * The ground truth is computed by mpmath. Mean absolute error:
1191 * float: 6.17992e-07
1192 * double: 4.60453e-12
1193 *
1194 * Reference:
1195 * R. Moore. "Algorithm AS 187: Derivatives of the incomplete gamma
1196 * integral". Journal of the Royal Statistical Society. 1982
1197 */
1198};
1199
1200template <typename Scalar>
1201struct gamma_sample_der_alpha_retval : igamma_retval<Scalar> {};
1202
1203template <typename Scalar>
1204struct gamma_sample_der_alpha_impl
1205 : igamma_generic_impl<Scalar, SAMPLE_DERIVATIVE> {
1206 /* Derivative of a Gamma random variable sample with respect to alpha.
1207 *
1208 * Consider a sample of a Gamma random variable with the concentration
1209 * parameter alpha: sample ~ Gamma(alpha, 1). The reparameterization
1210 * derivative that we want to compute is dsample / dalpha =
1211 * d igammainv(alpha, u) / dalpha, where u = igamma(alpha, sample).
1212 * However, this formula is numerically unstable and expensive, so instead
1213 * we use implicit differentiation:
1214 *
1215 * igamma(alpha, sample) = u, where u ~ Uniform(0, 1).
1216 * Apply d / dalpha to both sides:
1217 * d igamma(alpha, sample) / dalpha
1218 * + d igamma(alpha, sample) / dsample * dsample/dalpha = 0
1219 * d igamma(alpha, sample) / dalpha
1220 * + Gamma(sample | alpha, 1) dsample / dalpha = 0
1221 * dsample/dalpha = - (d igamma(alpha, sample) / dalpha)
1222 * / Gamma(sample | alpha, 1)
1223 *
1224 * Here Gamma(sample | alpha, 1) is the PDF of the Gamma distribution
1225 * (note that the derivative of the CDF w.r.t. sample is the PDF).
1226 * See the reference below for more details.
1227 *
1228 * The derivative of igamma(alpha, sample) is computed by forward
1229 * differentiation of the igamma code. Division by the Gamma PDF is performed
1230 * in the same code, increasing the accuracy and speed due to cancellation
1231 * of some terms.
1232 *
1233 * Accuracy estimation. For each alpha in [10^-2, 10^-1...10^3] we sample
1234 * 50 Gamma random variables sample ~ Gamma(sample | alpha, 1), a total of 300
1235 * points. The ground truth is computed by mpmath. Mean absolute error:
1236 * float: 2.1686e-06
1237 * double: 1.4774e-12
1238 *
1239 * Reference:
1240 * M. Figurnov, S. Mohamed, A. Mnih "Implicit Reparameterization Gradients".
1241 * 2018
1242 */
1243};
1244
1245/*****************************************************************************
1246 * Implementation of Riemann zeta function of two arguments, based on Cephes *
1247 *****************************************************************************/
1248
1249template <typename Scalar>
1250struct zeta_retval {
1251 typedef Scalar type;
1252};
1253
1254template <typename Scalar>
1255struct zeta_impl_series {
1256 EIGEN_DEVICE_FUNC
1257 static EIGEN_STRONG_INLINE Scalar run(const Scalar) {
1258 EIGEN_STATIC_ASSERT((internal::is_same<Scalar, Scalar>::value == false),
1259 THIS_TYPE_IS_NOT_SUPPORTED);
1260 return Scalar(0);
1261 }
1262};
1263
1264template <>
1265struct zeta_impl_series<float> {
1266 EIGEN_DEVICE_FUNC
1267 static EIGEN_STRONG_INLINE bool run(float& a, float& b, float& s, const float x, const float machep) {
1268 int i = 0;
1269 while(i < 9)
1270 {
1271 i += 1;
1272 a += 1.0f;
1273 b = numext::pow( a, -x );
1274 s += b;
1275 if( numext::abs(b/s) < machep )
1276 return true;
1277 }
1278
1279 //Return whether we are done
1280 return false;
1281 }
1282};
1283
1284template <>
1285struct zeta_impl_series<double> {
1286 EIGEN_DEVICE_FUNC
1287 static EIGEN_STRONG_INLINE bool run(double& a, double& b, double& s, const double x, const double machep) {
1288 int i = 0;
1289 while( (i < 9) || (a <= 9.0) )
1290 {
1291 i += 1;
1292 a += 1.0;
1293 b = numext::pow( a, -x );
1294 s += b;
1295 if( numext::abs(b/s) < machep )
1296 return true;
1297 }
1298
1299 //Return whether we are done
1300 return false;
1301 }
1302};
1303
1304template <typename Scalar>
1305struct zeta_impl {
1306 EIGEN_DEVICE_FUNC
1307 static Scalar run(Scalar x, Scalar q) {
1308 /* zeta.c
1309 *
1310 * Riemann zeta function of two arguments
1311 *
1312 *
1313 *
1314 * SYNOPSIS:
1315 *
1316 * double x, q, y, zeta();
1317 *
1318 * y = zeta( x, q );
1319 *
1320 *
1321 *
1322 * DESCRIPTION:
1323 *
1324 *
1325 *
1326 * inf.
1327 * - -x
1328 * zeta(x,q) = > (k+q)
1329 * -
1330 * k=0
1331 *
1332 * where x > 1 and q is not a negative integer or zero.
1333 * The Euler-Maclaurin summation formula is used to obtain
1334 * the expansion
1335 *
1336 * n
1337 * - -x
1338 * zeta(x,q) = > (k+q)
1339 * -
1340 * k=1
1341 *
1342 * 1-x inf. B x(x+1)...(x+2j)
1343 * (n+q) 1 - 2j
1344 * + --------- - ------- + > --------------------
1345 * x-1 x - x+2j+1
1346 * 2(n+q) j=1 (2j)! (n+q)
1347 *
1348 * where the B2j are Bernoulli numbers. Note that (see zetac.c)
1349 * zeta(x,1) = zetac(x) + 1.
1350 *
1351 *
1352 *
1353 * ACCURACY:
1354 *
1355 * Relative error for single precision:
1356 * arithmetic domain # trials peak rms
1357 * IEEE 0,25 10000 6.9e-7 1.0e-7
1358 *
1359 * Large arguments may produce underflow in powf(), in which
1360 * case the results are inaccurate.
1361 *
1362 * REFERENCE:
1363 *
1364 * Gradshteyn, I. S., and I. M. Ryzhik, Tables of Integrals,
1365 * Series, and Products, p. 1073; Academic Press, 1980.
1366 *
1367 */
1368
1369 int i;
1370 Scalar p, r, a, b, k, s, t, w;
1371
1372 const Scalar A[] = {
1373 Scalar(12.0),
1374 Scalar(-720.0),
1375 Scalar(30240.0),
1376 Scalar(-1209600.0),
1377 Scalar(47900160.0),
1378 Scalar(-1.8924375803183791606e9), /*1.307674368e12/691*/
1379 Scalar(7.47242496e10),
1380 Scalar(-2.950130727918164224e12), /*1.067062284288e16/3617*/
1381 Scalar(1.1646782814350067249e14), /*5.109094217170944e18/43867*/
1382 Scalar(-4.5979787224074726105e15), /*8.028576626982912e20/174611*/
1383 Scalar(1.8152105401943546773e17), /*1.5511210043330985984e23/854513*/
1384 Scalar(-7.1661652561756670113e18) /*1.6938241367317436694528e27/236364091*/
1385 };
1386
1387 const Scalar maxnum = NumTraits<Scalar>::infinity();
1388 const Scalar zero = Scalar(0.0), half = Scalar(0.5), one = Scalar(1.0);
1389 const Scalar machep = cephes_helper<Scalar>::machep();
1390 const Scalar nan = NumTraits<Scalar>::quiet_NaN();
1391
1392 if( x == one )
1393 return maxnum;
1394
1395 if( x < one )
1396 {
1397 return nan;
1398 }
1399
1400 if( q <= zero )
1401 {
1402 if(q == numext::floor(q))
1403 {
1404 if (x == numext::floor(x) && long(x) % 2 == 0) {
1405 return maxnum;
1406 }
1407 else {
1408 return nan;
1409 }
1410 }
1411 p = x;
1412 r = numext::floor(p);
1413 if (p != r)
1414 return nan;
1415 }
1416
1417 /* Permit negative q but continue sum until n+q > +9 .
1418 * This case should be handled by a reflection formula.
1419 * If q<0 and x is an integer, there is a relation to
1420 * the polygamma function.
1421 */
1422 s = numext::pow( q, -x );
1423 a = q;
1424 b = zero;
1425 // Run the summation in a helper function that is specific to the floating precision
1426 if (zeta_impl_series<Scalar>::run(a, b, s, x, machep)) {
1427 return s;
1428 }
1429
1430 // If b is zero, then the tail sum will also end up being zero.
1431 // Exiting early here can prevent NaNs for some large inputs, where
1432 // the tail sum computed below has term `a` which can overflow to `inf`.
1433 if (numext::equal_strict(b, zero)) {
1434 return s;
1435 }
1436
1437 w = a;
1438 s += b*w/(x-one);
1439 s -= half * b;
1440 a = one;
1441 k = zero;
1442
1443 for( i=0; i<12; i++ )
1444 {
1445 a *= x + k;
1446 b /= w;
1447 t = a*b/A[i];
1448 s = s + t;
1449 t = numext::abs(t/s);
1450 if( t < machep ) {
1451 break;
1452 }
1453 k += one;
1454 a *= x + k;
1455 b /= w;
1456 k += one;
1457 }
1458 return s;
1459 }
1460};
1461
1462/****************************************************************************
1463 * Implementation of polygamma function, requires C++11/C99 *
1464 ****************************************************************************/
1465
1466template <typename Scalar>
1467struct polygamma_retval {
1468 typedef Scalar type;
1469};
1470
1471#if !EIGEN_HAS_C99_MATH
1472
1473template <typename Scalar>
1474struct polygamma_impl {
1475 EIGEN_DEVICE_FUNC
1476 static EIGEN_STRONG_INLINE Scalar run(Scalar n, Scalar x) {
1477 EIGEN_STATIC_ASSERT((internal::is_same<Scalar, Scalar>::value == false),
1478 THIS_TYPE_IS_NOT_SUPPORTED);
1479 return Scalar(0);
1480 }
1481};
1482
1483#else
1484
1485template <typename Scalar>
1486struct polygamma_impl {
1487 EIGEN_DEVICE_FUNC
1488 static Scalar run(Scalar n, Scalar x) {
1489 Scalar zero = 0.0, one = 1.0;
1490 Scalar nplus = n + one;
1491 const Scalar nan = NumTraits<Scalar>::quiet_NaN();
1492
1493 // Check that n is a non-negative integer
1494 if (numext::floor(n) != n || n < zero) {
1495 return nan;
1496 }
1497 // Just return the digamma function for n = 0
1498 else if (n == zero) {
1499 return digamma_impl<Scalar>::run(x);
1500 }
1501 // Use the same implementation as scipy
1502 else {
1503 Scalar factorial = numext::exp(lgamma_impl<Scalar>::run(nplus));
1504 return numext::pow(-one, nplus) * factorial * zeta_impl<Scalar>::run(nplus, x);
1505 }
1506 }
1507};
1508
1509#endif // EIGEN_HAS_C99_MATH
1510
1511/************************************************************************************************
1512 * Implementation of betainc (incomplete beta integral), based on Cephes but requires C++11/C99 *
1513 ************************************************************************************************/
1514
1515template <typename Scalar>
1516struct betainc_retval {
1517 typedef Scalar type;
1518};
1519
1520#if !EIGEN_HAS_C99_MATH
1521
1522template <typename Scalar>
1523struct betainc_impl {
1524 EIGEN_DEVICE_FUNC
1525 static EIGEN_STRONG_INLINE Scalar run(Scalar a, Scalar b, Scalar x) {
1526 EIGEN_STATIC_ASSERT((internal::is_same<Scalar, Scalar>::value == false),
1527 THIS_TYPE_IS_NOT_SUPPORTED);
1528 return Scalar(0);
1529 }
1530};
1531
1532#else
1533
1534template <typename Scalar>
1535struct betainc_impl {
1536 EIGEN_DEVICE_FUNC
1537 static EIGEN_STRONG_INLINE Scalar run(Scalar, Scalar, Scalar) {
1538 /* betaincf.c
1539 *
1540 * Incomplete beta integral
1541 *
1542 *
1543 * SYNOPSIS:
1544 *
1545 * float a, b, x, y, betaincf();
1546 *
1547 * y = betaincf( a, b, x );
1548 *
1549 *
1550 * DESCRIPTION:
1551 *
1552 * Returns incomplete beta integral of the arguments, evaluated
1553 * from zero to x. The function is defined as
1554 *
1555 * x
1556 * - -
1557 * | (a+b) | | a-1 b-1
1558 * ----------- | t (1-t) dt.
1559 * - - | |
1560 * | (a) | (b) -
1561 * 0
1562 *
1563 * The domain of definition is 0 <= x <= 1. In this
1564 * implementation a and b are restricted to positive values.
1565 * The integral from x to 1 may be obtained by the symmetry
1566 * relation
1567 *
1568 * 1 - betainc( a, b, x ) = betainc( b, a, 1-x ).
1569 *
1570 * The integral is evaluated by a continued fraction expansion.
1571 * If a < 1, the function calls itself recursively after a
1572 * transformation to increase a to a+1.
1573 *
1574 * ACCURACY (float):
1575 *
1576 * Tested at random points (a,b,x) with a and b in the indicated
1577 * interval and x between 0 and 1.
1578 *
1579 * arithmetic domain # trials peak rms
1580 * Relative error:
1581 * IEEE 0,30 10000 3.7e-5 5.1e-6
1582 * IEEE 0,100 10000 1.7e-4 2.5e-5
1583 * The useful domain for relative error is limited by underflow
1584 * of the single precision exponential function.
1585 * Absolute error:
1586 * IEEE 0,30 100000 2.2e-5 9.6e-7
1587 * IEEE 0,100 10000 6.5e-5 3.7e-6
1588 *
1589 * Larger errors may occur for extreme ratios of a and b.
1590 *
1591 * ACCURACY (double):
1592 * arithmetic domain # trials peak rms
1593 * IEEE 0,5 10000 6.9e-15 4.5e-16
1594 * IEEE 0,85 250000 2.2e-13 1.7e-14
1595 * IEEE 0,1000 30000 5.3e-12 6.3e-13
1596 * IEEE 0,10000 250000 9.3e-11 7.1e-12
1597 * IEEE 0,100000 10000 8.7e-10 4.8e-11
1598 * Outputs smaller than the IEEE gradual underflow threshold
1599 * were excluded from these statistics.
1600 *
1601 * ERROR MESSAGES:
1602 * message condition value returned
1603 * incbet domain x<0, x>1 nan
1604 * incbet underflow nan
1605 */
1606
1607 EIGEN_STATIC_ASSERT((internal::is_same<Scalar, Scalar>::value == false),
1608 THIS_TYPE_IS_NOT_SUPPORTED);
1609 return Scalar(0);
1610 }
1611};
1612
1613/* Continued fraction expansion #1 for incomplete beta integral (small_branch = True)
1614 * Continued fraction expansion #2 for incomplete beta integral (small_branch = False)
1615 */
1616template <typename Scalar>
1617struct incbeta_cfe {
1618 EIGEN_DEVICE_FUNC
1619 static EIGEN_STRONG_INLINE Scalar run(Scalar a, Scalar b, Scalar x, bool small_branch) {
1620 EIGEN_STATIC_ASSERT((internal::is_same<Scalar, float>::value ||
1621 internal::is_same<Scalar, double>::value),
1622 THIS_TYPE_IS_NOT_SUPPORTED);
1623 const Scalar big = cephes_helper<Scalar>::big();
1624 const Scalar machep = cephes_helper<Scalar>::machep();
1625 const Scalar biginv = cephes_helper<Scalar>::biginv();
1626
1627 const Scalar zero = 0;
1628 const Scalar one = 1;
1629 const Scalar two = 2;
1630
1631 Scalar xk, pk, pkm1, pkm2, qk, qkm1, qkm2;
1632 Scalar k1, k2, k3, k4, k5, k6, k7, k8, k26update;
1633 Scalar ans;
1634 int n;
1635
1636 const int num_iters = (internal::is_same<Scalar, float>::value) ? 100 : 300;
1637 const Scalar thresh =
1638 (internal::is_same<Scalar, float>::value) ? machep : Scalar(3) * machep;
1639 Scalar r = (internal::is_same<Scalar, float>::value) ? zero : one;
1640
1641 if (small_branch) {
1642 k1 = a;
1643 k2 = a + b;
1644 k3 = a;
1645 k4 = a + one;
1646 k5 = one;
1647 k6 = b - one;
1648 k7 = k4;
1649 k8 = a + two;
1650 k26update = one;
1651 } else {
1652 k1 = a;
1653 k2 = b - one;
1654 k3 = a;
1655 k4 = a + one;
1656 k5 = one;
1657 k6 = a + b;
1658 k7 = a + one;
1659 k8 = a + two;
1660 k26update = -one;
1661 x = x / (one - x);
1662 }
1663
1664 pkm2 = zero;
1665 qkm2 = one;
1666 pkm1 = one;
1667 qkm1 = one;
1668 ans = one;
1669 n = 0;
1670
1671 do {
1672 xk = -(x * k1 * k2) / (k3 * k4);
1673 pk = pkm1 + pkm2 * xk;
1674 qk = qkm1 + qkm2 * xk;
1675 pkm2 = pkm1;
1676 pkm1 = pk;
1677 qkm2 = qkm1;
1678 qkm1 = qk;
1679
1680 xk = (x * k5 * k6) / (k7 * k8);
1681 pk = pkm1 + pkm2 * xk;
1682 qk = qkm1 + qkm2 * xk;
1683 pkm2 = pkm1;
1684 pkm1 = pk;
1685 qkm2 = qkm1;
1686 qkm1 = qk;
1687
1688 if (qk != zero) {
1689 r = pk / qk;
1690 if (numext::abs(ans - r) < numext::abs(r) * thresh) {
1691 return r;
1692 }
1693 ans = r;
1694 }
1695
1696 k1 += one;
1697 k2 += k26update;
1698 k3 += two;
1699 k4 += two;
1700 k5 += one;
1701 k6 -= k26update;
1702 k7 += two;
1703 k8 += two;
1704
1705 if ((numext::abs(qk) + numext::abs(pk)) > big) {
1706 pkm2 *= biginv;
1707 pkm1 *= biginv;
1708 qkm2 *= biginv;
1709 qkm1 *= biginv;
1710 }
1711 if ((numext::abs(qk) < biginv) || (numext::abs(pk) < biginv)) {
1712 pkm2 *= big;
1713 pkm1 *= big;
1714 qkm2 *= big;
1715 qkm1 *= big;
1716 }
1717 } while (++n < num_iters);
1718
1719 return ans;
1720 }
1721};
1722
1723/* Helper functions depending on the Scalar type */
1724template <typename Scalar>
1725struct betainc_helper {};
1726
1727template <>
1728struct betainc_helper<float> {
1729 /* Core implementation, assumes a large (> 1.0) */
1730 EIGEN_DEVICE_FUNC static EIGEN_STRONG_INLINE float incbsa(float aa, float bb,
1731 float xx) {
1732 float ans, a, b, t, x, onemx;
1733 bool reversed_a_b = false;
1734
1735 onemx = 1.0f - xx;
1736
1737 /* see if x is greater than the mean */
1738 if (xx > (aa / (aa + bb))) {
1739 reversed_a_b = true;
1740 a = bb;
1741 b = aa;
1742 t = xx;
1743 x = onemx;
1744 } else {
1745 a = aa;
1746 b = bb;
1747 t = onemx;
1748 x = xx;
1749 }
1750
1751 /* Choose expansion for optimal convergence */
1752 if (b > 10.0f) {
1753 if (numext::abs(b * x / a) < 0.3f) {
1754 t = betainc_helper<float>::incbps(a, b, x);
1755 if (reversed_a_b) t = 1.0f - t;
1756 return t;
1757 }
1758 }
1759
1760 ans = x * (a + b - 2.0f) / (a - 1.0f);
1761 if (ans < 1.0f) {
1762 ans = incbeta_cfe<float>::run(a, b, x, true /* small_branch */);
1763 t = b * numext::log(t);
1764 } else {
1765 ans = incbeta_cfe<float>::run(a, b, x, false /* small_branch */);
1766 t = (b - 1.0f) * numext::log(t);
1767 }
1768
1769 t += a * numext::log(x) + lgamma_impl<float>::run(a + b) -
1770 lgamma_impl<float>::run(a) - lgamma_impl<float>::run(b);
1771 t += numext::log(ans / a);
1772 t = numext::exp(t);
1773
1774 if (reversed_a_b) t = 1.0f - t;
1775 return t;
1776 }
1777
1778 EIGEN_DEVICE_FUNC
1779 static EIGEN_STRONG_INLINE float incbps(float a, float b, float x) {
1780 float t, u, y, s;
1781 const float machep = cephes_helper<float>::machep();
1782
1783 y = a * numext::log(x) + (b - 1.0f) * numext::log1p(-x) - numext::log(a);
1784 y -= lgamma_impl<float>::run(a) + lgamma_impl<float>::run(b);
1785 y += lgamma_impl<float>::run(a + b);
1786
1787 t = x / (1.0f - x);
1788 s = 0.0f;
1789 u = 1.0f;
1790 do {
1791 b -= 1.0f;
1792 if (b == 0.0f) {
1793 break;
1794 }
1795 a += 1.0f;
1796 u *= t * b / a;
1797 s += u;
1798 } while (numext::abs(u) > machep);
1799
1800 return numext::exp(y) * (1.0f + s);
1801 }
1802};
1803
1804template <>
1805struct betainc_impl<float> {
1806 EIGEN_DEVICE_FUNC
1807 static float run(float a, float b, float x) {
1808 const float nan = NumTraits<float>::quiet_NaN();
1809 float ans, t;
1810
1811 if (a <= 0.0f) return nan;
1812 if (b <= 0.0f) return nan;
1813 if ((x <= 0.0f) || (x >= 1.0f)) {
1814 if (x == 0.0f) return 0.0f;
1815 if (x == 1.0f) return 1.0f;
1816 // mtherr("betaincf", DOMAIN);
1817 return nan;
1818 }
1819
1820 /* transformation for small aa */
1821 if (a <= 1.0f) {
1822 ans = betainc_helper<float>::incbsa(a + 1.0f, b, x);
1823 t = a * numext::log(x) + b * numext::log1p(-x) +
1824 lgamma_impl<float>::run(a + b) - lgamma_impl<float>::run(a + 1.0f) -
1825 lgamma_impl<float>::run(b);
1826 return (ans + numext::exp(t));
1827 } else {
1828 return betainc_helper<float>::incbsa(a, b, x);
1829 }
1830 }
1831};
1832
1833template <>
1834struct betainc_helper<double> {
1835 EIGEN_DEVICE_FUNC
1836 static EIGEN_STRONG_INLINE double incbps(double a, double b, double x) {
1837 const double machep = cephes_helper<double>::machep();
1838
1839 double s, t, u, v, n, t1, z, ai;
1840
1841 ai = 1.0 / a;
1842 u = (1.0 - b) * x;
1843 v = u / (a + 1.0);
1844 t1 = v;
1845 t = u;
1846 n = 2.0;
1847 s = 0.0;
1848 z = machep * ai;
1849 while (numext::abs(v) > z) {
1850 u = (n - b) * x / n;
1851 t *= u;
1852 v = t / (a + n);
1853 s += v;
1854 n += 1.0;
1855 }
1856 s += t1;
1857 s += ai;
1858
1859 u = a * numext::log(x);
1860 // TODO: gamma() is not directly implemented in Eigen.
1861 /*
1862 if ((a + b) < maxgam && numext::abs(u) < maxlog) {
1863 t = gamma(a + b) / (gamma(a) * gamma(b));
1864 s = s * t * pow(x, a);
1865 }
1866 */
1867 t = lgamma_impl<double>::run(a + b) - lgamma_impl<double>::run(a) -
1868 lgamma_impl<double>::run(b) + u + numext::log(s);
1869 return s = numext::exp(t);
1870 }
1871};
1872
1873template <>
1874struct betainc_impl<double> {
1875 EIGEN_DEVICE_FUNC
1876 static double run(double aa, double bb, double xx) {
1877 const double nan = NumTraits<double>::quiet_NaN();
1878 const double machep = cephes_helper<double>::machep();
1879 // const double maxgam = 171.624376956302725;
1880
1881 double a, b, t, x, xc, w, y;
1882 bool reversed_a_b = false;
1883
1884 if (aa <= 0.0 || bb <= 0.0) {
1885 return nan; // goto domerr;
1886 }
1887
1888 if ((xx <= 0.0) || (xx >= 1.0)) {
1889 if (xx == 0.0) return (0.0);
1890 if (xx == 1.0) return (1.0);
1891 // mtherr("incbet", DOMAIN);
1892 return nan;
1893 }
1894
1895 if ((bb * xx) <= 1.0 && xx <= 0.95) {
1896 return betainc_helper<double>::incbps(aa, bb, xx);
1897 }
1898
1899 w = 1.0 - xx;
1900
1901 /* Reverse a and b if x is greater than the mean. */
1902 if (xx > (aa / (aa + bb))) {
1903 reversed_a_b = true;
1904 a = bb;
1905 b = aa;
1906 xc = xx;
1907 x = w;
1908 } else {
1909 a = aa;
1910 b = bb;
1911 xc = w;
1912 x = xx;
1913 }
1914
1915 if (reversed_a_b && (b * x) <= 1.0 && x <= 0.95) {
1916 t = betainc_helper<double>::incbps(a, b, x);
1917 if (t <= machep) {
1918 t = 1.0 - machep;
1919 } else {
1920 t = 1.0 - t;
1921 }
1922 return t;
1923 }
1924
1925 /* Choose expansion for better convergence. */
1926 y = x * (a + b - 2.0) - (a - 1.0);
1927 if (y < 0.0) {
1928 w = incbeta_cfe<double>::run(a, b, x, true /* small_branch */);
1929 } else {
1930 w = incbeta_cfe<double>::run(a, b, x, false /* small_branch */) / xc;
1931 }
1932
1933 /* Multiply w by the factor
1934 a b _ _ _
1935 x (1-x) | (a+b) / ( a | (a) | (b) ) . */
1936
1937 y = a * numext::log(x);
1938 t = b * numext::log(xc);
1939 // TODO: gamma is not directly implemented in Eigen.
1940 /*
1941 if ((a + b) < maxgam && numext::abs(y) < maxlog && numext::abs(t) < maxlog)
1942 {
1943 t = pow(xc, b);
1944 t *= pow(x, a);
1945 t /= a;
1946 t *= w;
1947 t *= gamma(a + b) / (gamma(a) * gamma(b));
1948 } else {
1949 */
1950 /* Resort to logarithms. */
1951 y += t + lgamma_impl<double>::run(a + b) - lgamma_impl<double>::run(a) -
1952 lgamma_impl<double>::run(b);
1953 y += numext::log(w / a);
1954 t = numext::exp(y);
1955
1956 /* } */
1957 // done:
1958
1959 if (reversed_a_b) {
1960 if (t <= machep) {
1961 t = 1.0 - machep;
1962 } else {
1963 t = 1.0 - t;
1964 }
1965 }
1966 return t;
1967 }
1968};
1969
1970#endif // EIGEN_HAS_C99_MATH
1971
1972} // end namespace internal
1973
1974namespace numext {
1975
1976template <typename Scalar>
1977EIGEN_DEVICE_FUNC inline EIGEN_MATHFUNC_RETVAL(lgamma, Scalar)
1978 lgamma(const Scalar& x) {
1979 return EIGEN_MATHFUNC_IMPL(lgamma, Scalar)::run(x);
1980}
1981
1982template <typename Scalar>
1983EIGEN_DEVICE_FUNC inline EIGEN_MATHFUNC_RETVAL(digamma, Scalar)
1984 digamma(const Scalar& x) {
1985 return EIGEN_MATHFUNC_IMPL(digamma, Scalar)::run(x);
1986}
1987
1988template <typename Scalar>
1989EIGEN_DEVICE_FUNC inline EIGEN_MATHFUNC_RETVAL(zeta, Scalar)
1990zeta(const Scalar& x, const Scalar& q) {
1991 return EIGEN_MATHFUNC_IMPL(zeta, Scalar)::run(x, q);
1992}
1993
1994template <typename Scalar>
1995EIGEN_DEVICE_FUNC inline EIGEN_MATHFUNC_RETVAL(polygamma, Scalar)
1996polygamma(const Scalar& n, const Scalar& x) {
1997 return EIGEN_MATHFUNC_IMPL(polygamma, Scalar)::run(n, x);
1998}
1999
2000template <typename Scalar>
2001EIGEN_DEVICE_FUNC inline EIGEN_MATHFUNC_RETVAL(erf, Scalar)
2002 erf(const Scalar& x) {
2003 return EIGEN_MATHFUNC_IMPL(erf, Scalar)::run(x);
2004}
2005
2006template <typename Scalar>
2007EIGEN_DEVICE_FUNC inline EIGEN_MATHFUNC_RETVAL(erfc, Scalar)
2008 erfc(const Scalar& x) {
2009 return EIGEN_MATHFUNC_IMPL(erfc, Scalar)::run(x);
2010}
2011
2012template <typename Scalar>
2013EIGEN_DEVICE_FUNC inline EIGEN_MATHFUNC_RETVAL(ndtri, Scalar)
2014 ndtri(const Scalar& x) {
2015 return EIGEN_MATHFUNC_IMPL(ndtri, Scalar)::run(x);
2016}
2017
2018template <typename Scalar>
2019EIGEN_DEVICE_FUNC inline EIGEN_MATHFUNC_RETVAL(igamma, Scalar)
2020 igamma(const Scalar& a, const Scalar& x) {
2021 return EIGEN_MATHFUNC_IMPL(igamma, Scalar)::run(a, x);
2022}
2023
2024template <typename Scalar>
2025EIGEN_DEVICE_FUNC inline EIGEN_MATHFUNC_RETVAL(igamma_der_a, Scalar)
2026 igamma_der_a(const Scalar& a, const Scalar& x) {
2027 return EIGEN_MATHFUNC_IMPL(igamma_der_a, Scalar)::run(a, x);
2028}
2029
2030template <typename Scalar>
2031EIGEN_DEVICE_FUNC inline EIGEN_MATHFUNC_RETVAL(gamma_sample_der_alpha, Scalar)
2032 gamma_sample_der_alpha(const Scalar& a, const Scalar& x) {
2033 return EIGEN_MATHFUNC_IMPL(gamma_sample_der_alpha, Scalar)::run(a, x);
2034}
2035
2036template <typename Scalar>
2037EIGEN_DEVICE_FUNC inline EIGEN_MATHFUNC_RETVAL(igammac, Scalar)
2038 igammac(const Scalar& a, const Scalar& x) {
2039 return EIGEN_MATHFUNC_IMPL(igammac, Scalar)::run(a, x);
2040}
2041
2042template <typename Scalar>
2043EIGEN_DEVICE_FUNC inline EIGEN_MATHFUNC_RETVAL(betainc, Scalar)
2044 betainc(const Scalar& a, const Scalar& b, const Scalar& x) {
2045 return EIGEN_MATHFUNC_IMPL(betainc, Scalar)::run(a, b, x);
2046}
2047
2048} // end namespace numext
2049} // end namespace Eigen
2050
2051#endif // EIGEN_SPECIAL_FUNCTIONS_H
Namespace containing all symbols from the Eigen library.
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)
Definition SpecialFunctionsArrayAPI.h:90
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)
Definition SpecialFunctionsArrayAPI.h:51
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sign_op< typename Derived::Scalar >, const Derived > sign(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 TensorCwiseTernaryOp< internal::scalar_betainc_op< typename XDerived::Scalar >, const ADerived, const BDerived, const XDerived > betainc(const ADerived &a, const BDerived &b, const XDerived &x)
Definition TensorGlobalFunctions.h:24
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)
Definition SpecialFunctionsArrayAPI.h:72
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_erfc_op< typename Derived::Scalar >, const Derived > erfc(const Eigen::ArrayBase< Derived > &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_digamma_op< typename Derived::Scalar >, const Derived > digamma(const Eigen::ArrayBase< Derived > &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)
Definition SpecialFunctionsArrayAPI.h:112
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)
Definition SpecialFunctionsArrayAPI.h:28
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)
Definition SpecialFunctionsArrayAPI.h:156