1 /** @file inifcns_nstdsums.cpp
3 * Implementation of some special functions that have a representation as nested sums.
5 * classical polylogarithm Li(n,x)
6 * multiple polylogarithm Li(lst(n_1,...,n_k),lst(x_1,...,x_k)
7 * nielsen's generalized polylogarithm S(n,p,x)
8 * harmonic polylogarithm H(lst(m_1,...,m_k),x)
9 * multiple zeta value mZeta(lst(m_1,...,m_k))
12 * - All formulae used can be looked up in the following publication:
13 * Nielsen's Generalized Polylogarithms, K.S.Kolbig, SIAM J.Math.Anal. 17 (1986), pp. 1232-1258.
14 * This document will be referenced as [Kol] throughout this source code.
15 * - Classical polylogarithms (Li) and nielsen's generalized polylogarithms (S) can be numerically
16 * evaluated in the whole complex plane. And of course, there is still room for speed optimizations ;-).
17 * - The calculation of classical polylogarithms is speed up by using Euler-MacLaurin summation.
18 * - The remaining functions can only be numerically evaluated with arguments lying in the unit sphere
19 * at the moment. Sorry. The evaluation especially for mZeta is very slow ... better not use it
21 * - The functions have no series expansion. To do it, you have to convert these functions
22 * into the appropriate objects from the nestedsums library, do the expansion and convert the
24 * - Numerical testing of this implementation has been performed by doing a comparison of results
25 * between this software and the commercial M.......... 4.1.
30 * GiNaC Copyright (C) 1999-2003 Johannes Gutenberg University Mainz, Germany
32 * This program is free software; you can redistribute it and/or modify
33 * it under the terms of the GNU General Public License as published by
34 * the Free Software Foundation; either version 2 of the License, or
35 * (at your option) any later version.
37 * This program is distributed in the hope that it will be useful,
38 * but WITHOUT ANY WARRANTY; without even the implied warranty of
39 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
40 * GNU General Public License for more details.
42 * You should have received a copy of the GNU General Public License
43 * along with this program; if not, write to the Free Software
44 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
54 #include "operators.h"
55 #include "relational.h"
62 // lookup table for Euler-MacLaurin optimization
64 std::vector<std::vector<cln::cl_N> > Xn;
68 // lookup table for Euler-Zagier-Sums (used for S_n,p(x))
70 std::vector<std::vector<cln::cl_N> > Yn;
71 int ynsize = 0; // number of Yn[]
72 int ynlength = 100; // initial length of all Yn[i]
75 //////////////////////
76 // helper functions //
77 //////////////////////
80 // This function calculates the X_n. The X_n are needed for the Euler-MacLaurin summation (EMS) of
81 // classical polylogarithms.
82 // With EMS the polylogs can be calculated as follows:
83 // Li_p (x) = \sum_{n=0}^\infty X_{p-2}(n) u^{n+1}/(n+1)! with u = -log(1-x)
84 // X_0(n) = B_n (Bernoulli numbers)
85 // X_p(n) = \sum_{k=0}^n binomial(n,k) B_{n-k} / (k+1) * X_{p-1}(k)
86 // The calculation of Xn depends on X0 and X{n-1}.
87 // X_0 is special, it holds only the non-zero Bernoulli numbers with index 2 or greater.
88 // This results in a slightly more complicated algorithm for the X_n.
89 // The first index in Xn corresponds to the index of the polylog minus 2.
90 // The second index in Xn corresponds to the index from the EMS.
91 static void fill_Xn(int n)
93 // rule of thumb. needs to be improved. TODO
94 const int initsize = Digits * 3 / 2;
97 // calculate X_2 and higher (corresponding to Li_4 and higher)
98 std::vector<cln::cl_N> buf(initsize);
99 std::vector<cln::cl_N>::iterator it = buf.begin();
101 *it = -(cln::expt(cln::cl_I(2),n+1) - 1) / cln::expt(cln::cl_I(2),n+1); // i == 1
103 for (int i=2; i<=initsize; i++) {
105 result = 0; // k == 0
107 result = Xn[0][i/2-1]; // k == 0
109 for (int k=1; k<i-1; k++) {
110 if ( !(((i-k) & 1) && ((i-k) > 1)) ) {
111 result = result + cln::binomial(i,k) * Xn[0][(i-k)/2-1] * Xn[n-1][k-1] / (k+1);
114 result = result - cln::binomial(i,i-1) * Xn[n-1][i-2] / 2 / i; // k == i-1
115 result = result + Xn[n-1][i-1] / (i+1); // k == i
122 // special case to handle the X_0 correct
123 std::vector<cln::cl_N> buf(initsize);
124 std::vector<cln::cl_N>::iterator it = buf.begin();
126 *it = cln::cl_I(-3)/cln::cl_I(4); // i == 1
128 *it = cln::cl_I(17)/cln::cl_I(36); // i == 2
130 for (int i=3; i<=initsize; i++) {
132 result = -Xn[0][(i-3)/2]/2;
133 *it = (cln::binomial(i,1)/cln::cl_I(2) + cln::binomial(i,i-1)/cln::cl_I(i))*result;
136 result = Xn[0][i/2-1] + Xn[0][i/2-1]/(i+1);
137 for (int k=1; k<i/2; k++) {
138 result = result + cln::binomial(i,k*2) * Xn[0][k-1] * Xn[0][i/2-k-1] / (k*2+1);
147 std::vector<cln::cl_N> buf(initsize/2);
148 std::vector<cln::cl_N>::iterator it = buf.begin();
149 for (int i=1; i<=initsize/2; i++) {
150 *it = bernoulli(i*2).to_cl_N();
159 // This function calculates the Y_n. The Y_n are needed for the evaluation of S_{n,p}(x).
160 // The Y_n are basically Euler-Zagier sums with all m_i=1. They are subsums in the Z-sum
161 // representing S_{n,p}(x).
162 // The first index in Y_n corresponds to the parameter p minus one, i.e. the depth of the
164 // The second index in Y_n corresponds to the running index of the outermost sum in the full Z-sum
165 // representing S_{n,p}(x).
166 // The calculation of Y_n uses the values from Y_{n-1}.
167 static void fill_Yn(int n, const cln::float_format_t& prec)
169 // TODO -> get rid of the magic number
170 const int initsize = ynlength;
171 //const int initsize = initsize_Yn;
172 cln::cl_N one = cln::cl_float(1, prec);
175 std::vector<cln::cl_N> buf(initsize);
176 std::vector<cln::cl_N>::iterator it = buf.begin();
177 std::vector<cln::cl_N>::iterator itprev = Yn[n-1].begin();
178 *it = (*itprev) / cln::cl_N(n+1) * one;
181 // sums with an index smaller than the depth are zero and need not to be calculated.
182 // calculation starts with depth, which is n+2)
183 for (int i=n+2; i<=initsize+n; i++) {
184 *it = *(it-1) + (*itprev) / cln::cl_N(i) * one;
190 std::vector<cln::cl_N> buf(initsize);
191 std::vector<cln::cl_N>::iterator it = buf.begin();
194 for (int i=2; i<=initsize; i++) {
195 *it = *(it-1) + 1 / cln::cl_N(i) * one;
203 // make Yn longer ...
204 static void make_Yn_longer(int newsize, const cln::float_format_t& prec)
207 cln::cl_N one = cln::cl_float(1, prec);
209 Yn[0].resize(newsize);
210 std::vector<cln::cl_N>::iterator it = Yn[0].begin();
212 for (int i=ynlength+1; i<=newsize; i++) {
213 *it = *(it-1) + 1 / cln::cl_N(i) * one;
217 for (int n=1; n<ynsize; n++) {
218 Yn[n].resize(newsize);
219 std::vector<cln::cl_N>::iterator it = Yn[n].begin();
220 std::vector<cln::cl_N>::iterator itprev = Yn[n-1].begin();
223 for (int i=ynlength+n+1; i<=newsize+n; i++) {
224 *it = *(it-1) + (*itprev) / cln::cl_N(i) * one;
234 static cln::cl_N Li_series(int n, const cln::cl_N& x, const cln::float_format_t& prec)
236 // check if precalculated values are sufficient
238 for (int i=xnsize; i<n-1; i++) {
243 // using Euler-MacLaurin summation
245 // Li_2. X_0 is special ...
246 std::vector<cln::cl_N>::const_iterator it = Xn[0].begin();
247 cln::cl_N u = -cln::log(cln::complex(cln::cl_float(1, prec), 0)-x);
248 cln::cl_N factor = u;
249 cln::cl_N res = u - u*u/4;
251 for (int i=1; true; i++) {
253 factor = factor * u*u / (2*i * (2*i+1));
254 res = res + (*it) * factor;
255 it++; // should we check it? or rely on initsize? ...
256 if (cln::zerop(res-resbuf))
264 std::vector<cln::cl_N>::const_iterator it = Xn[n-2].begin();
265 cln::cl_N u = -cln::log(cln::complex(cln::cl_float(1, prec), 0)-x);
266 cln::cl_N factor = u;
269 for (int i=1; true; i++) {
271 factor = factor * u / (i+1);
272 res = res + (*it) * factor;
273 it++; // should we check it? or rely on initsize? ...
274 if (cln::zerop(res-resbuf))
276 // should not be needed.
277 // if (!cln::zerop(*it)) {
287 // forward declaration needed by function C below
288 static numeric S_num(int n, int p, const numeric& x);
291 // helper function for classical polylog Li
292 static cln::cl_N Li_projection(int n, const cln::cl_N& x, const cln::float_format_t& prec)
294 if (cln::realpart(x) < 0.5) {
295 return Li_series(n, x, prec);
298 return -Li_series(2, 1-x, prec) - cln::log(x) * cln::log(1-x) + cln::zeta(2);
300 cln::cl_N result = -cln::expt(cln::log(x), n-1) * cln::log(1-x) / cln::factorial(n-1);
301 for (int j=0; j<n-1; j++) {
302 result = result + (S_num(n-j-1, 1, 1).to_cl_N() - S_num(1, n-j-1, 1-x).to_cl_N())
303 * cln::expt(cln::log(x), j) / cln::factorial(j) ;
311 // helper function for classical polylog Li
312 static numeric Li_num(int n, const numeric& x)
316 return -cln::log(1-x.to_cl_N());
327 return -(1-cln::expt(cln::cl_I(2),1-n)) * cln::zeta(n);
330 // what is the desired float format?
331 // first guess: default format
332 cln::float_format_t prec = cln::default_float_format;
333 const cln::cl_N value = x.to_cl_N();
334 // second guess: the argument's format
335 if (!x.real().is_rational())
336 prec = cln::float_format(cln::the<cln::cl_F>(cln::realpart(value)));
337 else if (!x.imag().is_rational())
338 prec = cln::float_format(cln::the<cln::cl_F>(cln::imagpart(value)));
341 if (cln::abs(value) > 1) {
342 cln::cl_N result = -cln::expt(cln::log(-value),n) / cln::factorial(n);
343 // check if argument is complex. if it is real, the new polylog has to be conjugated.
344 if (cln::zerop(cln::imagpart(value))) {
346 result = result + conjugate(Li_projection(n, cln::recip(value), prec));
349 result = result - conjugate(Li_projection(n, cln::recip(value), prec));
354 result = result + Li_projection(n, cln::recip(value), prec);
357 result = result - Li_projection(n, cln::recip(value), prec);
361 for (int j=0; j<n-1; j++) {
362 add = add + (1+cln::expt(cln::cl_I(-1),n-j)) * (1-cln::expt(cln::cl_I(2),1-n+j))
363 * Li_num(n-j,1).to_cl_N() * cln::expt(cln::log(-value),j) / cln::factorial(j);
365 result = result - add;
369 return Li_projection(n, value, prec);
374 // helper function for S(n,p,x)
375 static cln::cl_N numeric_nielsen(int n, int step)
379 for (int i=1; i<n; i++) {
380 res = res + numeric_nielsen(i, step-1) / cln::cl_I(i);
390 // helper function for S(n,p,x)
392 static cln::cl_N C(int n, int p)
396 for (int k=0; k<p; k++) {
397 for (int j=0; j<=(n+k-1)/2; j++) {
401 result = result - 2 * cln::expt(cln::pi(),2*j) * S_num(n-2*j,p,1).to_cl_N() / cln::factorial(2*j);
404 result = result + 2 * cln::expt(cln::pi(),2*j) * S_num(n-2*j,p,1).to_cl_N() / cln::factorial(2*j);
411 result = result + cln::factorial(n+k-1)
412 * cln::expt(cln::pi(),2*j) * S_num(n+k-2*j,p-k,1).to_cl_N()
413 / (cln::factorial(k) * cln::factorial(n-1) * cln::factorial(2*j));
416 result = result - cln::factorial(n+k-1)
417 * cln::expt(cln::pi(),2*j) * S_num(n+k-2*j,p-k,1).to_cl_N()
418 / (cln::factorial(k) * cln::factorial(n-1) * cln::factorial(2*j));
423 result = result - cln::factorial(n+k-1) * cln::expt(cln::pi(),2*j) * S_num(n+k-2*j,p-k,1).to_cl_N()
424 / (cln::factorial(k) * cln::factorial(n-1) * cln::factorial(2*j));
427 result = result + cln::factorial(n+k-1)
428 * cln::expt(cln::pi(),2*j) * S_num(n+k-2*j,p-k,1).to_cl_N()
429 / (cln::factorial(k) * cln::factorial(n-1) * cln::factorial(2*j));
437 if (((np)/2+n) & 1) {
438 result = -result - cln::expt(cln::pi(),np) / (np * cln::factorial(n-1) * cln::factorial(p));
441 result = -result + cln::expt(cln::pi(),np) / (np * cln::factorial(n-1) * cln::factorial(p));
449 // helper function for S(n,p,x)
450 // [Kol] remark to (9.1)
451 static cln::cl_N a_k(int k)
460 for (int m=2; m<=k; m++) {
461 result = result + cln::expt(cln::cl_N(-1),m) * cln::zeta(m) * a_k(k-m);
468 // helper function for S(n,p,x)
469 // [Kol] remark to (9.1)
470 static cln::cl_N b_k(int k)
479 for (int m=2; m<=k; m++) {
480 result = result + cln::expt(cln::cl_N(-1),m) * cln::zeta(m) * b_k(k-m);
487 // helper function for S(n,p,x)
488 static cln::cl_N S_series(int n, int p, const cln::cl_N& x, const cln::float_format_t& prec)
491 return Li_series(n+1, x, prec);
494 // TODO -> check for vector boundaries and do missing calculations
496 // check if precalculated values are sufficient
498 for (int i=ynsize; i<p-1; i++) {
503 // should be done otherwise
504 cln::cl_N xf = x * cln::cl_float(1, prec);
507 cln::cl_N resultbuffer;
509 for (i=p; true; i++) {
510 resultbuffer = result;
511 if (i-p >= ynlength) {
513 make_Yn_longer(ynlength*2, prec);
515 result = result + cln::expt(xf,i) / cln::expt(cln::cl_I(i),n+1) * Yn[p-2][i-p]; // should we check it? or rely on magic number? ...
516 if (cln::zerop(result-resultbuffer)) {
525 // helper function for S(n,p,x)
526 static cln::cl_N S_projection(int n, int p, const cln::cl_N& x, const cln::float_format_t& prec)
529 if (cln::abs(cln::realpart(x)) > cln::cl_F("0.5")) {
531 cln::cl_N result = cln::expt(cln::cl_I(-1),p) * cln::expt(cln::log(x),n)
532 * cln::expt(cln::log(1-x),p) / cln::factorial(n) / cln::factorial(p);
534 for (int s=0; s<n; s++) {
536 for (int r=0; r<p; r++) {
537 res2 = res2 + cln::expt(cln::cl_I(-1),r) * cln::expt(cln::log(1-x),r)
538 * S_series(p-r,n-s,1-x,prec) / cln::factorial(r);
540 result = result + cln::expt(cln::log(x),s) * (S_num(n-s,p,1).to_cl_N() - res2) / cln::factorial(s);
546 return S_series(n, p, x, prec);
550 // helper function for S(n,p,x)
551 static numeric S_num(int n, int p, const numeric& x)
555 // [Kol] (2.22) with (2.21)
556 return cln::zeta(p+1);
561 return cln::zeta(n+1);
566 for (int nu=0; nu<n; nu++) {
567 for (int rho=0; rho<=p; rho++) {
568 result = result + b_k(n-nu-1) * b_k(p-rho) * a_k(nu+rho+1)
569 * cln::factorial(nu+rho+1) / cln::factorial(rho) / cln::factorial(nu+1);
572 result = result * cln::expt(cln::cl_I(-1),n+p-1);
579 return -(1-cln::expt(cln::cl_I(2),-n)) * cln::zeta(n+1);
581 // throw std::runtime_error("don't know how to evaluate this function!");
584 // what is the desired float format?
585 // first guess: default format
586 cln::float_format_t prec = cln::default_float_format;
587 const cln::cl_N value = x.to_cl_N();
588 // second guess: the argument's format
589 if (!x.real().is_rational())
590 prec = cln::float_format(cln::the<cln::cl_F>(cln::realpart(value)));
591 else if (!x.imag().is_rational())
592 prec = cln::float_format(cln::the<cln::cl_F>(cln::imagpart(value)));
596 if (cln::realpart(value) < -0.5) {
598 cln::cl_N result = cln::expt(cln::cl_I(-1),p) * cln::expt(cln::log(value),n)
599 * cln::expt(cln::log(1-value),p) / cln::factorial(n) / cln::factorial(p);
601 for (int s=0; s<n; s++) {
603 for (int r=0; r<p; r++) {
604 res2 = res2 + cln::expt(cln::cl_I(-1),r) * cln::expt(cln::log(1-value),r)
605 * S_num(p-r,n-s,1-value).to_cl_N() / cln::factorial(r);
607 result = result + cln::expt(cln::log(value),s) * (S_num(n-s,p,1).to_cl_N() - res2) / cln::factorial(s);
614 if (cln::abs(value) > 1) {
618 for (int s=0; s<p; s++) {
619 for (int r=0; r<=s; r++) {
620 result = result + cln::expt(cln::cl_I(-1),s) * cln::expt(cln::log(-value),r) * cln::factorial(n+s-r-1)
621 / cln::factorial(r) / cln::factorial(s-r) / cln::factorial(n-1)
622 * S_num(n+s-r,p-s,cln::recip(value)).to_cl_N();
625 result = result * cln::expt(cln::cl_I(-1),n);
628 for (int r=0; r<n; r++) {
629 res2 = res2 + cln::expt(cln::log(-value),r) * C(n-r,p) / cln::factorial(r);
631 res2 = res2 + cln::expt(cln::log(-value),n+p) / cln::factorial(n+p);
633 result = result + cln::expt(cln::cl_I(-1),p) * res2;
638 return S_projection(n, p, value, prec);
643 // helper function for multiple polylogarithm
644 static cln::cl_N numeric_zsum(int n, std::vector<cln::cl_N>& x, std::vector<cln::cl_N>& m)
650 for (int i=1; i<n; i++) {
651 std::vector<cln::cl_N>::iterator be;
652 std::vector<cln::cl_N>::iterator en;
656 std::vector<cln::cl_N> xbuf(be, en);
660 std::vector<cln::cl_N> mbuf(be, en);
661 res = res + cln::expt(x[0],i) / cln::expt(i,m[0]) * numeric_zsum(i, xbuf, mbuf);
667 // helper function for harmonic polylogarithm
668 static cln::cl_N numeric_harmonic(int n, std::vector<cln::cl_N>& m)
674 for (int i=1; i<n; i++) {
675 std::vector<cln::cl_N>::iterator be;
676 std::vector<cln::cl_N>::iterator en;
680 std::vector<cln::cl_N> mbuf(be, en);
681 res = res + cln::recip(cln::expt(i,m[0])) * numeric_harmonic(i, mbuf);
687 /////////////////////////////
688 // end of helper functions //
689 /////////////////////////////
692 // Polylogarithm and multiple polylogarithm
694 static ex Li_eval(const ex& x1, const ex& x2)
700 return Li(x1,x2).hold();
704 static ex Li_evalf(const ex& x1, const ex& x2)
706 // classical polylogs
707 if (is_a<numeric>(x1) && is_a<numeric>(x2)) {
708 return Li_num(ex_to<numeric>(x1).to_int(), ex_to<numeric>(x2));
711 else if (is_a<lst>(x1) && is_a<lst>(x2)) {
712 for (int i=0; i<x1.nops(); i++) {
713 if (!is_a<numeric>(x1.op(i)))
714 return Li(x1,x2).hold();
715 if (!is_a<numeric>(x2.op(i)))
716 return Li(x1,x2).hold();
718 return Li(x1,x2).hold();
721 cln::cl_N m_1 = ex_to<numeric>(x1.op(x1.nops()-1)).to_cl_N();
722 cln::cl_N x_1 = ex_to<numeric>(x2.op(x2.nops()-1)).to_cl_N();
723 std::vector<cln::cl_N> x;
724 std::vector<cln::cl_N> m;
725 const int nops = ex_to<numeric>(x1.nops()).to_int();
726 for (int i=nops-2; i>=0; i--) {
727 m.push_back(ex_to<numeric>(x1.op(i)).to_cl_N());
728 x.push_back(ex_to<numeric>(x2.op(i)).to_cl_N());
733 for (int i=nops; true; i++) {
735 res = res + cln::expt(x_1,i) / cln::expt(i,m_1) * numeric_zsum(i, x, m);
736 if (cln::zerop(res-resbuf))
744 return Li(x1,x2).hold();
747 static ex Li_series(const ex& x1, const ex& x2, const relational& rel, int order, unsigned options)
750 seq.push_back(expair(Li(x1,x2), 0));
751 return pseries(rel,seq);
754 REGISTER_FUNCTION(Li, eval_func(Li_eval).evalf_func(Li_evalf).do_not_evalf_params().series_func(Li_series));
757 // Nielsen's generalized polylogarithm
759 static ex S_eval(const ex& x1, const ex& x2, const ex& x3)
764 return S(x1,x2,x3).hold();
767 static ex S_evalf(const ex& x1, const ex& x2, const ex& x3)
769 if (is_a<numeric>(x1) && is_a<numeric>(x2) && is_a<numeric>(x3)) {
770 if ((x3 == -1) && (x2 != 1)) {
771 // no formula to evaluate this ... sorry
772 // return S(x1,x2,x3).hold();
774 return S_num(ex_to<numeric>(x1).to_int(), ex_to<numeric>(x2).to_int(), ex_to<numeric>(x3));
776 return S(x1,x2,x3).hold();
779 static ex S_series(const ex& x1, const ex& x2, const ex& x3, const relational& rel, int order, unsigned options)
782 seq.push_back(expair(S(x1,x2,x3), 0));
783 return pseries(rel,seq);
786 REGISTER_FUNCTION(S, eval_func(S_eval).evalf_func(S_evalf).do_not_evalf_params().series_func(S_series));
789 // Harmonic polylogarithm
791 static ex H_eval(const ex& x1, const ex& x2)
793 return H(x1,x2).hold();
796 static ex H_evalf(const ex& x1, const ex& x2)
798 if (is_a<lst>(x1) && is_a<numeric>(x2)) {
799 for (int i=0; i<x1.nops(); i++) {
800 if (!is_a<numeric>(x1.op(i)))
801 return H(x1,x2).hold();
804 return H(x1,x2).hold();
807 cln::cl_N m_1 = ex_to<numeric>(x1.op(x1.nops()-1)).to_cl_N();
808 cln::cl_N x_1 = ex_to<numeric>(x2).to_cl_N();
809 std::vector<cln::cl_N> m;
810 const int nops = ex_to<numeric>(x1.nops()).to_int();
811 for (int i=nops-2; i>=0; i--) {
812 m.push_back(ex_to<numeric>(x1.op(i)).to_cl_N());
817 for (int i=nops; true; i++) {
819 res = res + cln::expt(x_1,i) / cln::expt(i,m_1) * numeric_harmonic(i, m);
820 if (cln::zerop(res-resbuf))
828 return H(x1,x2).hold();
831 static ex H_series(const ex& x1, const ex& x2, const relational& rel, int order, unsigned options)
834 seq.push_back(expair(H(x1,x2), 0));
835 return pseries(rel,seq);
838 REGISTER_FUNCTION(H, eval_func(H_eval).evalf_func(H_evalf).do_not_evalf_params().series_func(H_series));
841 // Multiple zeta value
843 static ex mZeta_eval(const ex& x1)
845 return mZeta(x1).hold();
848 static ex mZeta_evalf(const ex& x1)
851 for (int i=0; i<x1.nops(); i++) {
852 if (!is_a<numeric>(x1.op(i)))
853 return mZeta(x1).hold();
856 cln::cl_N m_1 = ex_to<numeric>(x1.op(x1.nops()-1)).to_cl_N();
857 std::vector<cln::cl_N> m;
858 const int nops = ex_to<numeric>(x1.nops()).to_int();
859 for (int i=nops-2; i>=0; i--) {
860 m.push_back(ex_to<numeric>(x1.op(i)).to_cl_N());
863 cln::float_format_t prec = cln::default_float_format;
864 cln::cl_N res = cln::complex(cln::cl_float(0, prec), 0);
866 for (int i=nops; true; i++) {
867 // to infinity and beyond ... timewise
869 res = res + cln::recip(cln::expt(i,m_1)) * numeric_harmonic(i, m);
870 if (cln::zerop(res-resbuf))
878 return mZeta(x1).hold();
881 static ex mZeta_series(const ex& x1, const relational& rel, int order, unsigned options)
884 seq.push_back(expair(mZeta(x1), 0));
885 return pseries(rel,seq);
888 REGISTER_FUNCTION(mZeta, eval_func(mZeta_eval).evalf_func(mZeta_evalf).do_not_evalf_params().series_func(mZeta_series));