|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
#pragma once |
|
|
|
#include "cephes/psi.h" |
|
#include "cephes/zeta.h" |
|
#include "config.h" |
|
#include "error.h" |
|
#include "trig.h" |
|
|
|
namespace xsf { |
|
namespace detail { |
|
|
|
|
|
constexpr double digamma_posroot = 1.4616321449683623; |
|
|
|
constexpr double digamma_posrootval = -9.2412655217294275e-17; |
|
|
|
constexpr double digamma_negroot = -0.504083008264455409; |
|
|
|
constexpr double digamma_negrootval = 7.2897639029768949e-17; |
|
|
|
template <typename T> |
|
XSF_HOST_DEVICE T digamma_zeta_series(T z, double root, double rootval) { |
|
T res = rootval; |
|
T coeff = -1.0; |
|
|
|
z = z - root; |
|
T term; |
|
for (int n = 1; n < 100; n++) { |
|
coeff *= -z; |
|
term = coeff * cephes::zeta(n + 1, root); |
|
res += term; |
|
if (std::abs(term) < std::numeric_limits<double>::epsilon() * std::abs(res)) { |
|
break; |
|
} |
|
} |
|
return res; |
|
} |
|
|
|
XSF_HOST_DEVICE inline std::complex<double> |
|
digamma_forward_recurrence(std::complex<double> z, std::complex<double> psiz, int n) { |
|
|
|
|
|
|
|
|
|
|
|
|
|
std::complex<double> res = psiz; |
|
|
|
for (int k = 0; k < n; k++) { |
|
res += 1.0 / (z + static_cast<double>(k)); |
|
} |
|
return res; |
|
} |
|
|
|
XSF_HOST_DEVICE inline std::complex<double> |
|
digamma_backward_recurrence(std::complex<double> z, std::complex<double> psiz, int n) { |
|
|
|
std::complex<double> res = psiz; |
|
|
|
for (int k = 1; k < n + 1; k++) { |
|
res -= 1.0 / (z - static_cast<double>(k)); |
|
} |
|
return res; |
|
} |
|
|
|
XSF_HOST_DEVICE inline std::complex<double> digamma_asymptotic_series(std::complex<double> z) { |
|
|
|
|
|
|
|
double bernoulli2k[] = {0.166666666666666667, -0.0333333333333333333, 0.0238095238095238095, |
|
-0.0333333333333333333, 0.0757575757575757576, -0.253113553113553114, |
|
1.16666666666666667, -7.09215686274509804, 54.9711779448621554, |
|
-529.124242424242424, 6192.12318840579710, -86580.2531135531136, |
|
1425517.16666666667, -27298231.0678160920, 601580873.900642368, |
|
-15116315767.0921569}; |
|
std::complex<double> rzz = 1.0 / z / z; |
|
std::complex<double> zfac = 1.0; |
|
std::complex<double> term; |
|
std::complex<double> res; |
|
|
|
if (!(std::isfinite(z.real()) && std::isfinite(z.imag()))) { |
|
|
|
|
|
|
|
|
|
return std::log(z); |
|
} |
|
|
|
res = std::log(z) - 0.5 / z; |
|
|
|
for (int k = 1; k < 17; k++) { |
|
zfac *= rzz; |
|
term = -bernoulli2k[k - 1] * zfac / (2 * static_cast<double>(k)); |
|
res += term; |
|
if (std::abs(term) < std::numeric_limits<double>::epsilon() * std::abs(res)) { |
|
break; |
|
} |
|
} |
|
return res; |
|
} |
|
|
|
} |
|
|
|
XSF_HOST_DEVICE inline double digamma(double z) { |
|
|
|
|
|
|
|
if (std::abs(z - detail::digamma_negroot) < 0.3) { |
|
return detail::digamma_zeta_series(z, detail::digamma_negroot, detail::digamma_negrootval); |
|
} |
|
return cephes::psi(z); |
|
} |
|
|
|
XSF_HOST_DEVICE inline float digamma(float z) { return static_cast<float>(digamma(static_cast<double>(z))); } |
|
|
|
XSF_HOST_DEVICE inline std::complex<double> digamma(std::complex<double> z) { |
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
double absz = std::abs(z); |
|
std::complex<double> res = 0; |
|
|
|
|
|
int smallabsz = 16; |
|
|
|
|
|
|
|
|
|
if (z.real() <= 0.0 && std::ceil(z.real()) == z) { |
|
|
|
set_error("digamma", SF_ERROR_SINGULAR, NULL); |
|
return {std::numeric_limits<double>::quiet_NaN(), std::numeric_limits<double>::quiet_NaN()}; |
|
} |
|
if (std::abs(z - detail::digamma_negroot) < 0.3) { |
|
|
|
return detail::digamma_zeta_series(z, detail::digamma_negroot, detail::digamma_negrootval); |
|
} |
|
|
|
if (z.real() < 0 and std::abs(z.imag()) < smallabsz) { |
|
|
|
|
|
|
|
|
|
res = -M_PI * cospi(z) / sinpi(z); |
|
z = 1.0 - z; |
|
absz = std::abs(z); |
|
} |
|
|
|
if (absz < 0.5) { |
|
|
|
|
|
res = -1.0 / z; |
|
z += 1.0; |
|
absz = std::abs(z); |
|
} |
|
|
|
if (std::abs(z - detail::digamma_posroot) < 0.5) { |
|
res += detail::digamma_zeta_series(z, detail::digamma_posroot, detail::digamma_posrootval); |
|
} else if (absz > smallabsz) { |
|
res += detail::digamma_asymptotic_series(z); |
|
} else if (z.real() >= 0.0) { |
|
double n = std::trunc(smallabsz - absz) + 1; |
|
std::complex<double> init = detail::digamma_asymptotic_series(z + n); |
|
res += detail::digamma_backward_recurrence(z + n, init, n); |
|
} else { |
|
|
|
double n = std::trunc(smallabsz - absz) - 1; |
|
std::complex<double> init = detail::digamma_asymptotic_series(z - n); |
|
res += detail::digamma_forward_recurrence(z - n, init, n); |
|
} |
|
return res; |
|
} |
|
|
|
XSF_HOST_DEVICE inline std::complex<float> digamma(std::complex<float> z) { |
|
return static_cast<std::complex<float>>(digamma(static_cast<std::complex<double>>(z))); |
|
} |
|
|
|
} |
|
|