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 n^s require log(x) ? */
329 : static long
330 50999 : get_needlog(GEN s)
331 : {
332 50999 : switch(typ(s))
333 : {
334 294 : case t_REAL: return 2; /* yes but not powcx */
335 47571 : case t_COMPLEX: return 1; /* yes using powcx */
336 3134 : default: return 0; /* no */
337 : }
338 : }
339 : /* [1^B,...,N^B] */
340 : GEN
341 12024 : vecpowug(long N, GEN B, long prec)
342 : {
343 12024 : GEN v, logp = NULL;
344 12024 : long gp[] = {evaltyp(t_INT)|_evallg(3), evalsigne(1)|evallgefint(3),0};
345 12024 : long p, precp = 2, prec0, prec1, needlog;
346 : forprime_t T;
347 12024 : if (N == 1) return mkvec(gen_1);
348 12003 : if (typ(B) == t_INT && lgefint(B) <= 3 && signe(B) >= 0)
349 168 : return vecpowuu(N, itou(B));
350 11835 : needlog = get_needlog(B);
351 11835 : prec1 = prec0 = prec;
352 11835 : if (needlog == 1) prec1 = powcx_prec(log2((double)N), B, prec);
353 11835 : u_forprime_init(&T, 2, N);
354 11835 : v = const_vec(N, NULL);
355 11835 : gel(v,1) = gen_1;
356 1565942 : while ((p = u_forprime_next(&T)))
357 : {
358 : long m, pk, oldpk;
359 : GEN u;
360 1554107 : gp[2] = p;
361 1554107 : if (needlog)
362 : {
363 96003 : if (!logp)
364 17444 : 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 78559 : ulong a = p >> 1, b = precp >> 1; /* p = 2a + 1, precp = 2b + 1 */
369 78559 : GEN z = atanhuu(a - b, a + b + 1, prec1); /* avoid overflow */
370 78559 : shiftr_inplace(z, 1); logp = addrr(logp, z);
371 : }
372 92479 : u = needlog == 1? powcx(gp, logp, B, prec0)
373 96003 : : mpexp(gmul(B, logp));
374 96003 : if (p == 2) logp = NULL; /* reset: precp must be odd */
375 : }
376 : else
377 1458104 : u = gpow(gp, B, prec0);
378 1554107 : precp = p;
379 1554107 : gel(v,p) = u; /* p^B */
380 1554107 : if (prec0 != prec) gel(v,p) = gprec_wtrunc(gel(v,p), prec);
381 3221888 : for (pk = p, oldpk = p; pk; oldpk = pk, pk = umuluu_le(pk,p,N))
382 : {
383 1667781 : if (pk != p) gel(v,pk) = gmul(gel(v,oldpk), gel(v,p));
384 46350452 : for (m = N/pk; m > 1; m--)
385 44682671 : if (gel(v,m) && m%p) gel(v, m*pk) = gmul(gel(v,m), gel(v,pk));
386 : }
387 : }
388 11835 : 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 : static GEN
404 252 : vecmulsqlv(GEN Q, GEN V)
405 : {
406 : long lq, i;
407 : GEN W;
408 252 : if (typ(V) != t_VEC) return RgV_Rg_mul(Q, V);
409 21 : lq = lg(Q); W = cgetg(lq, t_VEC);
410 672 : for (i = 1; i < lq; i++) gel(W, i) = vecmul(gel(Q, i), V);
411 21 : return W;
412 : }
413 :
414 : /* P = prime divisors of (squarefree) n, V[i] = i^s for i <= sq.
415 : * Return NULL if n is not sq-smooth, else f(n)n^s */
416 : static GEN
417 18241038 : smallfact(ulong n, GEN P, ulong sq, GEN V)
418 : {
419 : long i, l;
420 : ulong p, m, o;
421 : GEN c;
422 18241038 : if (n <= sq) return gel(V,n);
423 18192101 : l = lg(P); m = p = uel(P, l-1); if (p > sq) return NULL;
424 3855968 : for (i = l-2; i > 1; i--, m = o) { p = uel(P,i); o = m*p; if (o > sq) break; }
425 3822657 : c = gel(V,m); n /= m; /* m <= sq, o = m * p > sq */
426 3822657 : if (n > sq) { c = vecmul(c, gel(V,p)); n /= p; }
427 3813245 : return vecmul(c, gel(V,n));
428 : }
429 :
430 : static GEN
431 56 : _Qtor(GEN x, long prec)
432 56 : { return typ(x) == t_FRAC? fractor(x, prec): x; }
433 : static GEN
434 427 : Qtor(GEN x, long prec)
435 : {
436 427 : long tx = typ(x);
437 483 : if (tx == t_VEC || tx == t_COL) pari_APPLY_same(_Qtor(gel(x, i), prec));
438 392 : return tx == t_FRAC? fractor(x, prec): x;
439 : }
440 :
441 : static GEN
442 210 : vecf(long N, void *E, GEN (*f)(void *, ulong, long), long prec)
443 : {
444 : GEN v;
445 : long n;
446 210 : if (!f) return NULL;
447 112 : v = cgetg(N + 1, t_VEC);
448 32557 : for (n = 1; n <= N; n++) gel(v,n) = f(E, n, prec);
449 112 : return v;
450 : }
451 :
452 : /* Here #V = #F > 0 is small. Analogous to dotproduct, with following
453 : * semantic differences: uses V[1] = 1; V has scalar values but F may have
454 : * vector values */
455 : static GEN
456 259 : naivedirpowerssum(GEN V, GEN F, long prec)
457 : {
458 : GEN S;
459 259 : if (!F) S = RgV_sum(V);
460 : else
461 : {
462 147 : long n, l = lg(V);
463 147 : S = gel(F,1); /* V[1] = 1 */
464 33103 : for (n = 2; n < l; n++) S = gadd(S, gmul(gel(V, n), gel(F, n)));
465 : }
466 259 : return Qtor(S, prec);
467 : }
468 :
469 : static GEN
470 252 : smalldirpowerssum(long N, GEN s, void *E, GEN (*f)(void *, ulong, long),
471 : long both, long prec)
472 : {
473 : GEN F, V, S, SB, sb;
474 252 : if (!N)
475 : {
476 42 : if (!f) return both? mkvec2(gen_0, gen_0): gen_0;
477 28 : return gmul(gen_0, f(E, 1, prec));
478 : }
479 210 : V = vecpowug(N, s, prec);
480 210 : F = vecf(N, E, f, prec);
481 210 : S = naivedirpowerssum(V, F, prec);
482 210 : if (!both) return S;
483 :
484 70 : sb = conj_i(gsubsg(-1, s));
485 70 : if ((both==2 || !f) && gequal(s,sb))
486 21 : SB = S;
487 : else
488 : {
489 49 : GEN FB = (both == 1 && F)? conj_i(F): F;
490 49 : GEN VB = cgetg(N+1, t_VEC);
491 : long n;
492 49 : gel(VB, 1) = gel(V, 1); /* = 1 */
493 1428 : for (n = 2; n <= N; n++) gel(VB, n) = ginv(gmulsg(n, gel(V, n)));
494 49 : SB = naivedirpowerssum(VB, FB, prec);
495 : }
496 70 : return mkvec2(S, SB);
497 : }
498 :
499 : static void
500 87581 : v2unpack(GEN v, GEN *pV, GEN *pVB)
501 : {
502 87581 : if (typ(v) == t_COL) { *pV = gel(v,1); *pVB = gel(v,2); }
503 87476 : else { *pV = v; *pVB = NULL; }
504 87581 : }
505 : static GEN
506 71345 : v2pack(GEN V, GEN VB) { return VB? mkcol2(V,VB): V; }
507 :
508 : static GEN
509 63 : dirpowsuminit(GEN s, void *E, GEN (*f)(void *, ulong, long), GEN data,
510 : long both, long prec)
511 : {
512 63 : GEN onef = gel(data, 1), zervec = gel(data, 2), sqlpp = gel(data, 3);
513 63 : long n, sq = sqlpp[1], needlog = sqlpp[2], prec0 = sqlpp[3], prec1 = sqlpp[4];
514 63 : GEN V = cgetg(sq+1, t_VEC), W = cgetg(sq+1, t_VEC), Q = cgetg(sq+1, t_VEC);
515 63 : GEN VB = NULL, WB = NULL, QB = NULL;
516 63 : GEN c2, Q2, Q3, Q6, c2B = NULL, Q2B = NULL, Q3B = NULL, Q6B = NULL;
517 : GEN logp;
518 63 : long gp[] = {evaltyp(t_INT)|_evallg(3), evalsigne(1)|evallgefint(3),0};
519 :
520 63 : if (both == 1 || (both == 2 && !gequal(real_i(s), gneg(ghalf))))
521 21 : { VB = cgetg(sq+1, t_VEC); WB = cgetg(sq+1, t_VEC); QB = cgetg(sq+1, t_VEC);}
522 63 : gel(V, 1) = gel(W, 1) = gel(Q, 1) = onef;
523 63 : if (VB) { gel(VB, 1) = gel(WB, 1) = gel(QB, 1) = onef; }
524 63 : c2 = gpow(gen_2, s, prec0); if (VB) c2B = ginv(gmul2n(conj_i(c2), 1));
525 63 : if (f)
526 : {
527 42 : GEN tmp2 = f(E, 2, prec);
528 42 : c2 = gmul(c2, tmp2); if (VB) c2B = gmul(c2B, tmp2);
529 : }
530 63 : gel(V,2) = c2; /* f(2) 2^s */
531 63 : gel(W,2) = Qtor(gadd(c2, onef), prec0);
532 63 : gel(Q,2) = Qtor(gadd(vecsqr(c2), onef), prec0);
533 63 : if (VB)
534 : {
535 21 : gel(VB, 2) = c2B; gel(WB, 2) = Qtor(gadd(c2B, onef), prec0);
536 21 : gel(QB, 2) = Qtor(gadd(vecsqr(c2B), onef), prec0);
537 : }
538 63 : logp = NULL;
539 4074 : for (n = 3; n <= sq; n++)
540 : {
541 4011 : GEN u = NULL, uB = NULL, ks = f ? f(E, n, prec0) : gen_1;
542 4011 : long zks = !gequal0(ks);
543 4011 : if (odd(n))
544 : {
545 2023 : gp[2] = n;
546 2023 : if (needlog)
547 : {
548 476 : if (!logp)
549 42 : logp = logr_abs(utor(n, prec1));
550 : else
551 : { /* log n = log(n-2) + 2 atanh(1 / (n - 1)) */
552 434 : GEN z = atanhuu(1, n - 1, prec1);
553 434 : shiftr_inplace(z, 1); logp = addrr(logp, z);
554 : }
555 476 : if (zks)
556 476 : u = needlog == 1? powcx(gp, logp, s, prec0) : mpexp(gmul(s, logp));
557 : }
558 1547 : else if (zks) u = gpow(gp, s, prec0);
559 2023 : if (zks)
560 : {
561 2009 : if (VB) uB = gmul(ginv(gmulsg(n, conj_i(u))), ks);
562 2009 : u = gmul(u, ks); /* f(n) n^s */
563 : }
564 : }
565 : else
566 : {
567 1988 : u = vecmul(c2, gel(V, n >> 1));
568 1988 : if (VB) uB = vecmul(c2B, gel(VB, n >> 1));
569 : }
570 4011 : if (zks)
571 : { /* V[n]=f(n)n^s, W[n]=sum_{i<=n} f(i)i^s, Q[n]=sum_{i<=n} f(i^2)i^2s */
572 3983 : gel(V,n) = u;
573 3983 : gel(W,n) = gadd(gel(W, n-1), gel(V,n));
574 3983 : gel(Q,n) = gadd(gel(Q, n-1), vecsqr(gel(V,n)));
575 3983 : if (VB)
576 : {
577 462 : gel(VB,n) = uB;
578 462 : gel(WB,n) = gadd(gel(WB,n-1), gel(VB,n));
579 462 : gel(QB,n) = gadd(gel(QB,n-1), vecsqr(gel(VB,n)));
580 : }
581 : }
582 : else
583 : {
584 28 : gel(V,n) = zervec; gel(W,n) = gel(W, n-1); gel(Q,n) = gel(Q, n-1);
585 28 : if (VB)
586 : {
587 0 : gel(VB,n) = zervec; gel(WB,n) = gel(WB, n-1);
588 0 : gel(QB,n) = gel(QB, n-1);
589 : }
590 : }
591 : }
592 63 : Q2 = vecmulsqlv(Q, gel(V,2));
593 63 : Q3 = vecmulsqlv(Q, gel(V,3));
594 63 : Q6 = vecmulsqlv(Q, gel(V,6));
595 63 : if (VB)
596 : {
597 21 : Q2B = vecmulsqlv(QB, gel(VB,2));
598 21 : Q3B = vecmulsqlv(QB, gel(VB,3));
599 21 : Q6B = vecmulsqlv(QB, gel(VB,6));
600 : }
601 84 : return v2pack(mkvecn(6, V, W, Q, Q2, Q3, Q6),
602 21 : VB? mkvecn(6, VB, WB, QB, Q2B, Q3B, Q6B): NULL);
603 : }
604 :
605 : static GEN
606 63 : dirpowsumprimeloop(ulong N, GEN s, void *E, GEN (*f)(void *, ulong, long),
607 : GEN data, GEN W, GEN WB)
608 : {
609 : pari_sp av2;
610 63 : GEN zervec = gel(data, 2), S = zervec, SB = zervec, logp = NULL;
611 63 : GEN sqlpp = gel(data, 3);
612 : forprime_t T;
613 63 : long gp[] = {evaltyp(t_INT)|_evallg(3), evalsigne(1)|evallgefint(3),0};
614 63 : long p, precp = 0, sq = sqlpp[1], needlog = sqlpp[2];
615 63 : long prec0 = sqlpp[3], prec1 = sqlpp[4];
616 63 : u_forprime_init(&T, sq + 1, N);
617 63 : av2 = avma;
618 80969 : while ((p = u_forprime_next(&T)))
619 : {
620 80906 : GEN u = NULL, ks = f ? f(E, p, prec1) : gen_1;
621 80906 : long zks = !gequal0(ks);
622 80906 : gp[2] = p;
623 80906 : if (needlog)
624 : {
625 4690 : if (!logp)
626 42 : logp = logr_abs(utor(p, prec1));
627 : else
628 : { /* log p = log(precp) + 2 atanh((p - precp) / (p + precp)) */
629 4648 : ulong a = p >> 1, b = precp >> 1; /* p = 2a + 1, precp = 2b + 1 */
630 4648 : GEN z = atanhuu(a - b, a + b + 1, prec1); /* avoid overflow */
631 4648 : shiftr_inplace(z, 1); logp = addrr(logp, z);
632 : }
633 4690 : if (zks)
634 4690 : u = needlog == 1? powcx(gp, logp, s, prec0) : mpexp(gmul(s, logp));
635 : }
636 76216 : else { if (zks) u = gpow(gp, s, prec0); }
637 80906 : if (zks)
638 : {
639 80906 : S = gadd(S, vecmul(gel(W, N / p), gmul(ks, u)));
640 80906 : if (WB)
641 2345 : SB = gadd(SB, gdiv(vecmul(ks, gel(WB, N / p)), gmulsg(p, conj_i(u))));
642 : }
643 80906 : precp = p;
644 80906 : if ((p & 0x1ff) == 1)
645 : {
646 280 : if (!logp) gerepileall(av2, SB? 2: 1, &S, &SB);
647 0 : else gerepileall(av2, SB? 3: 2, &S, &logp, &SB);
648 : }
649 : }
650 63 : return v2pack(S, SB);
651 : }
652 :
653 : static GEN
654 1435 : add4(GEN a, GEN b, GEN c, GEN d) { return gadd(gadd(a,b), gadd(c,d)); }
655 :
656 : static void
657 413 : dirpowsumsqfloop(long N, long x1, long x2, long sq, GEN P,
658 : GEN *pS, GEN Z, GEN V, GEN Q, GEN Q2, GEN Q3, GEN Q6,
659 : GEN *pSB, GEN ZB, GEN VB, GEN QB, GEN Q2B, GEN Q3B, GEN Q6B)
660 : {
661 413 : GEN v = vecfactorsquarefreeu_coprime(x1, x2, P);
662 413 : long lv = lg(v), j;
663 806813 : for (j = 1; j < lv; j++)
664 806400 : if (gel(v,j))
665 : {
666 245126 : ulong d = x1 - 1 + j; /* squarefree, coprime to 6 */
667 245126 : GEN t = smallfact(d, gel(v,j), sq, V), u;
668 245126 : GEN tB = NULL, uB = NULL; /* = d^s */
669 : long a, b, c, e, q;
670 245126 : if (!t || gequal0(t)) continue;
671 48748 : if (VB) tB = vecinv(gmulsg(d, conj_i(t)));
672 : /* warning: gives 1/conj(f(d)) d^(-1-conj(s)), equal to
673 : f(d) d^(-1-conj(s)) only if |f(d)|=1. */
674 : /* S += f(d) * d^s * Z[q] */
675 48748 : q = N / d;
676 48748 : if (q == 1)
677 : {
678 17339 : *pS = gadd(*pS, t); if (VB) *pSB = gadd(*pSB, tB);
679 17339 : continue;
680 : }
681 31409 : if (q <= sq) { u = gel(Z, q); if (VB) uB = gel(ZB, q); }
682 : else
683 : {
684 1274 : a = usqrt(q); b = usqrt(q / 2); c = usqrt(q / 3); e = usqrt(q / 6);
685 1274 : u = add4(gel(Q,a), gel(Q2,b), gel(Q3,c), gel(Q6,e));
686 1274 : if (VB) uB = add4(gel(QB,a), gel(Q2B,b), gel(Q3B,c), gel(Q6B,e));
687 : }
688 31409 : *pS = gadd(*pS, vecmul(t, u)); if (VB) *pSB = gadd(*pSB, vecmul(tB, uB));
689 : }
690 413 : }
691 :
692 : static GEN
693 63 : dirpowsummakez(GEN V, GEN W, GEN VB, GEN WB, GEN onef, ulong sq)
694 : {
695 63 : GEN Z = cgetg(sq+1, t_VEC), ZB = NULL;
696 : ulong a, b, c, e, q;
697 : /* a,b,c,e = sqrt(q), sqrt(q/2), sqrt(q/3), sqrt(q/6)
698 : * Z[q] = Q[a] + 2^s Q[b] + 3^s Q[c] + 6^s Q[e], with Q[0] = 0 */
699 63 : gel(Z, 1) = onef;
700 63 : gel(Z, 2) = gel(W, 2);
701 63 : gel(Z, 3) = gel(W, 3);
702 63 : gel(Z, 4) = gel(Z, 5) = gel(W, 4);
703 63 : gel(Z, 6) = gel(Z, 7) = gadd(gel(W, 4), gel(V, 6));
704 63 : if (VB)
705 : {
706 21 : ZB = cgetg(sq+1, t_VEC);
707 21 : gel(ZB, 1) = onef;
708 21 : gel(ZB, 2) = gel(WB, 2);
709 21 : gel(ZB, 3) = gel(WB, 3);
710 21 : gel(ZB, 4) = gel(ZB, 5) = gel(WB, 4);
711 21 : gel(ZB, 6) = gel(ZB, 7) = gadd(gel(WB, 4), gel(VB, 6));
712 : }
713 63 : a = 2; b = c = e = 1;
714 3759 : for (q = 8; q <= sq; q++)
715 : { /* Gray code: at most one of a,b,c,d differs (by 1) from previous value */
716 3696 : GEN z = gel(Z, q - 1), zB = NULL;
717 : ulong na, nb, nc, ne, na2, nb2, nc2, ne2;
718 3696 : if (VB) zB = gel(ZB, q - 1);
719 3696 : if ((na = usqrt(q)) != a)
720 280 : { a = na; na2 = na * na; z = gadd(z, gel(V, na2));
721 280 : if (VB) zB = gadd(zB, gel(VB, na2)); }
722 3416 : else if ((nb = usqrt(q / 2)) != b)
723 203 : { b = nb; nb2 = 2 * nb * nb; z = gadd(z, gel(V, nb2));
724 203 : if (VB) zB = gadd(zB, gel(VB, nb2)); }
725 3213 : else if ((nc = usqrt(q / 3)) != c)
726 161 : { c = nc; nc2 = 3 * nc * nc; z = gadd(z, gel(V, nc2));
727 161 : if (VB) zB = gadd(zB, gel(VB, nc2)); }
728 3052 : else if ((ne = usqrt(q / 6)) != e)
729 98 : { e = ne; ne2 = 6 * ne * ne; z = gadd(z, gel(V, ne2));
730 98 : if (VB) zB = gadd(zB, gel(VB, ne2)); }
731 3696 : gel(Z, q) = z; if (VB) gel(ZB, q) = zB;
732 : }
733 63 : return v2pack(Z, ZB);
734 : }
735 :
736 : static const long step = 2048;
737 :
738 : /* both =
739 : * 0: sum_{n<=N}f(n)n^s
740 : * 1: sum for (f,s) and (conj(f),-1-s)
741 : * 2: sum for (f,s) and (f,-1-s), assuming |f(n)| in {0,1} */
742 : static GEN
743 203 : dirpowerssumfun_i(ulong N, GEN s, void *E, GEN (*f)(void *, ulong, long),
744 : long both, long prec)
745 : {
746 : pari_sp av;
747 : GEN P, V, W, Q, Q2, Q3, Q6, S, Z, onef, zervec;
748 203 : GEN VB = NULL, WB = NULL, QB = NULL;
749 203 : GEN Q2B = NULL, Q3B = NULL, Q6B = NULL, SB = NULL, ZB = NULL;
750 : GEN R, RB, data;
751 : ulong x1, sq;
752 : long prec0, prec1, needlog;
753 :
754 203 : if ((f && N < 49) || (!f && N < 1000))
755 140 : return smalldirpowerssum(N, s, E, f, both, prec);
756 63 : onef = f ? f(E, 1, prec) : gen_1;
757 63 : zervec = gmul(gen_0, onef);
758 63 : sq = usqrt(N);
759 63 : prec1 = prec0 = prec + EXTRAPREC64;
760 63 : s = gprec_w(s, prec0);
761 63 : needlog = get_needlog(s);
762 63 : if (needlog == 1) prec1 = powcx_prec(log2((double)N), s, prec);
763 63 : data = mkvec3(onef, zervec, mkvecsmall4(sq, needlog, prec0, prec1));
764 63 : v2unpack(dirpowsuminit(s, E, f, data, both, prec), &R, &RB);
765 63 : V = gel(R, 1); W = gel(R, 2); Q = gel(R, 3);
766 63 : Q2 = gel(R, 4); Q3 = gel(R, 5); Q6 = gel(R, 6);
767 63 : if (RB)
768 : {
769 21 : VB = gel(RB, 1); WB = gel(RB, 2); QB = gel(RB, 3);
770 21 : Q2B = gel(RB, 4); Q3B = gel(RB, 5); Q6B = gel(RB, 6);
771 : }
772 63 : v2unpack(dirpowsumprimeloop(N, s, E, f, data, W, WB), &S, &SB);
773 63 : v2unpack(dirpowsummakez(V, W, VB, WB, onef, sq), &Z, &ZB);
774 63 : P = mkvecsmall2(2, 3); av = avma;
775 63 : for(x1 = 1;; x1 += step)
776 350 : { /* beware overflow, fuse last two bins (avoid a tiny remainder) */
777 413 : ulong x2 = (N >= 2*step && N - 2*step >= x1)? x1-1 + step: N;
778 413 : dirpowsumsqfloop(N, x1, x2, sq, P, &S, Z, V, Q, Q2, Q3, Q6,
779 : &SB, ZB, VB, QB, Q2B, Q3B, Q6B);
780 413 : if (x2 == N) break;
781 350 : gerepileall(av, SB? 2: 1, &S, &SB);
782 : }
783 63 : return both? mkvec2(S, conj_i(VB? SB: S)): S;
784 : }
785 : GEN
786 203 : dirpowerssumfun(ulong N, GEN s, void *E, GEN (*f)(void *, ulong, long),
787 : long both, long prec)
788 : {
789 203 : pari_sp av = avma;
790 203 : return gerepilecopy(av, dirpowerssumfun_i(N, s, E, f, both, prec));
791 : }
792 :
793 : GEN
794 133 : dirpowerssum(ulong N, GEN s, long both, long prec)
795 133 : { return dirpowerssumfun(N, s, NULL, NULL, both, prec); }
796 :
797 : static GEN
798 13748 : gp_callUp(void *E, ulong x, long prec)
799 : {
800 13748 : long court[] = {evaltyp(t_INT)|_evallg(3), evalsigne(1)|evallgefint(3),0};
801 13748 : court[2] = x; return gp_callprec(E, court, prec);
802 : }
803 :
804 : GEN
805 154 : dirpowerssum0(GEN N, GEN s, GEN f, long both, long prec)
806 : {
807 154 : if (typ(N) != t_INT) pari_err_TYPE("dirpowerssum", N);
808 147 : if (signe(N) <= 0) N = gen_0;
809 147 : if (!f) return dirpowerssum(itou(N), s, both, prec);
810 70 : if (typ(f) != t_CLOSURE) pari_err_TYPE("dirpowerssum", f);
811 70 : return dirpowerssumfun(itou(N), s, (void*)f, gp_callUp, both, prec);
812 : }
813 :
814 : /*******************************************************************/
815 : /* Parallel dirpowerssumfun */
816 : /*******************************************************************/
817 : /* f is a totally multiplicative function of modulus 0 or 1
818 : * (essentially a Dirichlet character). Compute simultaneously
819 : * sum_{0 < n <= N} f(n)n^s and sum_{0 < n <= N} f(n)n^{-1-conj(s)}
820 : * Warning: s is conjugated, but not f. Main application for Riemann-Siegel,
821 : * where we need R(chi,s) and conj(R(chi,1-conj(s))). */
822 :
823 : static GEN
824 5389941 : mycallvec(void *f, ulong n, long prec)
825 : {
826 5389941 : GEN F = (GEN)f;
827 5389941 : if (!f) return gen_1;
828 405299 : if (typ(F) == t_CLOSURE) return gp_callUp(f, n, prec);
829 405299 : return gel(F, (n-1) % (lg(F)-1) + 1);
830 : }
831 :
832 : static GEN
833 2919 : gmulvecsqlv(GEN Q, GEN V)
834 : {
835 : long lq, i;
836 : GEN W;
837 2919 : if (typ(V) != t_VEC) return RgV_Rg_mul(Q, V);
838 1155 : lq = lg(Q); W = cgetg(lq, t_VEC);
839 61005 : for (i = 1; i < lq; i++) gel(W, i) = vecmul(gel(Q, i), V);
840 1155 : return W;
841 : }
842 :
843 : GEN
844 29132 : parsqfboth_worker(GEN gk, GEN vZ, GEN vVQ, GEN vV, GEN P, GEN Nsq)
845 : {
846 29132 : pari_sp av = avma;
847 : GEN Z, ZB, V, VB, vQ, Q, Q2, Q3, Q6, vQB, QB, Q2B, Q3B, Q6B, v, S, SB;
848 29132 : long k = itos(gk), N = Nsq[1];
849 29132 : long x1 = 1 + step * k, x2, j, lv;
850 29132 : ulong sq = Nsq[2];
851 29132 : v2unpack(vZ, &Z, &ZB);
852 29131 : v2unpack(vV, &V, &VB);
853 29131 : v2unpack(vVQ, &vQ, &vQB);
854 29130 : S = SB = is_vec_t(typ(gel(V,1))) ? zerovec(lg(gel(V,1)) - 1): gen_0;
855 29130 : Q = gel(vQ,1); Q2 = gel(vQ,2); Q3 = gel(vQ,3); Q6 = gel(vQ,4);
856 29130 : if (vQB)
857 0 : { QB = gel(vQB,1); Q2B = gel(vQB,2); Q3B = gel(vQB,3); Q6B = gel(vQB,4); }
858 : else
859 29130 : QB = Q2B = Q3B = Q6B = SB = NULL;
860 : /* beware overflow, fuse last two bins (avoid a tiny remainder) */
861 29130 : x2 = (N >= 2*step && N - 2*step >= x1)? x1 - 1 + step: N;
862 :
863 29130 : v = vecfactorsquarefreeu_coprime(x1, x2, P);
864 29132 : lv = lg(v);
865 59249664 : for (j = 1; j < lv; j++)
866 59221450 : if (gel(v,j))
867 : {
868 17996569 : ulong d = x1 - 1 + j; /* squarefree, coprime to 6 */
869 17996569 : GEN t = smallfact(d, gel(v,j), sq, V), u;
870 17956294 : GEN tB = NULL, uB = NULL; /* = f(d) d^s */
871 : ulong a, b, c, e, q;
872 17956294 : if (!t || gequal0(t)) continue;
873 3743282 : if (VB) tB = vecinv(gmulsg(d, conj_i(t)));
874 : /* warning: gives 1/conj(f(d)) d^(-1-conj(s)), equal to
875 : f(d) d^(-1-conj(s)) only if |f(d)|=1. */
876 : /* S += f(d) d^s * Z[q] */
877 3823390 : q = N / d;
878 3823390 : if (q == 1)
879 : {
880 1458993 : S = gadd(S, t); if (VB) SB = gadd(SB, tB);
881 1461818 : continue;
882 : }
883 2364397 : if (q <= sq) { u = gel(Z, q); if (VB) uB = gel(ZB, q); }
884 : else
885 : {
886 48421 : a = usqrt(q); b = usqrt(q / 2); c = usqrt(q / 3); e = usqrt(q / 6);
887 46977 : u = gadd(gadd(gel(Q,a), gel(Q2,b)), gadd(gel(Q3,c), gel(Q6,e)));
888 46977 : if (VB)
889 0 : uB = gadd(gadd(gel(QB,a), gel(Q2B,b)), gadd(gel(Q3B,c), gel(Q6B,e)));
890 : }
891 2362953 : S = gadd(S, vecmul(t, u)); if (VB) SB = gadd(SB, vecmul(tB, uB));
892 : }
893 28214 : return gerepilecopy(av, v2pack(S, SB));
894 : }
895 :
896 : GEN
897 38130 : parsumprimeWfunboth_worker(GEN gk, GEN s, GEN W, GEN WB, GEN f, GEN Nsqprec)
898 : {
899 : pari_sp av;
900 38130 : GEN S, SB = NULL, logp, tmp;
901 : forprime_t T;
902 38130 : long k = itou(gk), N = Nsqprec[1], sq = Nsqprec[2], precp;
903 38128 : long STEP = Nsqprec[3], prec0 = Nsqprec[4], prec1 = Nsqprec[5], p;
904 38128 : long gp[] = {evaltyp(t_INT)|_evallg(3), evalsigne(1)|evallgefint(3),0};
905 38128 : long needlog = get_needlog(s), nv;
906 :
907 38122 : if (isintzero(WB)) WB = NULL;
908 38112 : if (isintzero(f)) f = NULL;
909 38092 : tmp = mycallvec((void*)f, 1, prec1);
910 38091 : nv = typ(tmp) == t_VEC ? lg(tmp) - 1 : 0;
911 38091 : precp = 0; logp = NULL;
912 38091 : if (nv)
913 : {
914 15973 : S = const_vec(nv, real_0(prec1));
915 15964 : if (WB) SB = const_vec(nv, real_0(prec1));
916 : }
917 : else
918 : {
919 22118 : S = real_0(prec1);
920 22103 : if (WB) SB = real_0(prec1);
921 : }
922 38067 : u_forprime_init(&T, k * STEP + sq + 1, minss(N, (k + 1) * STEP + sq));
923 38075 : av = avma;
924 5201208 : while ((p = u_forprime_next(&T)))
925 : {
926 5165317 : GEN u = gen_0, ks = mycallvec((void*)f, p, prec1);
927 5165063 : long zks = !gequal0(ks);
928 5164990 : gp[2] = p;
929 5164990 : if (needlog)
930 : {
931 5164990 : if (!logp)
932 30268 : logp = logr_abs(utor(p, prec1));
933 : else
934 : { /* log p = log(precp) + 2 atanh((p - precp) / (p + precp)) */
935 5134722 : ulong a = p >> 1, b = precp >> 1; /* p = 2a + 1, precp = 2b + 1 */
936 5134722 : GEN z = atanhuu(a - b, a + b + 1, prec1); /* avoid overflow */
937 5132114 : shiftr_inplace(z, 1); logp = addrr(logp, z);
938 : }
939 5162782 : if (zks)
940 5162802 : u = needlog == 1? powcx(gp, logp, s, prec0) : mpexp(gmul(s, logp));
941 : }
942 0 : else { if (zks) u = gpow(gp, s, prec0); }
943 5165076 : if (zks)
944 : {
945 5165100 : S = gadd(S, vecmul(gel(W, N / p), gmul(ks, u)));
946 5163157 : if (WB)
947 0 : SB = gadd(SB, gdiv(vecmul(ks, gel(WB, N / p)), gmulsg(p, conj_i(u))));
948 : }
949 5163133 : precp = p;
950 5163133 : if ((p & 0x1ff) == 1)
951 : {
952 19285 : if (!logp) gerepileall(av, SB? 2: 1, &S, &SB);
953 19285 : else gerepileall(av, SB? 3: 2, &S, &logp, &SB);
954 : }
955 : }
956 37789 : return gcopy(v2pack(S, SB));
957 : }
958 :
959 : static GEN
960 1085 : pardirpowerssumfun_i(GEN f, ulong N, GEN s, long both, long prec)
961 : {
962 1085 : GEN P, V, W, Q, VB = NULL, WB = NULL, QB = NULL, c2, c2B = NULL;
963 1085 : GEN Q2, Q3, Q6, Q2B = NULL, Q3B = NULL, Q6B = NULL;
964 1085 : GEN S1, RES, Z, ZB = NULL, logp;
965 1085 : long gp[] = {evaltyp(t_INT)|_evallg(3), evalsigne(1)|evallgefint(3),0};
966 : ulong a, b, c, e, q, n, sq, fl;
967 1085 : long prec0, prec1, needlog, nv = 0;
968 1085 : GEN unvec = gen_1, zervec = gen_0, re0, re1, tmp2 = NULL;
969 :
970 1085 : if ((f && N < 49) || (!f && N < 10000UL))
971 112 : return smalldirpowerssum(N, s, (void*)f, mycallvec, both, prec);
972 973 : tmp2 = mycallvec((void*)f, 2, prec + EXTRAPRECWORD);
973 973 : if (is_vec_t(typ(tmp2)))
974 : {
975 385 : nv = lg(tmp2) - 1;
976 385 : if (!nv)
977 0 : return both? mkvec2(cgetg(1,t_VEC), cgetg(1,t_VEC)): cgetg(1,t_VEC);
978 385 : unvec = const_vec(nv, gen_1);
979 385 : zervec = const_vec(nv, gen_0);
980 : }
981 973 : fl = both && !gequal(real_i(s), gneg(ghalf));
982 973 : sq = usqrt(N);
983 973 : V = cgetg(sq+1, t_VEC); W = cgetg(sq+1, t_VEC); Q = cgetg(sq+1, t_VEC);
984 973 : if (fl)
985 0 : { VB = cgetg(sq+1, t_VEC); WB = cgetg(sq+1, t_VEC); QB = cgetg(sq+1, t_VEC); }
986 973 : prec1 = prec0 = prec + EXTRAPRECWORD;
987 973 : s = gprec_w(s, prec0);
988 973 : needlog = get_needlog(s);
989 973 : if (needlog == 1) prec1 = powcx_prec(log2((double)N), s, prec);
990 973 : gel(V,1) = gel(W,1) = gel(Q,1) = unvec;
991 973 : if (VB) { gel(VB,1) = gel(WB,1) = gel(QB,1) = unvec; }
992 973 : c2 = gpow(gen_2, s, prec0); if (VB) c2B = ginv(gmul2n(conj_i(c2), 1));
993 973 : re0 = real_0(prec0); re1 = real_1(prec0);
994 973 : if (f) { c2 = gmul(c2, tmp2); if (VB) c2B = gmul(c2B, tmp2); }
995 973 : gel(V,2) = c2; /* f(2) 2^s */
996 973 : gel(W,2) = gmul(re1, gadd(c2, unvec));
997 973 : gel(Q,2) = gmul(re1, gadd(vecsqr(c2), unvec));
998 973 : if (VB)
999 : {
1000 0 : gel(VB,2) = c2B; gel(WB,2) = gmul(re1, gadd(c2B, unvec));
1001 0 : gel(QB,2) = gmul(re1, gadd(vecsqr(c2B), unvec));
1002 : }
1003 973 : logp = NULL;
1004 154623 : for (n = 3; n <= sq; n++)
1005 : {
1006 153650 : GEN u = zervec, uB = zervec, ks = mycallvec((void*)f, n, prec);
1007 153650 : long zks = !gequal0(ks);
1008 153650 : if (odd(n))
1009 : {
1010 77287 : gp[2] = n;
1011 77287 : if (needlog)
1012 : {
1013 77287 : if (!logp)
1014 973 : logp = logr_abs(utor(n, prec1));
1015 : else
1016 : { /* log n = log(n-2) + 2 atanh(1 / (n - 1)) */
1017 76314 : GEN z = atanhuu(1, n - 1, prec1);
1018 76314 : shiftr_inplace(z, 1); logp = addrr(logp, z);
1019 : }
1020 77287 : if (zks)
1021 75712 : u = needlog == 1? powcx(gp, logp, s, prec0) : mpexp(gmul(s, logp));
1022 : }
1023 0 : else { if (zks) u = gpow(gp, s, prec0); }
1024 77287 : if (zks)
1025 : {
1026 75712 : if (VB) uB = gmul(ks, ginv(gmulsg(n, conj_i(u))));
1027 75712 : u = gmul(ks, u); /* f(n) n^s */
1028 : }
1029 : }
1030 : else
1031 : {
1032 76363 : u = vecmul(c2, gel(V, n >> 1));
1033 76363 : if (VB) uB = vecmul(c2B, gel(VB, n >> 1));
1034 : }
1035 153650 : gel(V,n) = u; /* f(n) n^s */
1036 153650 : gel(W,n) = gadd(gel(W,n-1), gel(V,n)); /* = sum_{i<=n} f(i)i^s */
1037 153650 : gel(Q,n) = gadd(gel(Q,n-1), vecsqr(gel(V,n)));/*= sum_{i<=n} f(i^2)i^2s*/
1038 153650 : if (VB)
1039 : {
1040 0 : gel(VB,n) = uB;
1041 0 : gel(WB,n) = gadd(gel(WB,n-1), gel(VB,n));
1042 0 : gel(QB,n) = gadd(gel(QB,n-1), vecsqr(gel(VB,n)));
1043 : }
1044 : }
1045 973 : Q2 = gmulvecsqlv(Q, gel(V,2));
1046 973 : Q3 = gmulvecsqlv(Q, gel(V,3));
1047 973 : Q6 = gmulvecsqlv(Q, gel(V,6));
1048 973 : if (VB)
1049 : {
1050 0 : Q2B = gmulvecsqlv(QB, gel(VB,2));
1051 0 : Q3B = gmulvecsqlv(QB, gel(VB,3));
1052 0 : Q6B = gmulvecsqlv(QB, gel(VB,6));
1053 : }
1054 973 : S1 = typ(zervec) == t_VEC? const_vec(nv, re0): re0;
1055 973 : RES = v2pack(S1, VB? S1: NULL);
1056 : {
1057 973 : long m = mt_nbthreads();
1058 973 : long STEP = maxss(N / (m * m), 1);
1059 973 : GEN VS = mkvecsmalln(5, N, sq, STEP, prec0, prec1);
1060 973 : GEN FUN = snm_closure(is_entry("_parsumprimeWfunboth_worker"),
1061 : mkvec5(s, W, WB? WB: gen_0, f? f: gen_0, VS));
1062 973 : RES = gadd(RES, parsum(gen_0, utoipos((N - 1) / STEP), FUN));
1063 : }
1064 973 : P = mkvecsmall2(2, 3);
1065 973 : Z = cgetg(sq+1, t_VEC);
1066 : /* a,b,c,e = sqrt(q), sqrt(q/2), sqrt(q/3), sqrt(q/6)
1067 : * Z[q] = Q[a] + 2^s Q[b] + 3^s Q[c] + 6^s Q[e], with Q[0] = 0 */
1068 973 : gel(Z, 1) = unvec;
1069 973 : gel(Z, 2) = gel(W, 2);
1070 973 : gel(Z, 3) = gel(W, 3);
1071 973 : gel(Z, 4) = gel(Z, 5) = gel(W, 4);
1072 973 : gel(Z, 6) = gel(Z, 7) = gadd(gel(W, 4), gel(V, 6));
1073 973 : if (VB)
1074 : {
1075 0 : ZB = cgetg(sq+1, t_VEC);
1076 0 : gel(ZB, 1) = unvec;
1077 0 : gel(ZB, 2) = gel(WB, 2);
1078 0 : gel(ZB, 3) = gel(WB, 3);
1079 0 : gel(ZB, 4) = gel(ZB, 5) = gel(WB, 4);
1080 0 : gel(ZB, 6) = gel(ZB, 7) = gadd(gel(WB, 4), gel(VB, 6));
1081 : }
1082 973 : a = 2; b = c = e = 1;
1083 149758 : for (q = 8; q <= sq; q++)
1084 : { /* Gray code: at most one of a,b,c,d differs (by 1) from previous value */
1085 148785 : GEN z = gel(Z, q - 1), zB = NULL;
1086 : ulong na, nb, nc, ne, na2, nb2, nc2, ne2;
1087 148785 : if (VB) zB = gel(ZB, q - 1);
1088 148785 : if ((na = usqrt(q)) != a)
1089 8911 : { a = na; na2 = na * na; z = gadd(z, gel(V, na2));
1090 8911 : if (VB) zB = gadd(zB, gel(VB, na2)); }
1091 139874 : else if ((nb = usqrt(q / 2)) != b)
1092 5999 : { b = nb; nb2 = 2 * nb * nb; z = gadd(z, gel(V, nb2));
1093 5999 : if (VB) zB = gadd(zB, gel(VB, nb2)); }
1094 133875 : else if ((nc = usqrt(q / 3)) != c)
1095 5257 : { c = nc; nc2 = 3 * nc * nc; z = gadd(z, gel(V, nc2));
1096 5257 : if (VB) zB = gadd(zB, gel(VB, nc2)); }
1097 128618 : else if ((ne = usqrt(q / 6)) != e)
1098 3010 : { e = ne; ne2 = 6 * ne * ne; z = gadd(z, gel(V, ne2));
1099 3010 : if (VB) zB = gadd(zB, gel(VB, ne2)); }
1100 148785 : gel(Z,q) = z; if (VB) gel(ZB,q) = zB;
1101 : }
1102 : {
1103 973 : GEN vQ = mkvec4(Q, Q2, Q3, Q6);
1104 973 : GEN vQB = VB? mkvec4(QB, Q2B, Q3B, Q6B): NULL;
1105 973 : GEN worker = snm_closure(is_entry("_parsqfboth_worker"),
1106 : mkvec5(v2pack(Z, ZB), v2pack(vQ, vQB), v2pack(V, VB),
1107 : P, mkvecsmall2(N, sq)));
1108 973 : RES = gadd(RES, parsum(gen_0, utoipos(maxss((N-1) / step - 1, 0)), worker));
1109 : }
1110 973 : return RES;
1111 : }
1112 :
1113 : GEN
1114 1085 : pardirpowerssumfun(GEN f, ulong N, GEN s, long both, long prec)
1115 : {
1116 1085 : pari_sp av = avma;
1117 1085 : return gerepilecopy(av, pardirpowerssumfun_i(f, N, s, both, prec));
1118 : }
1119 : GEN
1120 0 : pardirpowerssum0(GEN N, GEN s, GEN f, long both, long prec)
1121 : {
1122 0 : if (typ(N) != t_INT) pari_err_TYPE("pardirpowerssum", N);
1123 0 : return pardirpowerssumfun(f, itou(N), s, both, prec);
1124 : }
1125 : GEN
1126 0 : pardirpowerssum(ulong N, GEN s, long prec)
1127 0 : { return pardirpowerssumfun(NULL, N, s, 0, prec); }
|