flexiblesusy is hosted by Hepforge, IPPP Durham
FlexibleSUSY
Li.cpp
Go to the documentation of this file.
1// ====================================================================
2// This file is part of Polylogarithm.
3//
4// Polylogarithm is licenced under the MIT License.
5// ====================================================================
6
7#include "Li.hpp"
8#include "Li2.hpp"
9#include "Li3.hpp"
10#include "Li4.hpp"
11#include "Li5.hpp"
12#include "Li6.hpp"
13#include "eta.hpp"
14#include "factorial.hpp"
15#include "harmonic.hpp"
16#include "zeta.hpp"
17#include <cmath>
18#include <complex>
19#include <cstdint>
20#include <limits>
21
22namespace flexiblesusy {
23
24namespace {
25 constexpr double inf = std::numeric_limits<double>::infinity();
26 constexpr double nan = std::numeric_limits<double>::quiet_NaN();
27 constexpr double PI = 3.1415926535897932;
28
29 constexpr bool is_even(int64_t n) noexcept { return n % 2 == 0; }
30
31 constexpr bool is_finite(const std::complex<double>& z) noexcept
32 {
33 return std::isfinite(std::real(z)) && std::isfinite(std::imag(z));
34 }
35
37 std::complex<double> clog(const std::complex<double>& z) noexcept
38 {
39 const double n = std::hypot(std::real(z), std::imag(z));
40 double a = std::arg(z);
41
42 if (std::imag(z) == 0.0 && a < 0.0) {
43 a = -a;
44 }
45
46 return { std::log(n), a };
47 }
48
51 std::complex<double> Li_series(int64_t n, const std::complex<double>& z) noexcept
52 {
53 std::complex<double> sum(0.0, 0.0), sum_old(0.0, 0.0), p = z;
54 int64_t k = 0;
55
56 do {
57 k++;
58 sum_old = sum;
59 sum += p;
60 p *= z*std::pow(k/(1.0 + k), n);
61 if (!is_finite(p)) { break; }
62 } while (sum != sum_old &&
63 k < std::numeric_limits<int64_t>::max() - 2);
64
65 return sum;
66 }
67
69 std::complex<double> Li_unity_pos(int64_t n, const std::complex<double>& z) noexcept
70 {
71 const std::complex<double> lnz = clog(z);
72 const std::complex<double> lnz2 = lnz*lnz;
73 std::complex<double> sum(zeta(n), 0.0), p(1.0, 0.0);
74
75 for (int64_t j = 1; j < n - 1; ++j) {
76 p *= lnz/static_cast<double>(j);
77 sum += zeta(n - j)*p;
78 }
79
80 p *= lnz/static_cast<double>(n - 1);
81 sum += (harmonic(n - 1) - clog(-lnz))*p;
82
83 p *= lnz/static_cast<double>(n);
84 sum += zeta(0)*p;
85
86 p *= lnz/static_cast<double>(n + 1);
87 sum += zeta(-1)*p;
88
89 for (int64_t j = (n + 3); j < std::numeric_limits<int64_t>::max() - 2; j += 2) {
90 p *= lnz2/static_cast<double>((j - 1)*j);
91 const auto old_sum = sum;
92 sum += zeta(n - j)*p;
93 if (sum == old_sum) { break; }
94 }
95
96 return sum;
97 }
98
100 std::complex<double> stable_pow(const std::complex<double>& z, int64_t n) noexcept
101 {
102 if (std::imag(z) == 0) {
103 return { std::pow(std::real(z), n), 0.0 };
104 } else if (std::real(z) == 0) {
105 const double p = std::pow(std::imag(z), n);
106 if (n % 4 == 0) {
107 return { p, 0.0 };
108 } else if (n % 2 == 0) {
109 return { -p, 0.0 };
110 } else if ((n - 1) % 4 == 0) {
111 return { 0.0, p };
112 }
113 return { 0.0, -p };
114 }
115 return std::pow(z, n);
116 }
117
119 std::complex<double> Li_unity_neg(int64_t n, const std::complex<double>& z) noexcept
120 {
121 const std::complex<double> lnz = clog(z);
122 const std::complex<double> lnz2 = lnz*lnz;
123 std::complex<double> sum = std::tgamma(1 - n)*stable_pow(-lnz, n - 1);
124 std::complex<double> lnzk, sum_old, term;
125 int64_t k;
126
127 if (is_even(n)) {
128 lnzk = lnz;
129 k = 1;
130 } else {
131 lnzk = lnz2;
132 sum += zeta(n);
133 k = 2;
134 }
135
136 do {
137 term = zeta(n - k)*inv_fac(k)*lnzk;
138 if (!is_finite(term)) { break; }
139 sum_old = sum;
140 sum += term;
141 lnzk *= lnz2;
142 k += 2;
143 } while (sum != sum_old);
144
145 return sum;
146 }
147
149 std::complex<double> Li_rest(int64_t n, const std::complex<double>& z) noexcept
150 {
151 const std::complex<double> lnz = clog(-z);
152 const std::complex<double> lnz2 = lnz*lnz;
153 const int64_t kmax = is_even(n) ? n/2 : (n - 1)/2;
154 std::complex<double> p = is_even(n) ? 1.0 : lnz;
155 std::complex<double> sum(0.0, 0.0), old_sum;
156
157 for (int64_t k = kmax; k != 0; --k) {
158 const double ifac = inv_fac(n - 2*k);
159 if (ifac == 0) { return 2.0*sum; }
160 sum += neg_eta(2*k)*ifac*p;
161 p *= lnz2;
162 if (sum == old_sum) { break; }
163 }
164
165 return 2.0*sum - p*inv_fac(n);
166 }
167
168} // anonymous namespace
169
181std::complex<double> Li(int64_t n, const std::complex<double>& z) noexcept
182{
183 if (std::isnan(std::real(z)) || std::isnan(std::imag(z))) {
184 return {nan, nan};
185 } else if (std::isinf(std::real(z)) || std::isinf(std::imag(z))) {
186 return {-inf, 0.0};
187 } else if (z == 0.0) {
188 return {0.0, 0.0};
189 } else if (z == 1.0) {
190 if (n <= 0) {
191 return {inf, inf};
192 }
193 return {zeta(n), 0.0};
194 } else if (z == -1.0) {
195 return {neg_eta(n), 0.0};
196 } else if (n < -1) {
197 // arXiv:2010.09860
198 const double nz = std::norm(z);
199 const double nl = std::norm(clog(z));
200 if (4*PI*PI*nz < nl) {
201 return Li_series(n, z);
202 } else if (nl < 0.512*0.512*4*PI*PI) {
203 return Li_unity_neg(n, z);
204 }
205 const auto sqrtz = std::sqrt(z);
206 return std::exp2(n - 1)*(Li(n, sqrtz) + Li(n, -sqrtz));
207 } else if (n == -1) {
208 return z/((1.0 - z)*(1.0 - z));
209 } else if (n == 0) {
210 return z/(1.0 - z);
211 } else if (n == 1) {
212 return -clog(1.0 - z);
213 } else if (n == 2) {
214 return Li2(z);
215 } else if (n == 3) {
216 return Li3(z);
217 } else if (n == 4) {
218 return Li4(z);
219 } else if (n == 5) {
220 return Li5(z);
221 } else if (n == 6) {
222 return Li6(z);
223 } else if (std::norm(z) <= 0.75*0.75) {
224 return Li_series(n, z);
225 } else if (std::norm(z) >= 1.4*1.4) {
226 const double sgn = is_even(n) ? -1.0 : 1.0;
227 return sgn*Li_series(n, 1.0/z) + Li_rest(n, z);
228 }
229 return Li_unity_pos(n, z);
230}
231
232} // namespace flexiblesusy
std::complex< double > Li_series(int64_t n, const std::complex< double > &z) noexcept
Definition: Li.cpp:51
std::complex< double > clog(const std::complex< double > &z) noexcept
complex logarithm, converts -0.0 to 0.0
Definition: Li.cpp:37
std::complex< double > stable_pow(const std::complex< double > &z, int64_t n) noexcept
returns z^n, treating Re(z) == 0 and Im(z) == 0 in a stable way
Definition: Li.cpp:100
constexpr bool is_even(int64_t n) noexcept
Definition: Li.cpp:29
std::complex< double > Li_unity_neg(int64_t n, const std::complex< double > &z) noexcept
Series expansion of Li_n(z) around z ~ 1, n < 0.
Definition: Li.cpp:119
std::complex< double > Li_rest(int64_t n, const std::complex< double > &z) noexcept
returns remainder from inversion formula
Definition: Li.cpp:149
std::complex< double > Li_unity_pos(int64_t n, const std::complex< double > &z) noexcept
Series expansion of Li_n(z) around z ~ 1, n > 0.
Definition: Li.cpp:69
constexpr T norm(const Complex< T > &z) noexcept
Definition: complex.hpp:66
double neg_eta(int64_t n) noexcept
negative Dirichlet eta function
Definition: eta.cpp:81
std::complex< double > Li(int64_t n, const std::complex< double > &z) noexcept
Complex polylogarithm .
Definition: Li.cpp:181
constexpr T arg(const Complex< T > &z) noexcept
Definition: complex.hpp:42
double zeta(int64_t n) noexcept
Riemann zeta function for integer arguments.
Definition: zeta.cpp:83
std::complex< double > Li5(const std::complex< double > &z_) noexcept
Complex polylogarithm .
Definition: Li5.cpp:42
double Li4(double x) noexcept
Real 4-th order polylogarithm .
Definition: Li4.cpp:136
double Li2(double x) noexcept
Real dilogarithm .
Definition: Li2.cpp:68
bool is_finite(const gsl_vector *x)
Returns true if GSL vector contains only finite elements, false otherwise.
Definition: gsl_utils.cpp:32
double Li3(double x) noexcept
Real trilogarithm .
Definition: Li3.cpp:104
Complex< T > log(const Complex< T > &z) noexcept
Definition: complex.hpp:54
double harmonic(int64_t n) noexcept
harmonic number n
Definition: harmonic.cpp:54
double inv_fac(int64_t n) noexcept
negative Dirichlet eta function
Definition: factorial.cpp:80
std::complex< double > Li6(const std::complex< double > &z_) noexcept
Complex polylogarithm .
Definition: Li6.cpp:42
auto sum(Idx ini, Idx fin, Function f) -> decltype(EvalEigenXpr< Idx >(ini, f))
Definition: sum.hpp:91