Code coverage tests

This page documents the degree to which the PARI/GP source code is tested by our public test suite, distributed with the source distribution in directory src/test/. This is measured by the gcov utility; we then process gcov output using the lcov frond-end.

We test a few variants depending on Configure flags on the pari.math.u-bordeaux.fr machine (x86_64 architecture), and agregate them in the final report:

The target is to exceed 90% coverage for all mathematical modules (given that branches depending on DEBUGLEVEL or DEBUGMEM are not covered). This script is run to produce the results below.

LCOV - code coverage report
Current view: top level - basemath - dirichlet.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.18.0 lcov report (development 29806-4d001396c7) Lines: 481 504 95.4 %
Date: 2024-12-21 09:08:57 Functions: 42 45 93.3 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2015  The PARI group.
       2             : 
       3             : This file is part of the PARI/GP package.
       4             : 
       5             : PARI/GP is free software; you can redistribute it and/or modify it under the
       6             : terms of the GNU General Public License as published by the Free Software
       7             : Foundation; either version 2 of the License, or (at your option) any later
       8             : version. It is distributed in the hope that it will be useful, but WITHOUT
       9             : ANY WARRANTY WHATSOEVER.
      10             : 
      11             : Check the License for details. You should have received a copy of it, along
      12             : with the package; see the file 'COPYING'. If not, write to the Free Software
      13             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      14             : 
      15             : /********************************************************************/
      16             : /**                                                                **/
      17             : /**           Dirichlet series through Euler product               **/
      18             : /**                                                                **/
      19             : /********************************************************************/
      20             : #include "pari.h"
      21             : #include "paripriv.h"
      22             : 
      23             : static void
      24          28 : err_direuler(GEN x)
      25          28 : { pari_err_DOMAIN("direuler","constant term","!=", gen_1,x); }
      26             : 
      27             : /* s = t_POL (tolerate t_SER of valuation 0) of constant term = 1
      28             :  * d = minimal such that p^d > X
      29             :  * V indexed by 1..X will contain the a_n
      30             :  * v[1..n] contains the indices nj such that V[nj] != 0 */
      31             : static long
      32       28700 : dirmuleuler_small(GEN V, GEN v, long n, ulong p, GEN s, long d)
      33             : {
      34       28700 :   long i, j, m = n, D = minss(d+2, lg(s));
      35       28700 :   ulong q = 1, X = lg(V)-1;
      36             : 
      37       94724 :   for (i = 3, q = p; i < D; i++, q *= p) /* q*p does not overflow */
      38             :   {
      39       66024 :     GEN aq = gel(s,i);
      40       66024 :     if (gequal0(aq)) continue;
      41             :     /* j = 1 */
      42       53753 :     gel(V,q) = aq;
      43       53753 :     v[++n] = q;
      44     3268013 :     for (j = 2; j <= m; j++)
      45             :     {
      46     3214260 :       ulong nj = umuluu_le(uel(v,j), q, X);
      47     3214260 :       if (!nj) continue;
      48      192017 :       gel(V,nj) = gmul(aq, gel(V,v[j]));
      49      192017 :       v[++n] = nj;
      50             :     }
      51             :   }
      52       28700 :   return n;
      53             : }
      54             : 
      55             : /* ap != 0 for efficiency, p > sqrt(X) */
      56             : static void
      57      308798 : dirmuleuler_large(GEN V, ulong p, GEN ap)
      58             : {
      59      308798 :   long j, jp, X = lg(V)-1;
      60      308798 :   gel(V,p) = ap;
      61     1506547 :   for (j = 2, jp = 2*p; jp <= X; j++, jp += p) gel(V,jp) = gmul(ap, gel(V,j));
      62      308798 : }
      63             : 
      64             : static ulong
      65       10269 : direulertou(GEN a, GEN fl(GEN))
      66             : {
      67       10269 :   if (typ(a) != t_INT)
      68             :   {
      69          49 :     a = fl(a);
      70          28 :     if (typ(a) != t_INT) pari_err_TYPE("direuler", a);
      71             :   }
      72       10248 :   return signe(a)<=0 ? 0: itou(a);
      73             : }
      74             : 
      75             : static GEN
      76        3724 : direuler_Sbad(GEN V, GEN v, GEN Sbad, ulong *n)
      77             : {
      78        3724 :   long i, l = lg(Sbad);
      79        3724 :   ulong X = lg(V)-1;
      80        3724 :   GEN pbad = gen_1;
      81        9646 :   for (i = 1; i < l; i++)
      82             :   {
      83        5957 :     GEN ai = gel(Sbad,i);
      84             :     ulong q;
      85        5957 :     if (typ(ai) != t_VEC || lg(ai) != 3)
      86          14 :       pari_err_TYPE("direuler [bad primes]",ai);
      87        5943 :     q = gtou(gel(ai,1));
      88        5936 :     if (q <= X)
      89             :     {
      90        4809 :       long d = ulogint(X, q) + 1;
      91        4809 :       GEN s = direuler_factor(gel(ai,2), d);
      92        4795 :       *n = dirmuleuler_small(V, v, *n, q, s, d);
      93        4795 :       pbad = muliu(pbad, q);
      94             :     }
      95             :   }
      96        3689 :   return pbad;
      97             : }
      98             : 
      99             : GEN
     100         672 : direuler_bad(void *E, GEN (*eval)(void *,GEN,long), GEN a,GEN b,GEN c, GEN Sbad)
     101             : {
     102             :   ulong au, bu, X, sqrtX, n, p;
     103         672 :   pari_sp av0 = avma;
     104             :   GEN gp, v, V;
     105             :   forprime_t T;
     106         672 :   au = direulertou(a, gceil);
     107         665 :   bu = direulertou(b, gfloor);
     108         658 :   X = c ? direulertou(c, gfloor): bu;
     109         651 :   if (X == 0) return cgetg(1,t_VEC);
     110         644 :   if (bu > X) bu = X;
     111         644 :   if (!u_forprime_init(&T, au, bu)) { set_avma(av0); return mkvec(gen_1); }
     112         630 :   v = vecsmall_ei(X, 1);
     113         630 :   V = vec_ei(X, 1);
     114         630 :   n = 1;
     115         630 :   if (Sbad) Sbad = direuler_Sbad(V, v, Sbad, &n);
     116         595 :   p = 1; gp = cgetipos(3); sqrtX = usqrt(X);
     117        8316 :   while (p <= sqrtX && (p = u_forprime_next(&T)))
     118        7742 :     if (!Sbad || umodiu(Sbad, p))
     119             :     {
     120        7637 :       long d = ulogint(X, p) + 1; /* minimal d such that p^d > X */
     121             :       GEN s;
     122        7637 :       gp[2] = p; s = eval(E, gp, d);
     123        7616 :       n = dirmuleuler_small(V, v, n, p, s, d);
     124             :     }
     125      740201 :   while ((p = u_forprime_next(&T))) /* sqrt(X) < p <= X */
     126      739627 :     if (!Sbad || umodiu(Sbad, p))
     127             :     {
     128             :       GEN s;
     129      739620 :       gp[2] = p; s = eval(E, gp, 2); /* s either t_POL or t_SER of val 0 */
     130      739620 :       if (lg(s) > 3 && !gequal0(gel(s,3)))
     131      139153 :         dirmuleuler_large(V, p, gel(s,3));
     132             :     }
     133         574 :   return gerepilecopy(av0,V);
     134             : }
     135             : 
     136             : /* return a t_SER or a truncated t_POL to precision n */
     137             : GEN
     138      752066 : direuler_factor(GEN s, long n)
     139             : {
     140      752066 :   long t = typ(s);
     141      752066 :   if (is_scalar_t(t))
     142             :   {
     143       33194 :     if (!gequal1(s)) err_direuler(s);
     144       33180 :     return scalarpol_shallow(s,0);
     145             :   }
     146      718872 :   switch(t)
     147             :   {
     148        5712 :     case t_POL: break; /* no need to RgXn_red */
     149      712845 :     case t_RFRAC:
     150             :     {
     151      712845 :       GEN p = gel(s,1), q = gel(s,2);
     152      712845 :       q = RgXn_red_shallow(q,n);
     153      712845 :       s = RgXn_inv(q, n);
     154      712845 :       if (typ(p) == t_POL && varn(p) == varn(q))
     155             :       {
     156          28 :         p = RgXn_red_shallow(p, n);
     157          28 :         s = RgXn_mul(s, p, n);
     158             :       }
     159             :       else
     160      712817 :         if (!gequal1(p)) s = RgX_Rg_mul(s, p);
     161      712845 :       if (!signe(s) || !gequal1(gel(s,2))) err_direuler(s);
     162      712831 :       break;
     163             :     }
     164         308 :     case t_SER:
     165         308 :       if (!signe(s) || valser(s) || !gequal1(gel(s,2))) err_direuler(s);
     166         308 :       break;
     167           7 :     default: pari_err_TYPE("direuler", s);
     168             :   }
     169      718851 :   return s;
     170             : }
     171             : 
     172             : struct eval_bad
     173             : {
     174             :   void *E;
     175             :   GEN (*eval)(void *, GEN);
     176             : };
     177             : static GEN
     178      688303 : eval_bad(void *E, GEN p, long n)
     179             : {
     180      688303 :   struct eval_bad *d = (struct eval_bad*) E;
     181      688303 :   return direuler_factor(d->eval(d->E, p), n);
     182             : }
     183             : GEN
     184         301 : direuler(void *E, GEN (*eval)(void *, GEN), GEN a, GEN b, GEN c)
     185             : {
     186             :   struct eval_bad d;
     187         301 :   d.E= E; d.eval = eval;
     188         301 :   return direuler_bad((void*)&d, eval_bad, a, b, c, NULL);
     189             : }
     190             : 
     191             : static GEN
     192       31157 : primelist(forprime_t *T, GEN Sbad, long n, long *running)
     193             : {
     194       31157 :   GEN P = cgetg(n+1, t_VECSMALL);
     195             :   long i, j;
     196      302631 :   for (i = 1, j = 1; i <= n; i++)
     197             :   {
     198      275548 :     ulong p = u_forprime_next(T);
     199      275548 :     if (!p) { *running = 0; break; }
     200      271474 :     if (Sbad && umodiu(Sbad, p)==0) continue;
     201      266791 :     uel(P,j++) = p;
     202             :   }
     203       31157 :   setlg(P, j);
     204       31157 :   return P;
     205             : }
     206             : 
     207             : GEN
     208        4081 : pardireuler(GEN worker, GEN a, GEN b, GEN c, GEN Sbad)
     209             : {
     210             :   ulong au, bu, X, sqrtX, n, snX, nX;
     211        4081 :   pari_sp av0 = avma;
     212             :   GEN v, V;
     213             :   forprime_t T;
     214             :   struct pari_mt pt;
     215        4081 :   long running = 1, pending = 0;
     216        4081 :   au = direulertou(a, gceil);
     217        4081 :   bu = direulertou(b, gfloor);
     218        4081 :   X = c ? direulertou(c, gfloor): bu;
     219        4081 :   if (X == 0) return cgetg(1,t_VEC);
     220        4081 :   if (bu > X) bu = X;
     221        4081 :   if (!u_forprime_init(&T, au, bu)) { set_avma(av0); return mkvec(gen_1); }
     222        4074 :   v = vecsmall_ei(X, 1);
     223        4074 :   V = vec_ei(X, 1);
     224        4074 :   n = 1;
     225        4074 :   if (Sbad) Sbad = direuler_Sbad(V, v, Sbad, &n);
     226        4074 :   sqrtX = usqrt(X); snX = uprimepi(sqrtX); nX = uprimepi(X);
     227        4074 :   if (snX)
     228             :   {
     229        4060 :     GEN P = primelist(&T, Sbad, snX, &running);
     230        4060 :     GEN R = gel(closure_callgenvec(worker, mkvec2(P, utoi(X))), 2);
     231        4060 :     long i, l = lg(P);
     232       20349 :     for (i = 1; i < l; i++)
     233             :     {
     234       16289 :       GEN s = gel(R,i);
     235       16289 :       n = dirmuleuler_small(V, v, n, uel(P,i), s, lg(s));
     236             :     }
     237          14 :   } else snX = 1;
     238        4074 :   mt_queue_start_lim(&pt, worker, (nX+snX-1)/snX);
     239       34212 :   while (running || pending)
     240             :   {
     241             :     GEN done;
     242       30138 :     GEN P = running? primelist(&T, Sbad, snX, &running): NULL;
     243       30138 :     mt_queue_submit(&pt, 0, P ? mkvec2(P, utoi(X)): NULL);
     244       30138 :     done = mt_queue_get(&pt, NULL, &pending);
     245       30138 :     if (done)
     246             :     {
     247       27097 :       GEN P = gel(done,1), R = gel(done,2);
     248       27097 :       long j, l = lg(P);
     249      277599 :       for (j=1; j<l; j++)
     250             :       {
     251      250502 :         GEN F = gel(R,j);
     252      250502 :         if (degpol(F) && !gequal0(gel(F,3)))
     253      169645 :           dirmuleuler_large(V, uel(P,j), gel(F,3));
     254             :       }
     255             :     }
     256             :   }
     257        4074 :   mt_queue_end(&pt);
     258        4074 :   return gerepilecopy(av0,V);
     259             : }
     260             : 
     261             : /********************************************************************/
     262             : /**                                                                **/
     263             : /**                 DIRPOWERS and DIRPOWERSSUM                     **/
     264             : /**                                                                **/
     265             : /********************************************************************/
     266             : 
     267             : /* [1^B,...,N^B] */
     268             : GEN
     269         686 : vecpowuu(long N, ulong B)
     270             : {
     271             :   GEN v;
     272             :   long p, i;
     273             :   forprime_t T;
     274             : 
     275         686 :   if (B <= 8000)
     276             :   {
     277         686 :     if (!B) return const_vec(N,gen_1);
     278         679 :     v = cgetg(N+1, t_VEC); if (N == 0) return v;
     279         679 :     gel(v,1) = gen_1;
     280         679 :     if (B == 1)
     281       92736 :       for (i = 2; i <= N; i++) gel(v,i) = utoipos(i);
     282         469 :     else if (B == 2)
     283             :     {
     284             :       ulong o, s;
     285         273 :       if (N & HIGHMASK)
     286           0 :         for (i = 2, o = 3; i <= N; i++, o += 2)
     287           0 :           gel(v,i) = addiu(gel(v,i-1), o);
     288             :       else
     289       31073 :         for (i = 2, s = 1, o = 3; i <= N; i++, s += o, o += 2)
     290       30800 :           gel(v,i) = utoipos(s + o);
     291             :     }
     292         196 :     else if (B == 3)
     293         840 :       for (i = 2; i <= N; i++) gel(v,i) = powuu(i, B);
     294             :     else
     295             :     {
     296         182 :       long k, Bk, e = expu(N);
     297        7553 :       for (i = 3; i <= N; i += 2) gel(v,i) = powuu(i, B);
     298        1239 :       for (k = 1; k <= e; k++)
     299             :       {
     300        1057 :         N >>= 1; Bk = B * k;
     301        8498 :         for (i = 1; i <= N; i += 2) gel(v, i << k) = shifti(gel(v, i), Bk);
     302             :       }
     303             :     }
     304         679 :     return v;
     305             :   }
     306           0 :   v = const_vec(N, NULL);
     307           0 :   u_forprime_init(&T, 3, N);
     308           0 :   while ((p = u_forprime_next(&T)))
     309             :   {
     310             :     long m, pk, oldpk;
     311           0 :     gel(v,p) = powuu(p, B);
     312           0 :     for (pk = p, oldpk = p; pk; oldpk = pk, pk = umuluu_le(pk,p,N))
     313             :     {
     314           0 :       if (pk != p) gel(v,pk) = mulii(gel(v,oldpk), gel(v,p));
     315           0 :       for (m = N/pk; m > 1; m--)
     316           0 :         if (gel(v,m) && m%p) gel(v, m*pk) = mulii(gel(v,m), gel(v,pk));
     317             :     }
     318             :   }
     319           0 :   gel(v,1) = gen_1;
     320           0 :   for (i = 2; i <= N; i+=2)
     321             :   {
     322           0 :     long vi = vals(i);
     323           0 :     gel(v,i) = shifti(gel(v,i >> vi), B * vi);
     324             :   }
     325           0 :   return v;
     326             : }
     327             : 
     328             : /* does x^s require log(x) ? */
     329             : static long
     330       12913 : get_needlog(GEN s)
     331             : {
     332       12913 :   switch(typ(s))
     333             :   {
     334         294 :     case t_REAL: return 2; /* yes but not powcx */
     335        9485 :     case t_COMPLEX: return 1; /* yes using powcx */
     336        3134 :     default: return 0; /* no */
     337             :   }
     338             : }
     339             : /* [1^B,...,N^B] */
     340             : GEN
     341       12052 : vecpowug(long N, GEN B, long prec)
     342             : {
     343       12052 :   GEN v, logp = NULL;
     344       12052 :   long gp[] = {evaltyp(t_INT)|_evallg(3), evalsigne(1)|evallgefint(3),0};
     345       12052 :   long p, precp = 2, prec0, prec1, needlog;
     346             :   forprime_t T;
     347       12052 :   if (N == 1) return mkvec(gen_1);
     348       12031 :   if (typ(B) == t_INT && lgefint(B) <= 3 && signe(B) >= 0)
     349         168 :     return vecpowuu(N, itou(B));
     350       11863 :   needlog = get_needlog(B);
     351       11863 :   prec1 = prec0 = prec;
     352       11863 :   if (needlog == 1) prec1 = powcx_prec(log2((double)N), B, prec);
     353       11863 :   u_forprime_init(&T, 2, N);
     354       11863 :   v = const_vec(N, NULL);
     355       11863 :   gel(v,1) = gen_1;
     356     1566278 :   while ((p = u_forprime_next(&T)))
     357             :   {
     358             :     long m, pk, oldpk;
     359             :     GEN u;
     360     1554415 :     gp[2] = p;
     361     1554415 :     if (needlog)
     362             :     {
     363       96311 :       if (!logp)
     364       17500 :         logp = logr_abs(utor(p, prec1));
     365             :       else
     366             :       { /* Assuming p and precp are odd,
     367             :          * log p = log(precp) + 2 atanh((p - precp) / (p + precp)) */
     368       78811 :         ulong a = p >> 1, b = precp >> 1; /* p = 2a + 1, precp = 2b + 1 */
     369       78811 :         GEN z = atanhuu(a - b, a + b + 1, prec1); /* avoid overflow */
     370       78811 :         shiftr_inplace(z, 1); logp = addrr(logp, z);
     371             :       }
     372       92787 :       u = needlog == 1? powcx(gp, logp, B, prec0)
     373       96311 :                       : mpexp(gmul(B, logp));
     374       96311 :       if (p == 2) logp = NULL; /* reset: precp must be odd */
     375             :     }
     376             :     else
     377     1458104 :       u = gpow(gp, B, prec0);
     378     1554415 :     precp = p;
     379     1554415 :     gel(v,p) = u; /* p^B */
     380     1554415 :     if (prec0 != prec) gel(v,p) = gprec_wtrunc(gel(v,p), prec);
     381     3222658 :     for (pk = p, oldpk = p; pk; oldpk = pk, pk = umuluu_le(pk,p,N))
     382             :     {
     383     1668243 :       if (pk != p) gel(v,pk) = gmul(gel(v,oldpk), gel(v,p));
     384    46352230 :       for (m = N/pk; m > 1; m--)
     385    44683987 :         if (gel(v,m) && m%p) gel(v, m*pk) = gmul(gel(v,m), gel(v,pk));
     386             :     }
     387             :   }
     388       11863 :   return v;
     389             : }
     390             : 
     391             : GEN
     392         665 : dirpowers(long n, GEN x, long prec)
     393             : {
     394             :   pari_sp av;
     395             :   GEN v;
     396         665 :   if (n <= 0) return cgetg(1, t_VEC);
     397         651 :   av = avma; v = vecpowug(n, x, prec);
     398         651 :   if (typ(x) == t_INT && lgefint(x) <= 3 && signe(x) >= 0 && cmpiu(x, 2) <= 0)
     399         133 :     return v;
     400         518 :   return gerepilecopy(av, v);
     401             : }
     402             : 
     403             : /* f is a totally multiplicative function of modulus 0 or 1
     404             :  * (essentially a Dirichlet character). Compute simultaneously
     405             :  * sum_{0 < n <= N} f(n)n^s and sum_{0 < n <= N} f(n)n^{-1-conj(s)}
     406             :  * Warning: s is conjugated, but not f. Main application for Riemann-Siegel,
     407             :  * where we need R(chi,s) and conj(R(chi,1-conj(s))). */
     408             : 
     409             : static GEN
     410        3234 : vecmulsqlv(GEN Q, GEN V)
     411             : {
     412             :   long l, i;
     413             :   GEN W;
     414        3234 :   if (typ(V) != t_VEC) return RgV_Rg_mul(Q, V);
     415        1239 :   l = lg(Q); W = cgetg(l, t_VEC);
     416       62244 :   for (i = 1; i < l; i++) gel(W, i) = vecmul(gel(Q, i), V);
     417        1239 :   return W;
     418             : }
     419             : 
     420             : /* P = prime divisors of (squarefree) n, V[i] = i^s for i <= sq.
     421             :  * Return NULL if n is not sq-smooth, else f(n)n^s */
     422             : static GEN
     423    18209536 : smallfact(ulong n, GEN P, ulong sq, GEN V)
     424             : {
     425             :   long i, l;
     426             :   ulong p, m, o;
     427             :   GEN c;
     428    18209536 :   if (n <= sq) return gel(V,n);
     429    18160557 :   l = lg(P); m = p = uel(P, l-1); if (p > sq) return NULL;
     430     3846395 :   for (i = l-2; i > 1; i--, m = o) { p = uel(P,i); o = m*p; if (o > sq) break; }
     431     3813089 :   c = gel(V,m); n /= m; /* m <= sq, o = m * p > sq */
     432     3813089 :   if (n > sq) { c = vecmul(c, gel(V,p)); n /= p; }
     433     3805450 :   return vecmul(c, gel(V,n));
     434             : }
     435             : 
     436             : static GEN
     437         973 : _Qtor(GEN x, long prec)
     438         973 : { return typ(x) == t_FRAC? fractor(x, prec): x; }
     439             : static GEN
     440        2457 : Qtor(GEN x, long prec)
     441             : {
     442        2457 :   long tx = typ(x);
     443        3430 :   if (tx == t_VEC || tx == t_COL) pari_APPLY_same(_Qtor(gel(x, i), prec));
     444        1589 :   return tx == t_FRAC? fractor(x, prec): x;
     445             : }
     446             : 
     447             : static GEN
     448         238 : vecf(long N, void *E, GEN (*f)(void *, ulong, long), long prec)
     449             : {
     450             :   GEN v;
     451             :   long n;
     452         238 :   if (!f) return NULL;
     453         140 :   v = cgetg(N + 1, t_VEC);
     454       33481 :   for (n = 1; n <= N; n++) gel(v,n) = f(E, n, prec);
     455         140 :   return v;
     456             : }
     457             : 
     458             : /* Here #V = #F > 0 is small. Analogous to dotproduct, with following
     459             :  * semantic differences: uses V[1] = 1; V has scalar values but F may have
     460             :  * vector values */
     461             : static GEN
     462         301 : naivedirpowerssum(GEN V, GEN F, long prec)
     463             : {
     464             :   GEN S;
     465         301 :   if (!F) S = RgV_sum(V);
     466             :   else
     467             :   {
     468         189 :     long n, l = lg(V);
     469         189 :     S = gel(F,1); /* V[1] = 1 */
     470       34447 :     for (n = 2; n < l; n++) S = gadd(S, gmul(gel(V, n), gel(F, n)));
     471             :   }
     472         301 :   return Qtor(S, prec);
     473             : }
     474             : 
     475             : static GEN
     476         280 : smalldirpowerssum(long N, GEN s, void *E, GEN (*f)(void *, ulong, long),
     477             :                   long both, long prec)
     478             : {
     479             :   GEN F, V, S, SB;
     480         280 :   if (!N)
     481             :   {
     482          42 :     if (!f) return both? mkvec2(gen_0, gen_0): gen_0;
     483          28 :     return gmul(gen_0, f(E, 1, prec));
     484             :   }
     485         238 :   V = vecpowug(N, s, prec);
     486         238 :   F = vecf(N, E, f, prec);
     487         238 :   S = naivedirpowerssum(V, F, prec);
     488         238 :   if (!both) return S;
     489          84 :   if ((both==2 || !f) && gequalm1(gmul2n(real_i(s), 1)))
     490          21 :     SB = S;
     491             :   else
     492             :   {
     493          63 :     GEN FB = (both == 1 && F)? conj_i(F): F;
     494             :     long n;
     495        1876 :     for (n = 2; n <= N; n++) gel(V, n) = ginv(gmulsg(n, gel(V, n)));
     496          63 :     SB = naivedirpowerssum(V, FB, prec);
     497             :   }
     498          84 :   return mkvec2(S, SB);
     499             : }
     500             : 
     501             : static void
     502       68526 : v2unpack(GEN v, GEN *pV, GEN *pVB)
     503             : {
     504       68526 :   if (typ(v) == t_COL) { *pV = gel(v,1); *pVB = gel(v,2); }
     505       68393 :   else { *pV = v; *pVB = NULL; }
     506       68526 : }
     507             : static GEN
     508       69514 : v2pack(GEN V, GEN VB) { return VB? mkcol2(V,VB): V; }
     509             : 
     510             : static GEN
     511        1050 : dirpowsuminit(GEN s, GEN onef, GEN zerf, void *E, GEN (*f)(void *, ulong, long),              GEN data, long both)
     512             : {
     513        1050 :   long needlog = data[1], prec0 = data[2], prec1 = data[3];
     514        1050 :   ulong a, b, c, e, q, n, sq = usqrt(data[4]);
     515        1050 :   GEN V = cgetg(sq+1, t_VEC), W = cgetg(sq+1, t_VEC), VB = NULL, WB = NULL;
     516        1050 :   GEN Q = cgetg(sq+1, t_VEC), Z = cgetg(sq+1, t_VEC), QB = NULL, ZB = NULL;
     517        1050 :   GEN logp, c2, Q2, Q3, Q6, c2B = NULL, Q2B = NULL, Q3B = NULL, Q6B = NULL;
     518        1050 :   long gp[] = {evaltyp(t_INT)|_evallg(3), evalsigne(1)|evallgefint(3),0};
     519             : 
     520        1050 :   if (both == 1 || (both == 2 && !gequal(real_i(s), gneg(ghalf))))
     521          28 :   { VB = cgetg(sq+1, t_VEC); WB = cgetg(sq+1, t_VEC); QB = cgetg(sq+1, t_VEC);}
     522        1050 :   gel(V, 1) = gel(W, 1) = gel(Q, 1) = onef;
     523        1050 :   if (VB) { gel(VB, 1) = gel(WB, 1) = gel(QB, 1) = onef; }
     524        1050 :   c2 = gpow(gen_2, s, prec0); if (VB) c2B = ginv(gmul2n(conj_i(c2), 1));
     525        1050 :   if (f)
     526             :   {
     527        1029 :     GEN tmp2 = f(E, 2, prec0);
     528        1029 :     c2 = gmul(c2, tmp2); if (VB) c2B = gmul(c2B, tmp2);
     529             :   }
     530        1050 :   gel(V,2) = c2; /* f(2) 2^s */
     531        1050 :   gel(W,2) = Qtor(gadd(c2, onef), prec0);
     532        1050 :   gel(Q,2) = Qtor(gadd(vecsqr(c2), onef), prec0);
     533        1050 :   if (VB)
     534             :   {
     535          28 :     gel(VB, 2) = c2B; gel(WB, 2) = Qtor(gadd(c2B, onef), prec0);
     536          28 :     gel(QB, 2) = Qtor(gadd(vecsqr(c2B), onef), prec0);
     537             :   }
     538        1050 :   logp = NULL;
     539      158795 :   for (n = 3; n <= sq; n++)
     540             :   {
     541      157745 :     GEN u = NULL, uB = NULL, ks = f ? f(E, n, prec0) : gen_1;
     542      157745 :     if (gequal0(ks)) ks = NULL;
     543      157745 :     if (odd(n))
     544             :     {
     545       79352 :       gp[2] = n;
     546       79352 :       if (needlog)
     547             :       {
     548       77805 :         if (!logp)
     549        1029 :           logp = logr_abs(utor(n, prec1));
     550             :         else
     551             :         { /* log n = log(n-2) + 2 atanh(1 / (n - 1)) */
     552       76776 :           GEN z = atanhuu(1, n - 1, prec1);
     553       76776 :           shiftr_inplace(z, 1); logp = addrr(logp, z);
     554             :         }
     555       77805 :         if (ks)
     556       76230 :           u = needlog == 1? powcx(gp, logp, s, prec0) : mpexp(gmul(s, logp));
     557             :       }
     558        1547 :       else if (ks) u = gpow(gp, s, prec0);
     559       79352 :       if (ks)
     560             :       {
     561       77763 :         if (VB) uB = gmul(ginv(gmulsg(n, conj_i(u))), ks);
     562       77763 :         u = gmul(u, ks); /* f(n) n^s */
     563             :       }
     564             :     }
     565             :     else
     566             :     {
     567       78393 :       u = vecmul(c2, gel(V, n >> 1));
     568       78393 :       if (VB) uB = vecmul(c2B, gel(VB, n >> 1));
     569             :     }
     570      157745 :     if (ks)
     571             :     {
     572      154693 :       gel(V,n) = u; /* f(n) n^s */
     573      154693 :       gel(W,n) = gadd(gel(W, n-1), gel(V,n));        /*= sum_{i<=n} f(i)i^s */
     574      154693 :       gel(Q,n) = gadd(gel(Q, n-1), vecsqr(gel(V,n)));/*= sum_{i<=n} f(i^2)i^2s*/
     575      154693 :       if (VB)
     576             :       {
     577         483 :         gel(VB,n) = uB;
     578         483 :         gel(WB,n) = gadd(gel(WB,n-1), gel(VB,n));
     579         483 :         gel(QB,n) = gadd(gel(QB,n-1), vecsqr(gel(VB,n)));
     580             :       }
     581             :     }
     582             :     else
     583             :     {
     584        3052 :       gel(V,n) = zerf; gel(W,n) = gel(W, n-1); gel(Q,n) = gel(Q, n-1);
     585        3052 :       if (VB)
     586          21 :       { gel(VB,n) = zerf; gel(WB,n) = gel(WB, n-1); gel(QB,n) = gel(QB, n-1); }
     587             :     }
     588             :   }
     589        1050 :   Q2 = vecmulsqlv(Q, gel(V,2));
     590        1050 :   Q3 = vecmulsqlv(Q, gel(V,3));
     591        1050 :   Q6 = vecmulsqlv(Q, gel(V,6));
     592        1050 :   if (VB)
     593             :   {
     594          28 :     Q2B = vecmulsqlv(QB, gel(VB,2));
     595          28 :     Q3B = vecmulsqlv(QB, gel(VB,3));
     596          28 :     Q6B = vecmulsqlv(QB, gel(VB,6));
     597             :   }
     598             :   /* a,b,c,e = sqrt(q), sqrt(q/2), sqrt(q/3), sqrt(q/6)
     599             :    * Z[q] = Q[a] + 2^s Q[b] + 3^s Q[c] + 6^s Q[e], with Q[0] = 0 */
     600        1050 :   gel(Z, 1) = onef;
     601        1050 :   gel(Z, 2) = gel(W, 2);
     602        1050 :   gel(Z, 3) = gel(W, 3);
     603        1050 :   gel(Z, 4) = gel(Z, 5) = gel(W, 4);
     604        1050 :   gel(Z, 6) = gel(Z, 7) = gadd(gel(W, 4), gel(V, 6));
     605        1050 :   if (VB)
     606             :   {
     607          28 :     ZB = cgetg(sq+1, t_VEC);
     608          28 :     gel(ZB, 1) = onef;
     609          28 :     gel(ZB, 2) = gel(WB, 2);
     610          28 :     gel(ZB, 3) = gel(WB, 3);
     611          28 :     gel(ZB, 4) = gel(ZB, 5) = gel(WB, 4);
     612          28 :     gel(ZB, 6) = gel(ZB, 7) = gadd(gel(WB, 4), gel(VB, 6));
     613             :   }
     614        1050 :   a = 2; b = c = e = 1;
     615      153545 :   for (q = 8; q <= sq; q++)
     616             :   { /* Gray code: at most one of a,b,c,d differs (by 1) from previous value */
     617      152495 :     GEN z = gel(Z, q - 1), zB = NULL;
     618             :     ulong na, nb, nc, ne, na2, nb2, nc2, ne2;
     619      152495 :     if (VB) zB = gel(ZB, q - 1);
     620      152495 :     if ((na = usqrt(q)) != a)
     621        9191 :     { a = na; na2 = na * na; z = gadd(z, gel(V, na2));
     622        9191 :       if (VB) zB = gadd(zB, gel(VB, na2)); }
     623      143304 :     else if ((nb = usqrt(q / 2)) != b)
     624        6216 :     { b = nb; nb2 = 2 * nb * nb; z = gadd(z, gel(V, nb2));
     625        6216 :       if (VB) zB = gadd(zB, gel(VB, nb2)); }
     626      137088 :     else if ((nc = usqrt(q / 3)) != c)
     627        5418 :     { c = nc; nc2 = 3 * nc * nc; z = gadd(z, gel(V, nc2));
     628        5418 :       if (VB) zB = gadd(zB, gel(VB, nc2)); }
     629      131670 :     else if ((ne = usqrt(q / 6)) != e)
     630        3108 :     { e = ne; ne2 = 6 * ne * ne; z = gadd(z, gel(V, ne2));
     631        3108 :       if (VB) zB = gadd(zB, gel(VB, ne2)); }
     632      152495 :     gel(Z, q) = z; if (VB) gel(ZB, q) = zB;
     633             :   }
     634        1078 :   return v2pack(mkvecn(7, V, W, Q, Q2, Q3, Q6, Z),
     635          28 :                 VB? mkvecn(7, VB, WB, QB, Q2B, Q3B, Q6B, ZB): NULL);
     636             : }
     637             : 
     638             : static GEN
     639       38313 : sumprimeloop(forprime_t *pT, GEN s, long N, GEN data, GEN S, GEN W, GEN WB,
     640             :              void *E, GEN (*f)(void *, ulong, long))
     641             : {
     642       38313 :   pari_sp av = avma;
     643       38313 :   long needlog = data[1], prec0 = data[2], prec1 = data[3];
     644       38313 :   long gp[] = {evaltyp(t_INT)|_evallg(3), evalsigne(1)|evallgefint(3),0};
     645       38313 :   ulong p, precp = 0;
     646       38313 :   GEN logp = NULL, SB = WB? S: NULL;
     647             : 
     648     5281003 :   while ((p = u_forprime_next(pT)))
     649             :   {
     650     5242319 :     GEN u = NULL, ks;
     651     5242319 :     if (!f) ks = gen_1;
     652             :     else
     653             :     {
     654     5173432 :       ks = f(E, p, prec1);
     655     5173259 :       if (gequal0(ks)) ks = NULL;
     656     5173656 :       if (ks && gequal1(ks)) ks = gen_1;
     657             :     }
     658     5242519 :     gp[2] = p;
     659     5242519 :     if (needlog)
     660             :     {
     661     5166303 :       if (!logp)
     662       30409 :         logp = logr_abs(utor(p, prec1));
     663             :       else
     664             :       { /* Assuming p and precp are odd,
     665             :          * log p = log(precp) + 2 atanh((p - precp) / (p + precp)) */
     666     5135894 :         ulong a = p >> 1, b = precp >> 1; /* p = 2a + 1, precp = 2b + 1 */
     667     5135894 :         GEN z = atanhuu(a - b, a + b + 1, prec1); /* avoid overflow */
     668     5138688 :         shiftr_inplace(z, 1); logp = addrr(logp, z);
     669             :       }
     670    10331045 :       if (ks) u = needlog == 1? powcx(gp, logp, s, prec0)
     671     5165538 :                               : mpexp(gmul(s, logp));
     672             :     }
     673       76216 :     else if (ks) u = gpow(gp, s, prec0);
     674     5245047 :     if (ks)
     675             :     {
     676     5244998 :       S = gadd(S, vecmul(gel(W, N / p), ks == gen_1? u: gmul(ks, u)));
     677     5242640 :       if (WB)
     678             :       {
     679        2457 :         GEN w = gel(WB, N / p);
     680        2457 :         if (ks != gen_1) w = vecmul(ks, w);
     681        2457 :         SB = gadd(SB, gdiv(w, gmulsg(p, conj_i(u))));
     682             :       }
     683             :     }
     684     5242689 :     precp = p;
     685     5242689 :     if ((p & 0x1ff) == 1)
     686             :     {
     687       19564 :       if (!logp) gerepileall(av, SB? 2: 1, &S, &SB);
     688       19284 :       else gerepileall(av, SB? 3: 2, &S, &logp, &SB);
     689             :     }
     690             :   }
     691       38074 :   return v2pack(S, SB);
     692             : }
     693             : 
     694             : static GEN
     695       48475 : add4(GEN a, GEN b, GEN c, GEN d) { return gadd(gadd(a,b), gadd(c,d)); }
     696             : 
     697             : static const long step = 2048;
     698             : static int
     699       29561 : mksqfloop(long N, long x1, GEN R, GEN RB, GEN *pS, GEN *pSB)
     700             : {
     701       29561 :   GEN V = gel(R,1), Q = gel(R,3), Q2 = gel(R,4);
     702       29561 :   GEN Q3 = gel(R,5), Q6 = gel(R,6), Z = gel(R,7);
     703       29561 :   GEN v, VB = NULL, QB = NULL, Q2B = NULL, Q3B = NULL, Q6B = NULL, ZB = NULL;
     704       29561 :   long x2, j, lv, sq = lg(V)-1;
     705             : 
     706       29561 :   if (RB)
     707          28 :   { VB = gel(RB,1); QB = gel(RB,3); Q2B = gel(RB,4);
     708          28 :     Q3B = gel(RB,5), Q6B = gel(RB,6); ZB = gel(RB,7); }
     709             :   /* beware overflow, fuse last two bins (avoid a tiny remainder) */
     710       29561 :   x2 = (N >= 2*step && N - 2*step >= x1)? x1-1 + step: N;
     711       29561 :   v = vecfactorsquarefreeu_coprime(x1, x2, mkvecsmall2(2, 3));
     712       28767 :   lv = lg(v);
     713    59940816 :   for (j = 1; j < lv; j++)
     714    59911256 :     if (gel(v,j))
     715             :     {
     716    18210168 :       ulong d = x1 - 1 + j; /* squarefree, coprime to 6 */
     717    18210168 :       GEN t = smallfact(d, gel(v,j), sq, V), u;
     718    18182520 :       GEN tB = NULL, uB = NULL; /* = f(d) d^s */
     719             :       long a, b, c, e, q;
     720    18182520 :       if (!t || gequal0(t)) continue;
     721     3795599 :       if (VB) tB = vecinv(gmulsg(d, conj_i(t)));
     722             :       /* warning: gives 1/conj(f(d)) d^(-1-conj(s)), equal to
     723             :          f(d) d^(-1-conj(s)) only if |f(d)|=1. */
     724             :       /* S += f(d) * d^s * Z[q] */
     725     3866633 :       q = N / d;
     726     3866633 :       if (q == 1)
     727             :       {
     728     1474440 :         *pS = gadd(*pS, t); if (VB) *pSB = gadd(*pSB, tB);
     729     1477329 :         continue;
     730             :       }
     731     2392193 :       if (q <= sq) { u = gel(Z, q); if (VB) uB = gel(ZB, q); }
     732             :       else
     733             :       {
     734       48446 :         a = usqrt(q); b = usqrt(q / 2); c = usqrt(q / 3); e = usqrt(q / 6);
     735       48293 :         u = add4(gel(Q,a), gel(Q2,b), gel(Q3,c), gel(Q6,e));
     736       48293 :         if (VB) uB = add4(gel(QB,a), gel(Q2B,b), gel(Q3B,c), gel(Q6B,e));
     737             :       }
     738     2392040 :       *pS = gadd(*pS, vecmul(t, u)); if (VB) *pSB = gadd(*pSB, vecmul(tB, uB));
     739             :     }
     740       29560 :   return x2 == N;
     741             : }
     742             : 
     743             : static GEN
     744        1050 : mkdata(long N, GEN s, long prec)
     745             : {
     746        1050 :   long needlog, prec0, prec1, m = mt_nbthreads(), STEP = maxss(N / (m * m), 1);
     747        1050 :   prec1 = prec0 = prec + EXTRAPRECWORD;
     748        1050 :   needlog = get_needlog(s);
     749        1050 :   if (needlog == 1) prec1 = powcx_prec(log2((double)N), s, prec);
     750        1050 :   return mkvecsmalln(5, needlog, prec0, prec1, N, STEP);
     751             : }
     752             : 
     753             : static GEN
     754       13748 : gp_callUp(void *E, ulong x, long prec)
     755             : {
     756       13748 :   long n[] = {evaltyp(t_INT)|_evallg(3), evalsigne(1)|evallgefint(3),0};
     757       13748 :   n[2] = x; return gp_callprec(E, n, prec);
     758             : }
     759             : 
     760             : /* set *p0 and *p1 to 0 and 1 in the algebra where f takes its values */
     761             : static int
     762        1050 : mk01(void *E, GEN (*f)(void*,ulong,long), long prec, GEN *p0, GEN *p1)
     763             : {
     764        1050 :   *p0 = gen_0; *p1 = gen_1;
     765        1050 :   if (!f) return 1;
     766        1029 :   *p1 = f(E, 1, prec);
     767        1029 :   if (is_vec_t(typ(*p1)))
     768             :   {
     769         406 :     long l = lg(*p1);
     770         406 :     if (l == 1) { *p0 = *p1 = NULL; return 0; }
     771         406 :     *p0 = const_vec(l-1, gen_0);
     772             :   }
     773        1029 :   return 1;
     774             : }
     775             : static GEN
     776           0 : mktrivial(long both)
     777             : {
     778           0 :   if (!both) return cgetg(1, t_VEC);
     779           0 :   return mkvec2(cgetg(1,t_VEC), cgetg(1,t_VEC));
     780             : }
     781             : 
     782             : /* both is
     783             :  * 0: sum_{n<=N}f(n)n^s
     784             :  * 1: sum for (f,s) and (conj(f),-1-s)
     785             :  * 2: sum for (f,s) and (f,-1-s), assuming |f(n)| in {0,1} */
     786             : static GEN
     787         203 : dirpowerssumfun_i(ulong N, GEN s, void *E, GEN (*f)(void *, ulong, long),
     788             :                   long both, long prec)
     789             : {
     790             :   forprime_t T;
     791             :   pari_sp av;
     792             :   GEN onef, zerf, R, RB, W, WB, S, SB, data;
     793             :   ulong x1;
     794             : 
     795         203 :   if ((f && N < 49) || (!f && N < 1000))
     796         140 :     return smalldirpowerssum(N, s, E, f, both, prec);
     797          63 :   if (!mk01(E, f, prec, &zerf, &onef)) return mktrivial(both);
     798          63 :   data = mkdata(N, s, prec); s = gprec_w(s, prec + EXTRAPRECWORD);
     799          63 :   v2unpack(dirpowsuminit(s, onef, zerf, E, f, data, both), &R, &RB);
     800          63 :   W = gel(R,2); WB = RB? gel(RB,2): NULL;
     801          63 :   av = avma; u_forprime_init(&T, lg(W), N);
     802          63 :   v2unpack(sumprimeloop(&T, s, N, data, zerf, W, WB, E, f), &S, &SB);
     803         413 :   for(x1 = 1;; x1 += step)
     804             :   {
     805         413 :     if (mksqfloop(N, x1, R, RB, &S, &SB))
     806          63 :       return both? mkvec2(S, conj_i(SB? SB: S)): S;
     807         350 :     gerepileall(av, SB? 2: 1, &S, &SB);
     808             :   }
     809             : }
     810             : GEN
     811         203 : dirpowerssumfun(ulong N, GEN s, void *E, GEN (*f)(void *, ulong, long),
     812             :                 long both, long prec)
     813             : {
     814         203 :   pari_sp av = avma;
     815         203 :   return gerepilecopy(av, dirpowerssumfun_i(N, s, E, f, both, prec));
     816             : }
     817             : 
     818             : GEN
     819         133 : dirpowerssum(ulong N, GEN s, long both, long prec)
     820         133 : { return dirpowerssumfun(N, s, NULL, NULL, both, prec); }
     821             : 
     822             : GEN
     823         154 : dirpowerssum0(GEN N, GEN s, GEN f, long both, long prec)
     824             : {
     825         154 :   if (typ(N) != t_INT) pari_err_TYPE("dirpowerssum", N);
     826         147 :   if (signe(N) <= 0) N = gen_0;
     827         147 :   if (!f) return dirpowerssum(itou(N), s, both, prec);
     828          70 :   if (typ(f) != t_CLOSURE) pari_err_TYPE("dirpowerssum", f);
     829          70 :   return dirpowerssumfun(itou(N), s, (void*)f, gp_callUp, both, prec);
     830             : }
     831             : 
     832             : GEN
     833       29148 : parsqf_worker(GEN gk, GEN vR, GEN gN)
     834             : {
     835       29148 :   pari_sp av = avma;
     836             :   GEN R, RB, onef, S, SB;
     837       29148 :   long k = itou(gk), N = itou(gN), x1 = 1 + step * k;
     838       29148 :   v2unpack(vR, &R, &RB); onef = gmael(R,1,1);
     839       29148 :   S = SB = is_vec_t(typ(onef)) ? zerovec(lg(onef) - 1): gen_0;
     840       29148 :   (void)mksqfloop(N, x1, R, RB, &S, &SB);
     841       29147 :   return gerepilecopy(av, v2pack(S, RB? SB: NULL));
     842             : }
     843             : 
     844             : static GEN
     845     5350137 : mycallvec(void *f, ulong n, long prec)
     846             : {
     847     5350137 :   GEN F = (GEN)f;
     848     5350137 :   if (!f) return gen_1;
     849      390617 :   if (typ(F) == t_CLOSURE) return gp_callUp(f, n, prec);
     850      390617 :   return gel(F, (n-1) % (lg(F)-1) + 1);
     851             : }
     852             : 
     853             : GEN
     854       38274 : parsumprimefun_worker(GEN gk, GEN s, GEN zerf, GEN data, GEN vW, GEN f)
     855             : {
     856       38274 :   pari_sp av = avma;
     857             :   forprime_t T;
     858             :   GEN W, WB, S;
     859       38274 :   long k = itou(gk), sq, N = data[4], STEP = data[5];
     860             : 
     861       38272 :   v2unpack(vW, &W, &WB); sq = lg(W)-1;
     862       38265 :   if (isintzero(f)) f = NULL;
     863       38265 :   u_forprime_init(&T, k * STEP + sq + 1, minss(N, (k + 1) * STEP + sq));
     864       38246 :   S = sumprimeloop(&T, s, N, data, zerf, W, WB, (void*)f, mycallvec);
     865       38267 :   return gerepilecopy(av, S);
     866             : }
     867             : 
     868             : static GEN
     869         987 : vR_get_vW(GEN vR)
     870             : {
     871             :   GEN R, RB, W, WB;
     872         987 :   v2unpack(vR, &R, &RB); W = gel(R,2); WB = RB? gel(RB,2): NULL;
     873         987 :   return v2pack(W, WB);
     874             : }
     875             : 
     876             : static GEN
     877         987 : halfconj(long both, GEN V)
     878         987 : { return both ? mkvec2(gel(V, 1), gconj(gel(V, 2))) : V; }
     879             : 
     880             : static GEN
     881        1127 : pardirpowerssumfun_i(GEN f, ulong N, GEN s, long both, long prec)
     882             : {
     883             :   GEN worker, worker2, data, vR, onef, zerf;
     884             : 
     885        1127 :   if ((f && N < 49) || (!f && N < 10000UL))
     886         140 :     return smalldirpowerssum(N, s, (void*)f, mycallvec, both, prec);
     887         987 :   if (!mk01((void*)f, mycallvec, prec, &zerf, &onef)) return mktrivial(both);
     888         987 :   data = mkdata(N, s, prec); s = gprec_w(s, prec + EXTRAPRECWORD);
     889         987 :   vR = dirpowsuminit(s, onef, zerf, (void*)f, mycallvec, data, both);
     890         987 :   worker = snm_closure(is_entry("_parsumprimefun_worker"),
     891             :                        mkvecn(5, s, zerf, data, vR_get_vW(vR), f? f: gen_0));
     892         987 :   worker2 = snm_closure(is_entry("_parsqf_worker"), mkvec2(vR, utoi(N)));
     893        1974 :   return halfconj(both,
     894         987 :                   gadd(parsum(gen_0, utoipos((N-1) / data[5]), worker),
     895         987 :                        parsum(gen_0, utoipos(maxss((N-1) / step - 1, 0)), worker2)));
     896             : }
     897             : 
     898             : GEN
     899        1127 : pardirpowerssumfun(GEN f, ulong N, GEN s, long both, long prec)
     900             : {
     901        1127 :   pari_sp av = avma;
     902        1127 :   return gerepilecopy(av, pardirpowerssumfun_i(f, N, s, both, prec));
     903             : }
     904             : GEN
     905           0 : pardirpowerssum0(GEN N, GEN s, GEN f, long both, long prec)
     906             : {
     907           0 :   if (typ(N) != t_INT) pari_err_TYPE("pardirpowerssum", N);
     908           0 :   return pardirpowerssumfun(f, itou(N), s, both, prec);
     909             : }
     910             : GEN
     911           0 : pardirpowerssum(ulong N, GEN s, long prec)
     912           0 : { return pardirpowerssumfun(NULL, N, s, 0, prec); }

Generated by: LCOV version 1.16