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 : /** L-functions **/
18 : /** **/
19 : /********************************************************************/
20 :
21 : #include "pari.h"
22 : #include "paripriv.h"
23 :
24 : #define DEBUGLEVEL DEBUGLEVEL_lfun
25 :
26 : /*******************************************************************/
27 : /* Accessors */
28 : /*******************************************************************/
29 :
30 : static GEN
31 11749 : mysercoeff(GEN x, long n)
32 : {
33 11749 : long N = n - valser(x);
34 11749 : return (N < 0)? gen_0: gel(x, N+2);
35 : }
36 :
37 : long
38 12733 : ldata_get_type(GEN ldata) { return mael3(ldata, 1, 1, 1); }
39 :
40 : GEN
41 42021 : ldata_get_an(GEN ldata) { return gel(ldata, 1); }
42 :
43 : GEN
44 35413 : ldata_get_dual(GEN ldata) { return gel(ldata, 2); }
45 :
46 : long
47 2280 : ldata_isreal(GEN ldata) { return isintzero(gel(ldata, 2)); }
48 :
49 : GEN
50 238344 : ldata_get_gammavec(GEN ldata) { return gel(ldata, 3); }
51 :
52 : long
53 13809 : ldata_get_degree(GEN ldata) { return lg(gel(ldata, 3))-1; }
54 :
55 : GEN
56 99421 : ldata_get_k(GEN ldata)
57 : {
58 99421 : GEN w = gel(ldata,4);
59 99421 : if (typ(w) == t_VEC) w = gel(w,1);
60 99421 : return w;
61 : }
62 :
63 : /* a_n = O(n^{k1 + epsilon}) */
64 : GEN
65 98 : ldata_get_k1(GEN ldata)
66 : {
67 98 : GEN w = gel(ldata,4);
68 98 : if (typ(w) == t_VEC) return gel(w,2);
69 : /* by default, assume that k1 = k-1 and even (k-1)/2 for entire functions */
70 98 : w = gaddgs(w,-1);
71 98 : return ldata_get_residue(ldata)? w: gmul2n(w, -1);
72 : }
73 :
74 : /* a_n = O(n^{k1 + epsilon}) */
75 : static double
76 74464 : ldata_get_k1_dbl(GEN ldata)
77 : {
78 74464 : GEN w = gel(ldata,4);
79 : double k;
80 74464 : if (typ(w) == t_VEC) return gtodouble(gel(w,2));
81 : /* by default, assume that k1 = k-1 and even (k-1)/2 for entire functions */
82 73043 : k = gtodouble(w);
83 73043 : return ldata_get_residue(ldata)? k-1: (k-1)/2.;
84 : }
85 :
86 : GEN
87 170651 : ldata_get_conductor(GEN ldata) { return gel(ldata, 5); }
88 :
89 : GEN
90 68164 : ldata_get_rootno(GEN ldata) { return gel(ldata, 6); }
91 :
92 : GEN
93 128549 : ldata_get_residue(GEN ldata) { return lg(ldata) == 7 ? NULL: gel(ldata, 7); }
94 :
95 : long
96 93074 : linit_get_type(GEN linit) { return mael(linit, 1, 1); }
97 :
98 : GEN
99 139330 : linit_get_ldata(GEN linit) { return gel(linit, 2); }
100 :
101 : GEN
102 175415 : linit_get_tech(GEN linit) { return gel(linit, 3); }
103 :
104 : long
105 167565 : is_linit(GEN data)
106 : {
107 112170 : return lg(data) == 4 && typ(data) == t_VEC
108 279735 : && typ(gel(data, 1)) == t_VECSMALL;
109 : }
110 :
111 : GEN
112 21598 : lfun_get_step(GEN tech) { return gmael(tech, 2, 1);}
113 :
114 : GEN
115 21598 : lfun_get_pol(GEN tech) { return gmael(tech, 2, 2);}
116 :
117 : GEN
118 5066 : lfun_get_Residue(GEN tech) { return gmael(tech, 2, 3);}
119 :
120 : GEN
121 34329 : lfun_get_k2(GEN tech) { return gmael(tech, 3, 1);}
122 :
123 : GEN
124 13130 : lfun_get_w2(GEN tech) { return gmael(tech, 3, 2);}
125 :
126 : GEN
127 13130 : lfun_get_expot(GEN tech) { return gmael(tech, 3, 3);}
128 :
129 : GEN
130 5033 : lfun_get_factgammavec(GEN tech) { return gmael(tech, 3, 4); }
131 :
132 : /* Handle complex Vga whose sum is real */
133 : static GEN
134 80946 : sumVga(GEN Vga) { return real_i(vecsum(Vga)); }
135 : /* sum_i max (Im v[i],0) */
136 : static double
137 17768 : sumVgaimpos(GEN v)
138 : {
139 17768 : double d = 0.;
140 17768 : long i, l = lg(v);
141 47557 : for (i = 1; i < l; i++)
142 : {
143 29789 : GEN c = imag_i(gel(v,i));
144 29789 : if (gsigne(c) > 0) d += gtodouble(c);
145 : }
146 17768 : return d;
147 : }
148 :
149 : static long
150 30047 : vgaell(GEN Vga)
151 : {
152 30047 : if (lg(Vga) == 3)
153 19071 : { GEN c = gsub(gel(Vga,1), gel(Vga,2)); return gequal1(c) || gequalm1(c); }
154 10976 : return 0;
155 : }
156 : int
157 71774 : Vgaeasytheta(GEN Vga) { return lg(Vga)-1 == 1 || vgaell(Vga); }
158 : /* return b(n) := a(n) * n^c, when Vgaeasytheta(Vga) is set */
159 : static GEN
160 11004 : antwist(GEN an, GEN Vga, long prec)
161 : {
162 : long l, i;
163 11004 : GEN b, c = vecmin(Vga);
164 11004 : if (gequal0(c)) return an;
165 1764 : l = lg(an); b = cgetg(l, t_VEC);
166 1764 : if (gequal1(c))
167 : {
168 1036 : if (typ(an) == t_VECSMALL)
169 8631 : for (i = 1; i < l; i++) gel(b,i) = mulss(an[i], i);
170 : else
171 19551 : for (i = 1; i < l; i++) gel(b,i) = gmulgu(gel(an,i), i);
172 : }
173 : else
174 : {
175 728 : GEN v = vecpowug(l-1, c, prec);
176 728 : if (typ(an) == t_VECSMALL)
177 0 : for (i = 1; i < l; i++) gel(b,i) = gmulsg(an[i], gel(v,i));
178 : else
179 20965 : for (i = 1; i < l; i++) gel(b,i) = gmul(gel(an,i), gel(v,i));
180 : }
181 1764 : return b;
182 : }
183 :
184 : static GEN
185 6902 : theta_dual(GEN theta, GEN bn)
186 : {
187 6902 : if (typ(bn)==t_INT) return NULL;
188 : else
189 : {
190 77 : GEN thetad = shallowcopy(theta), ldata = linit_get_ldata(theta);
191 77 : GEN Vga = ldata_get_gammavec(ldata);
192 77 : GEN tech = shallowcopy(linit_get_tech(theta));
193 77 : GEN an = theta_get_an(tech);
194 77 : long prec = nbits2prec(theta_get_bitprec(tech));
195 77 : GEN vb = ldata_vecan(bn, lg(an)-1, prec);
196 77 : if (!theta_get_m(tech) && Vgaeasytheta(Vga)) vb = antwist(vb, Vga, prec);
197 77 : gel(tech,1) = vb;
198 77 : gel(thetad,3) = tech; return thetad;
199 : }
200 : }
201 :
202 : static GEN
203 37710 : domain_get_dom(GEN domain) { return gel(domain,1); }
204 : static long
205 16560 : domain_get_der(GEN domain) { return mael2(domain, 2, 1); }
206 : static long
207 22524 : domain_get_bitprec(GEN domain) { return mael2(domain, 2, 2); }
208 : GEN
209 38235 : lfun_get_domain(GEN tech) { return gel(tech,1); }
210 : long
211 77 : lfun_get_bitprec(GEN tech){ return domain_get_bitprec(lfun_get_domain(tech)); }
212 : GEN
213 0 : lfun_get_dom(GEN tech) { return domain_get_dom(lfun_get_domain(tech)); }
214 :
215 : GEN
216 2035 : lfunprod_get_fact(GEN tech) { return gel(tech, 2); }
217 :
218 : GEN
219 47239 : theta_get_an(GEN tdata) { return gel(tdata, 1);}
220 : GEN
221 7133 : theta_get_K(GEN tdata) { return gel(tdata, 2);}
222 : GEN
223 5935 : theta_get_R(GEN tdata) { return gel(tdata, 3);}
224 : long
225 60608 : theta_get_bitprec(GEN tdata) { return itos(gel(tdata, 4));}
226 : long
227 91552 : theta_get_m(GEN tdata) { return itos(gel(tdata, 5));}
228 : GEN
229 48611 : theta_get_tdom(GEN tdata) { return gel(tdata, 6);}
230 : GEN
231 52601 : theta_get_isqrtN(GEN tdata) { return gel(tdata, 7);}
232 :
233 : /*******************************************************************/
234 : /* Helper functions related to Gamma products */
235 : /*******************************************************************/
236 : /* x != 0 */
237 : static int
238 5747 : serisscalar(GEN x)
239 : {
240 : long i;
241 5747 : if (valser(x)) return 0;
242 7448 : for (i = lg(x)-1; i > 3; i--) if (!gequal0(gel(x,i))) return 0;
243 5509 : return 1;
244 : }
245 :
246 : /* return -itos(s) >= 0 if scalar s is (approximately) equal to a nonpositive
247 : * integer, and -1 otherwise */
248 : static long
249 14728 : isnegint(GEN s)
250 : {
251 14728 : GEN r = ground(real_i(s));
252 14728 : if (signe(r) <= 0 && gequal(s, r)) return -itos(r);
253 14602 : return -1;
254 : }
255 : /* if s = a + O(x^n), a <= 0 integer, replace by a + b*x^n + O(x^(n+1)) */
256 : static GEN
257 5768 : serextendifnegint(GEN s, GEN b, long *ext)
258 : {
259 5768 : if (!signe(s) || (serisscalar(s) && isnegint(gel(s,2)) >= 0))
260 : {
261 112 : long l = lg(s);
262 112 : GEN t = cgetg(l+1, t_SER);
263 301 : gel(t, l) = b; while (--l > 1) gel(t,l) = gel(s,l);
264 112 : if (gequal0(gel(t,2))) gel(t,2) = gen_0;
265 112 : t[1] = s[1]; s = normalizeser(t); *ext = 1;
266 : }
267 5768 : return s;
268 : }
269 :
270 : /* r/x + O(1), r != 0 */
271 : static GEN
272 3878 : serpole(GEN r)
273 : {
274 3878 : GEN s = cgetg(3, t_SER);
275 3878 : s[1] = evalsigne(1)|evalvalser(-1)|evalvarn(0);
276 3878 : gel(s,2) = r; return s;
277 : }
278 : /* a0 + a1 x + O(x^e), e >= 0 */
279 : static GEN
280 6706 : deg1ser_shallow(GEN a1, GEN a0, long v, long e)
281 6706 : { return RgX_to_ser(deg1pol_shallow(a1, a0, v), e+2); }
282 :
283 : /* pi^(-s/2) Gamma(s/2) */
284 : static GEN
285 8624 : gamma_R(GEN s, long *ext, long prec)
286 : {
287 8624 : GEN s2 = gmul2n(s, -1);
288 : long ms;
289 :
290 8624 : if (typ(s) == t_SER)
291 4018 : s2 = serextendifnegint(s2, ghalf, ext);
292 4606 : else if ((ms = isnegint(s2)) >= 0)
293 : {
294 35 : GEN r = gmul(powPis(stoi(ms),prec), gdivsg(odd(ms)? -2: 2, mpfact(ms)));
295 35 : return serpole(r);
296 : }
297 8589 : return gdiv(ggamma(s2,prec), powPis(s2,prec));
298 : }
299 : /* gamma_R(s)gamma_R(s+1) = 2 (2pi)^(-s) Gamma(s) */
300 : static GEN
301 6363 : gamma_C(GEN s, long *ext, long prec)
302 : {
303 : long ms;
304 6363 : if (typ(s) == t_SER)
305 1750 : s = serextendifnegint(s, gen_1, ext);
306 4613 : else if ((ms = isnegint(s)) >= 0)
307 : {
308 0 : GEN r = gmul(pow2Pis(stoi(ms),prec), gdivsg(odd(ms)? -2: 2, mpfact(ms)));
309 0 : return serpole(r);
310 : }
311 6363 : return gmul2n(gdiv(ggamma(s,prec), pow2Pis(s,prec)), 1);
312 : }
313 :
314 : static GEN
315 1708 : gammafrac(GEN r, long d)
316 : {
317 1708 : long i, l = labs(d) + 1, j = (d > 0)? 0: 2*d;
318 1708 : GEN T, v = cgetg(l, t_COL);
319 4466 : for (i = 1; i < l; i++, j += 2)
320 2758 : gel(v,i) = deg1pol_shallow(gen_1, gaddgs(r, j), 0);
321 1708 : T = RgV_prod(v); return d > 0? T: mkrfrac(gen_1, T);
322 : }
323 :
324 : /*
325 : GR(s)=Pi^-(s/2)*gamma(s/2);
326 : GC(s)=2*(2*Pi)^-s*gamma(s)
327 : gdirect(F,s)=prod(i=1,#F,GR(s+F[i]))
328 : gfact(F,s)=
329 : { my([R,A,B]=gammafactor(F), [a,e]=A, [b,f]=B, p=poldegree(R));
330 : subst(R,x,s) * (2*Pi)^-p * prod(i=1,#a,GR(s+a[i])^e[i])
331 : * prod(i=1,#b,GC(s+b[i])^f[i]); }
332 : */
333 : static GEN
334 15225 : gammafactor(GEN Vga)
335 : {
336 15225 : long i, r, c, l = lg(Vga);
337 15225 : GEN v, P, a, b, e, f, E, F = cgetg(l, t_VEC), R = gen_1;
338 43001 : for (i = 1; i < l; ++i)
339 : {
340 27776 : GEN a = gel(Vga,i), r = gmul2n(real_i(a), -1);
341 27776 : long q = itos(gfloor(r)); /* [Re a/2] */
342 27776 : r = gmul2n(gsubgs(r, q), 1);
343 27776 : gel(F,i) = gequal0(imag_i(a)) ? r : mkcomplex(r, imag_i(a)); /* 2{Re a/2} + I*(Im a) */
344 27776 : if (q) R = gmul(R, gammafrac(gel(F,i), q));
345 : }
346 15225 : F = vec_reduce(F, &E); l = lg(E);
347 15225 : v = cgetg(l, t_VEC);
348 37772 : for (i = 1; i < l; i++)
349 22547 : gel(v,i) = mkvec2(gsub(gel(F,i),gfloor(real_i(gel(F,i)))), stoi(E[i]));
350 15225 : gen_sort_inplace(v, (void*)cmp_universal, cmp_nodata, &P);
351 15225 : a = cgetg(l, t_VEC); e = cgetg(l, t_VECSMALL);
352 15225 : b = cgetg(l, t_VEC); f = cgetg(l, t_VECSMALL);
353 31395 : for (i = r = c = 1; i < l;)
354 16170 : if (i==l-1 || cmp_universal(gel(v,i), gel(v,i+1)))
355 9793 : { gel(a, r) = gel(F, P[i]); e[r++] = E[P[i]]; i++; }
356 : else
357 6377 : { gel(b, c) = gel(F, P[i]); f[c++] = E[P[i]]; i+=2; }
358 15225 : setlg(a, r); setlg(e, r);
359 15225 : setlg(b, c); setlg(f, c); return mkvec3(R, mkvec2(a,e), mkvec2(b,f));
360 : }
361 :
362 : static GEN
363 3640 : polgammaeval(GEN F, GEN s)
364 : {
365 3640 : GEN r = poleval(F, s);
366 3640 : if (typ(s) != t_SER && gequal0(r))
367 : { /* here typ(F) = t_POL */
368 : long e;
369 7 : for (e = 1;; e++)
370 : {
371 7 : F = RgX_deriv(F); r = poleval(F,s);
372 7 : if (!gequal0(r)) break;
373 : }
374 7 : if (e > 1) r = gdiv(r, mpfact(e));
375 7 : r = serpole(r); setvalser(r, e);
376 : }
377 3640 : return r;
378 : }
379 : static long
380 1799 : rfrac_degree(GEN R)
381 : {
382 1799 : GEN a = gel(R,1), b = gel(R,2);
383 1799 : return ((typ(a) == t_POL)? degpol(a): 0) - degpol(b);
384 : }
385 : static GEN
386 14000 : fracgammaeval(GEN F, GEN s, long prec)
387 : {
388 14000 : GEN R = gel(F,1);
389 : long d;
390 14000 : switch(typ(R))
391 : {
392 42 : case t_POL:
393 42 : d = degpol(R);
394 42 : R = polgammaeval(R, s); break;
395 1799 : case t_RFRAC:
396 1799 : d = rfrac_degree(R);
397 1799 : R = gdiv(polgammaeval(gel(R,1), s), polgammaeval(gel(R,2), s)); break;
398 12159 : default: return R;
399 : }
400 1841 : return gmul(R, powrs(Pi2n(1,prec), -d));
401 : }
402 :
403 : static GEN
404 14000 : gammafactproduct(GEN F, GEN s, long *ext, long prec)
405 : {
406 14000 : pari_sp av = avma;
407 14000 : GEN R = gel(F,2), Rw = gel(R,1), Re = gel(R,2);
408 14000 : GEN C = gel(F,3), Cw = gel(C,1), Ce = gel(C,2), z = fracgammaeval(F,s,prec);
409 14000 : long i, lR = lg(Rw), lC = lg(Cw);
410 14000 : *ext = 0;
411 22624 : for (i = 1; i < lR; i++)
412 8624 : z = gmul(z, gpowgs(gamma_R(gadd(s,gel(Rw, i)), ext, prec), Re[i]));
413 20363 : for (i = 1; i < lC; i++)
414 6363 : z = gmul(z, gpowgs(gamma_C(gadd(s,gel(Cw, i)), ext, prec), Ce[i]));
415 14000 : return gerepileupto(av, z);
416 : }
417 :
418 : static int
419 4739 : gammaordinary(GEN Vga, GEN s)
420 : {
421 4739 : long i, d = lg(Vga)-1;
422 12985 : for (i = 1; i <= d; i++)
423 : {
424 8337 : GEN z = gadd(s, gel(Vga,i));
425 : long e;
426 8337 : if (gexpo(imag_i(z)) < -10)
427 : {
428 8337 : z = real_i(z);
429 8337 : if (gsigne(z) <= 0) { (void)grndtoi(z, &e); if (e < -10) return 0; }
430 : }
431 : }
432 4648 : return 1;
433 : }
434 :
435 : /* Exponent A of t in asymptotic expansion; K(t) ~ C t^A exp(-pi d t^(2/d)).
436 : * suma = vecsum(Vga)*/
437 : static double
438 74457 : gammavec_expo(long d, double suma) { return (1 - d + suma) / d; }
439 :
440 : /*******************************************************************/
441 : /* First part: computations only involving Theta(t) */
442 : /*******************************************************************/
443 :
444 : static void
445 118747 : get_cone(GEN t, double *r, double *a)
446 : {
447 118747 : const long prec = LOWDEFAULTPREC;
448 118747 : if (typ(t) == t_COMPLEX)
449 : {
450 21098 : t = gprec_w(t, prec);
451 21098 : *r = gtodouble(gabs(t, prec));
452 21098 : *a = fabs(gtodouble(garg(t, prec)));
453 : }
454 : else
455 : {
456 97649 : *r = fabs(gtodouble(t));
457 97649 : *a = 0.;
458 : }
459 118747 : if (!*r && !*a) pari_err_DOMAIN("lfunthetainit","t","=",gen_0,t);
460 118740 : }
461 : /* slightly larger cone than necessary, to avoid round-off problems */
462 : static void
463 70136 : get_cone_fuzz(GEN t, double *r, double *a)
464 70136 : { get_cone(t, r, a); *r -= 1e-10; if (*a) *a += 1e-10; }
465 :
466 : /* Initialization m-th Theta derivative. tdom is either
467 : * - [rho,alpha]: assume |t| >= rho and |arg(t)| <= alpha
468 : * - a positive real scalar: assume t real, t >= tdom;
469 : * - a complex number t: compute at t;
470 : * N is the conductor (either the true one from ldata or a guess from
471 : * lfunconductor) */
472 : long
473 56696 : lfunthetacost(GEN ldata, GEN tdom, long m, long bitprec)
474 : {
475 56696 : pari_sp av = avma;
476 56696 : GEN Vga = ldata_get_gammavec(ldata);
477 56696 : long d = lg(Vga)-1;
478 56696 : double k1 = maxdd(ldata_get_k1_dbl(ldata), 0.);
479 56696 : double c = d/2., a, A, B, logC, al, rho, T;
480 56696 : double N = gtodouble(ldata_get_conductor(ldata));
481 :
482 56696 : if (!N) pari_err_TYPE("lfunthetaneed [missing conductor]", ldata);
483 56696 : if (typ(tdom) == t_VEC && lg(tdom) == 3)
484 : {
485 7 : rho= gtodouble(gel(tdom,1));
486 7 : al = gtodouble(gel(tdom,2));
487 : }
488 : else
489 56689 : get_cone_fuzz(tdom, &rho, &al);
490 56689 : A = gammavec_expo(d, gtodouble(sumVga(Vga))); set_avma(av);
491 56689 : a = (A+k1+1) + (m-1)/c;
492 56689 : if (fabs(a) < 1e-10) a = 0.;
493 56689 : logC = c*M_LN2 - log(c)/2;
494 : /* +1: fudge factor */
495 56689 : B = M_LN2*bitprec+logC+m*log(2*M_PI) + 1 + (k1+1)*log(N)/2 - (k1+m+1)*log(rho);
496 56689 : if (al)
497 : { /* t = rho e^(i*al), T^(1/c) = Re(t^(1/c)) > 0, T = rho cos^c(al/c) */
498 10556 : double z = cos(al/c);
499 10556 : if (z <= 0)
500 7 : pari_err_DOMAIN("lfunthetaneed", "arg t", ">", dbltor(c*M_PI/2), tdom);
501 10549 : T = (d == 2 && typ(tdom) != t_VEC)? gtodouble(real_i(tdom)): rho*pow(z,c);
502 10549 : B -= log(z) * (c * (k1+A+1) + m);
503 : }
504 : else
505 46133 : T = rho;
506 56682 : if (B <= 0) return 0;
507 56682 : A = floor(0.9 + dblcoro526(a,c,B) / T * sqrt(N));
508 56682 : if (dblexpo(A) >= BITS_IN_LONG-1) pari_err_OVERFLOW("lfunthetacost");
509 56675 : return (long)A;
510 : }
511 : long
512 21 : lfunthetacost0(GEN L, GEN tdom, long m, long bitprec)
513 : {
514 : long n;
515 21 : if (is_linit(L) && linit_get_type(L)==t_LDESC_THETA)
516 7 : {
517 7 : GEN tech = linit_get_tech(L);
518 7 : n = lg(theta_get_an(tech))-1;
519 : }
520 : else
521 : {
522 14 : pari_sp av = avma;
523 14 : GEN ldata = lfunmisc_to_ldata_shallow(L);
524 14 : n = lfunthetacost(ldata, tdom? tdom: gen_1, m, bitprec);
525 7 : set_avma(av);
526 : }
527 14 : return n;
528 : }
529 :
530 : static long
531 5544 : fracgammadegree(GEN FVga)
532 5544 : { GEN F = gel(FVga,1); return (typ(F)==t_RFRAC)? degpol(gel(F,2)): 0; }
533 :
534 : /* Poles of a L-function can be represented in the following ways:
535 : * 1) Nothing (ldata has only 6 components, ldata_get_residue = NULL).
536 : * 2) a complex number (single pole at s = k with given residue, unknown if 0).
537 : * 3) A vector (possibly empty) of 2-component vectors [a, ra], where a is the
538 : * pole, ra a t_SER: its Taylor expansion at a. A t_VEC encodes the polar
539 : * part of L, a t_COL, the polar part of Lambda */
540 :
541 : /* 'a' a complex number (pole), 'r' the polar part of L at 'a';
542 : * return 'R' the polar part of Lambda at 'a' */
543 : static GEN
544 3920 : rtoR(GEN a, GEN r, GEN FVga, GEN N, long prec)
545 : {
546 3920 : long v = lg(r)-2, d = fracgammadegree(FVga), ext;
547 3920 : GEN Na, as = deg1ser_shallow(gen_1, a, varn(r), v);
548 3920 : Na = gpow(N, gdivgu(as, 2), prec);
549 : /* make up for a possible loss of accuracy */
550 3920 : if (d) as = deg1ser_shallow(gen_1, a, varn(r), v + d);
551 3920 : return gmul(gmul(r, Na), gammafactproduct(FVga, as, &ext, prec));
552 : }
553 :
554 : /* assume r in normalized form: t_VEC of pairs [be,re] */
555 : GEN
556 3689 : lfunrtopoles(GEN r)
557 : {
558 3689 : long j, l = lg(r);
559 3689 : GEN v = cgetg(l, t_VEC);
560 7609 : for (j = 1; j < l; j++)
561 : {
562 3920 : GEN rj = gel(r,j), a = gel(rj,1);
563 3920 : gel(v,j) = a;
564 : }
565 3689 : gen_sort_inplace(v, (void*)&cmp_universal, cmp_nodata, NULL);
566 3689 : return v;
567 : }
568 :
569 : /* r / x + O(1) */
570 : static GEN
571 3969 : simple_pole(GEN r)
572 3969 : { return isintzero(r)? gen_0: serpole(r); }
573 : static GEN
574 4809 : normalize_simple_pole(GEN r, GEN k)
575 : {
576 4809 : long tx = typ(r);
577 4809 : if (is_vec_t(tx)) return r;
578 3969 : if (!is_scalar_t(tx)) pari_err_TYPE("lfunrootres [poles]", r);
579 3969 : return mkvec(mkvec2(k, simple_pole(r)));
580 : }
581 : /* normalize the description of a polar part */
582 : static GEN
583 4536 : normalizepoles(GEN r, GEN k)
584 : {
585 : long iv, j, l;
586 : GEN v;
587 4536 : if (!is_vec_t(typ(r))) return normalize_simple_pole(r, k);
588 2163 : v = cgetg_copy(r, &l);
589 5425 : for (j = iv = 1; j < l; j++)
590 : {
591 3262 : GEN rj = gel(r,j), a = gel(rj,1), ra = gel(rj,2);
592 3262 : if (!is_scalar_t(typ(a)) || typ(ra) != t_SER)
593 0 : pari_err_TYPE("lfunrootres [poles]",r);
594 3262 : if (valser(ra) >= 0) continue;
595 3248 : gel(v,iv++) = rj;
596 : }
597 2163 : setlg(v, iv); return v;
598 : }
599 : static int
600 7371 : residues_known(GEN r)
601 : {
602 7371 : long i, l = lg(r);
603 7371 : if (isintzero(r)) return 0;
604 7070 : if (!is_vec_t(typ(r))) return 1;
605 9968 : for (i = 1; i < l; i++)
606 : {
607 6167 : GEN ri = gel(r,i);
608 6167 : if (!is_vec_t(typ(ri)) || lg(ri)!=3)
609 0 : pari_err_TYPE("lfunrootres [poles]",r);
610 6167 : if (isintzero(gel(ri, 2))) return 0;
611 : }
612 3801 : return 1;
613 : }
614 :
615 : /* Compute R's from r's (r = Taylor devts of L(s), R of Lambda(s)).
616 : * 'r/eno' passed to override the one from ldata */
617 : static GEN
618 15568 : lfunrtoR_i(GEN ldata, GEN r, GEN eno, long prec)
619 : {
620 15568 : GEN Vga = ldata_get_gammavec(ldata), N = ldata_get_conductor(ldata);
621 : GEN R, vr, FVga;
622 15568 : pari_sp av = avma;
623 : long lr, j, jR;
624 15568 : GEN k = ldata_get_k(ldata);
625 :
626 15568 : if (!r || isintzero(eno) || !residues_known(r))
627 11032 : return gen_0;
628 4536 : r = normalizepoles(r, k);
629 4536 : if (typ(r) == t_COL) return gerepilecopy(av, r);
630 3689 : if (typ(ldata_get_dual(ldata)) != t_INT)
631 0 : pari_err(e_MISC,"please give the Taylor development of Lambda");
632 3689 : vr = lfunrtopoles(r); lr = lg(vr);
633 3689 : FVga = gammafactor(Vga);
634 3689 : R = cgetg(2*lr, t_COL);
635 7609 : for (j = jR = 1; j < lr; j++)
636 : {
637 3920 : GEN rj = gel(r,j), a = gel(rj,1), ra = gel(rj,2);
638 3920 : GEN Ra = rtoR(a, ra, FVga, N, prec);
639 3920 : GEN b = gsub(k, conj_i(a));
640 3920 : if (lg(Ra)-2 < -valser(Ra))
641 0 : pari_err(e_MISC,
642 : "please give more terms in L function's Taylor development at %Ps", a);
643 3920 : gel(R,jR++) = mkvec2(a, Ra);
644 3920 : if (!tablesearch(vr, b, (int (*)(GEN,GEN))&cmp_universal))
645 : {
646 3731 : GEN mX = gneg(pol_x(varn(Ra)));
647 3731 : GEN Rb = gmul(eno, gsubst(conj_i(Ra), varn(Ra), mX));
648 3731 : gel(R,jR++) = mkvec2(b, Rb);
649 : }
650 : }
651 3689 : setlg(R, jR); return gerepilecopy(av, R);
652 : }
653 : static GEN
654 15106 : lfunrtoR_eno(GEN ldata, GEN eno, long prec)
655 15106 : { return lfunrtoR_i(ldata, ldata_get_residue(ldata), eno, prec); }
656 : static GEN
657 13454 : lfunrtoR(GEN ldata, long prec)
658 13454 : { return lfunrtoR_eno(ldata, ldata_get_rootno(ldata), prec); }
659 :
660 : static long
661 13454 : prec_fix(long prec)
662 : {
663 : #ifndef LONG_IS_64BIT
664 : /* make sure that default accuracy is the same on 32/64bit */
665 1922 : if (odd(prec)) prec += EXTRAPREC64;
666 : #endif
667 13454 : return prec;
668 : }
669 :
670 : /* thetainit using {an: n <= L}; if (m = 0 && easytheta), an2 is an * n^al */
671 : static GEN
672 13454 : lfunthetainit0(GEN ldata, GEN tdom, GEN an2, long m,
673 : long bitprec, long extrabit)
674 : {
675 13454 : long prec = nbits2prec(bitprec);
676 13454 : GEN tech, N = ldata_get_conductor(ldata);
677 13454 : GEN K = gammamellininvinit(ldata, m, bitprec + extrabit);
678 13454 : GEN R = lfunrtoR(ldata, prec);
679 13454 : if (!tdom) tdom = gen_1;
680 13454 : if (typ(tdom) != t_VEC)
681 : {
682 : double r, a;
683 13447 : get_cone_fuzz(tdom, &r, &a);
684 13447 : tdom = mkvec2(dbltor(r), a? dbltor(a): gen_0);
685 : }
686 13454 : prec += maxss(EXTRAPREC64, nbits2extraprec(extrabit));
687 13454 : tech = mkvecn(7, an2,K,R, stoi(bitprec), stoi(m), tdom,
688 : gsqrt(ginv(N), prec_fix(prec)));
689 13454 : return mkvec3(mkvecsmall(t_LDESC_THETA), ldata, tech);
690 : }
691 :
692 : /* tdom: 1) positive real number r, t real, t >= r; or
693 : * 2) [r,a], describing the cone |t| >= r, |arg(t)| <= a */
694 : static GEN
695 7812 : lfunthetainit_i(GEN data, GEN tdom, long m, long bit)
696 : {
697 7812 : GEN ldata = lfunmisc_to_ldata_shallow(data);
698 7812 : long b = 32, L = lfunthetacost(ldata, tdom, m, bit), prec = nbits2prec(bit);
699 7798 : GEN ldatan = ldata_newprec(ldata, prec);
700 7798 : GEN an = ldata_vecan(ldata_get_an(ldatan), L, prec);
701 7798 : GEN Vga = ldata_get_gammavec(ldatan);
702 7798 : if (m == 0 && Vgaeasytheta(Vga)) an = antwist(an, Vga, prec);
703 7798 : if (typ(an) != t_VECSMALL) b = maxss(b, gexpo(an));
704 7798 : return lfunthetainit0(ldatan, tdom, an, m, bit, b);
705 : }
706 :
707 : GEN
708 322 : lfunthetainit(GEN ldata, GEN tdom, long m, long bitprec)
709 : {
710 322 : pari_sp av = avma;
711 322 : GEN S = lfunthetainit_i(ldata, tdom? tdom: gen_1, m, bitprec);
712 322 : return gerepilecopy(av, S);
713 : }
714 :
715 : GEN
716 1274 : lfunan(GEN ldata, long L, long prec)
717 : {
718 1274 : pari_sp av = avma;
719 : GEN an ;
720 1274 : ldata = ldata_newprec(lfunmisc_to_ldata_shallow(ldata), prec);
721 1274 : an = gerepilecopy(av, ldata_vecan(ldata_get_an(ldata), L, prec));
722 1218 : if (typ(an) != t_VEC) an = vecsmall_to_vec_inplace(an);
723 1218 : return an;
724 : }
725 :
726 : static GEN
727 10773 : mulrealvec(GEN x, GEN y)
728 : {
729 10773 : if (is_vec_t(typ(x)) && is_vec_t(typ(y)))
730 84 : pari_APPLY_same(mulreal(gel(x,i),gel(y,i)))
731 : else
732 10745 : return mulreal(x,y);
733 : }
734 : static GEN
735 21199 : gmulvec(GEN x, GEN y)
736 : {
737 21199 : if (is_vec_t(typ(x)) && is_vec_t(typ(y)))
738 1974 : pari_APPLY_same(gmul(gel(x,i),gel(y,i)))
739 : else
740 20674 : return gmul(x,y);
741 : }
742 : static GEN
743 6881 : gdivvec(GEN x, GEN y)
744 : {
745 6881 : if (is_vec_t(typ(x)) && is_vec_t(typ(y)))
746 1498 : pari_APPLY_same(gdiv(gel(x,i),gel(y,i)))
747 : else
748 6454 : return gdiv(x,y);
749 : }
750 :
751 : static GEN
752 3486 : gsubvec(GEN x, GEN y)
753 : {
754 3486 : if (is_vec_t(typ(x)) && !is_vec_t(typ(y)))
755 0 : pari_APPLY_same(gsub(gel(x,i),y))
756 : else
757 3486 : return gsub(x,y);
758 : }
759 :
760 : /* return [1^(2/d), 2^(2/d),...,lim^(2/d)] */
761 : static GEN
762 7133 : mkvroots(long d, long lim, long prec)
763 : {
764 7133 : if (d <= 4)
765 : {
766 6783 : GEN v = cgetg(lim+1,t_VEC);
767 : long n;
768 6783 : switch(d)
769 : {
770 2072 : case 1:
771 38204 : for (n=1; n <= lim; n++) gel(v,n) = sqru(n);
772 2072 : return v;
773 1267 : case 2:
774 227213 : for (n=1; n <= lim; n++) gel(v,n) = utoipos(n);
775 1267 : return v;
776 1983 : case 4:
777 6079017 : for (n=1; n <= lim; n++) gel(v,n) = sqrtr(utor(n, prec));
778 1983 : return v;
779 : }
780 : }
781 1811 : return vecpowug(lim, gdivgu(gen_2,d), prec);
782 : }
783 :
784 : GEN
785 54169 : lfunthetacheckinit(GEN data, GEN t, long m, long bitprec)
786 : {
787 54169 : if (is_linit(data) && linit_get_type(data)==t_LDESC_THETA)
788 : {
789 48611 : GEN tdom, thetainit = linit_get_tech(data);
790 48611 : long bitprecnew = theta_get_bitprec(thetainit);
791 48611 : long m0 = theta_get_m(thetainit);
792 : double r, al, rt, alt;
793 48611 : if (m0 != m)
794 0 : pari_err_DOMAIN("lfuntheta","derivative order","!=", stoi(m),stoi(m0));
795 48611 : if (bitprec > bitprecnew) goto INIT;
796 48611 : get_cone(t, &rt, &alt);
797 48611 : tdom = theta_get_tdom(thetainit);
798 48611 : r = gtodouble(gel(tdom,1));
799 48611 : al= gtodouble(gel(tdom,2)); if (rt >= r && alt <= al) return data;
800 : }
801 5558 : INIT:
802 7399 : return lfunthetainit_i(data, t, m, bitprec);
803 : }
804 :
805 : static GEN
806 14711879 : get_an(GEN an, long n)
807 : {
808 14711879 : if (typ(an) == t_VECSMALL) { long a = an[n]; if (a) return stoi(a); }
809 14711879 : else { GEN a = gel(an,n); if (a && !gequal0(a)) return a; }
810 12686150 : return NULL;
811 : }
812 : /* x * an[n] */
813 : static GEN
814 11134889 : mul_an(GEN an, long n, GEN x)
815 : {
816 11134889 : if (typ(an) == t_VECSMALL) { long a = an[n]; if (a) return gmulsg(a,x); }
817 5748279 : else { GEN a = gel(an,n); if (a && !gequal0(a)) return gmul(a,x); }
818 2293982 : return NULL;
819 : }
820 : /* 2*t^a * x **/
821 : static GEN
822 191210 : mulT(GEN t, GEN a, GEN x, long prec)
823 : {
824 191210 : if (gequal0(a)) return gmul2n(x,1);
825 17898 : return gmul(x, gmul2n(gequal1(a)? t: gpow(t,a,prec), 1));
826 : }
827 :
828 : static GEN
829 28308410 : vecan_cmul(void *E, GEN P, long a, GEN x)
830 : {
831 : (void)E;
832 28308410 : if (typ(P) == t_VECSMALL)
833 23778184 : return (a==0 || !P[a])? NULL: gmulsg(P[a], x);
834 : else
835 4530226 : return (a==0 || !gel(P,a))? NULL: gmul(gel(P,a), x);
836 : }
837 : /* d=2, 2 sum_{n <= N} a(n) (n t)^al q^n, q = exp(-2pi t),
838 : * an2[n] = a(n) * n^al */
839 : static GEN
840 154042 : theta2_i(GEN an2, long N, GEN t, GEN al, long prec)
841 : {
842 154042 : GEN S, q, pi2 = Pi2n(1,prec);
843 154038 : const struct bb_algebra *alg = get_Rg_algebra();
844 154038 : setsigne(pi2,-1); q = gexp(gmul(pi2, t), prec);
845 : /* Brent-Kung in case the a_n are small integers */
846 154041 : S = gen_bkeval(an2, N, q, 1, NULL, alg, vecan_cmul);
847 154038 : return mulT(t, al, S, prec);
848 : }
849 : static GEN
850 148348 : theta2(GEN an2, long N, GEN t, GEN al, long prec)
851 : {
852 148348 : pari_sp av = avma;
853 148348 : return gerepileupto(av, theta2_i(an2, N, t, al, prec));
854 : }
855 :
856 : /* d=1, 2 sum_{n <= N} a_n (n t)^al q^(n^2), q = exp(-pi t^2),
857 : * an2[n] is a_n n^al */
858 : static GEN
859 37170 : theta1(GEN an2, long N, GEN t, GEN al, long prec)
860 : {
861 37170 : GEN q = gexp(gmul(negr(mppi(prec)), gsqr(t)), prec);
862 37170 : GEN vexp = gsqrpowers(q, N), S = gen_0;
863 37170 : pari_sp av = avma;
864 : long n;
865 6748474 : for (n = 1; n <= N; n++)
866 : {
867 6711304 : GEN c = mul_an(an2, n, gel(vexp,n));
868 6711304 : if (c)
869 : {
870 5644015 : S = gadd(S, c);
871 5644015 : if (gc_needed(av, 3)) S = gerepileupto(av, S);
872 : }
873 : }
874 37170 : return mulT(t, al, S, prec);
875 : }
876 :
877 : /* If m > 0, compute m-th derivative of theta(t) = theta0(t/sqrt(N))
878 : * with absolute error 2^-bitprec; theta(t)=\sum_{n\ge1}a(n)K(nt/N^(1/2)) */
879 : GEN
880 46952 : lfuntheta(GEN data, GEN t, long m, long bitprec)
881 : {
882 46952 : pari_sp ltop = avma;
883 : long limt, d;
884 : GEN isqN, vecan, Vga, ldata, theta, thetainit, S;
885 : long n, prec;
886 :
887 46952 : theta = lfunthetacheckinit(data, t, m, bitprec);
888 46945 : ldata = linit_get_ldata(theta);
889 46945 : thetainit = linit_get_tech(theta);
890 46945 : vecan = theta_get_an(thetainit);
891 46945 : isqN = theta_get_isqrtN(thetainit);
892 46945 : prec = maxss(realprec(isqN), nbits2prec(bitprec));
893 46945 : t = gprec_w(t, prec);
894 46945 : limt = lg(vecan)-1;
895 46945 : if (theta == data)
896 45111 : limt = minss(limt, lfunthetacost(ldata, t, m, bitprec));
897 46945 : if (!limt)
898 : {
899 14 : set_avma(ltop); S = real_0_bit(-bitprec);
900 14 : if (!is_real_t(typ(t)) || !ldata_isreal(ldata))
901 7 : S = gerepilecopy(ltop, mkcomplex(S,S));
902 14 : return S;
903 : }
904 46931 : t = gmul(t, isqN);
905 46931 : Vga = ldata_get_gammavec(ldata);
906 46931 : d = lg(Vga)-1;
907 46931 : if (m == 0 && Vgaeasytheta(Vga))
908 : {
909 42864 : if (theta_get_m(thetainit) > 0) vecan = antwist(vecan, Vga, prec);
910 42864 : if (d == 1) S = theta1(vecan, limt, t, gel(Vga,1), prec);
911 5694 : else S = theta2_i(vecan, limt, t, vecmin(Vga), prec);
912 : }
913 : else
914 : {
915 4067 : GEN K = theta_get_K(thetainit);
916 4067 : GEN vroots = mkvroots(d, limt, prec);
917 : pari_sp av;
918 4067 : t = gpow(t, gdivgu(gen_2,d), prec);
919 4067 : S = gen_0; av = avma;
920 14715946 : for (n = 1; n <= limt; ++n)
921 : {
922 14711879 : GEN nt, an = get_an(vecan, n);
923 14711879 : if (!an) continue;
924 2025729 : nt = gmul(gel(vroots,n), t);
925 2025729 : if (m) an = gmul(an, powuu(n, m));
926 2025729 : S = gadd(S, gmul(an, gammamellininvrt(K, nt, bitprec)));
927 2025729 : if ((n & 0x1ff) == 0) S = gerepileupto(av, S);
928 : }
929 4067 : if (m) S = gmul(S, gpowgs(isqN, m));
930 : }
931 46931 : return gerepileupto(ltop, S);
932 : }
933 :
934 : /*******************************************************************/
935 : /* Second part: Computation of L-Functions. */
936 : /*******************************************************************/
937 :
938 : struct lfunp {
939 : long precmax, Dmax, D, M, m0, nmax, d, vgaell;
940 : double k1, dc, dw, dh, MAXs, sub;
941 : GEN L, an, bn;
942 : };
943 :
944 : static void
945 17768 : lfunp_set(GEN ldata, long der, long bitprec, struct lfunp *S)
946 : {
947 17768 : const long derprec = (der > 1)? dbllog2(mpfact(der)): 0; /* log2(der!) */
948 : GEN Vga, N, L, k;
949 : long k1, d, m, M, flag, nmax;
950 : double a, A, E, hd, Ep, d2, suma, maxs, mins, sub, B0,B1;
951 : double logN2, logC, Lestimate, Mestimate;
952 :
953 17768 : Vga = ldata_get_gammavec(ldata);
954 17768 : S->d = d = lg(Vga)-1; d2 = d/2.;
955 :
956 17768 : suma = gtodouble(sumVga(Vga));
957 17768 : k = ldata_get_k(ldata);
958 17768 : N = ldata_get_conductor(ldata);
959 17768 : logN2 = log(gtodouble(N)) / 2;
960 17768 : maxs = S->dc + S->dw;
961 17768 : mins = S->dc - S->dw;
962 17768 : S->MAXs = maxdd(maxs, gtodouble(k)-mins);
963 :
964 : /* we compute Lambda^(der)(s) / der!; need to compensate for L^(der)(s)
965 : * ln |gamma(s)| ~ -(pi/4) \sum_i |Im(s + a_i)|; max with 1: fudge factor */
966 17768 : a = (M_PI/(4*M_LN2))*(d*S->dh + sumVgaimpos(Vga));
967 17768 : S->D = (long)ceil(bitprec + derprec + maxdd(a, 1));
968 17768 : E = M_LN2*S->D; /* D:= required absolute bitprec */
969 :
970 17768 : Ep = E + maxdd(M_PI * S->dh * d2, (d*S->MAXs + suma - 1) * log(E));
971 17768 : hd = d2*M_PI*M_PI / Ep;
972 17768 : S->m0 = (long)ceil(M_LN2/hd);
973 17768 : hd = M_LN2/S->m0;
974 :
975 17768 : logC = d2*M_LN2 - log(d2)/2;
976 17768 : k1 = maxdd(ldata_get_k1_dbl(ldata), 0.);
977 17768 : S->k1 = k1; /* assume |a_n| << n^k1 with small implied constant */
978 17768 : A = gammavec_expo(d, suma);
979 :
980 17768 : sub = 0.;
981 17768 : if (mins > 1)
982 : {
983 4739 : GEN sig = dbltor(mins);
984 4739 : sub += logN2*mins;
985 4739 : if (gammaordinary(Vga, sig))
986 : {
987 : long ext;
988 4648 : GEN gas = gammafactproduct(gammafactor(Vga), sig, &ext, LOWDEFAULTPREC);
989 4648 : if (typ(gas) != t_SER)
990 : {
991 4648 : double dg = dbllog2(gas);
992 4648 : if (dg > 0) sub += dg * M_LN2;
993 : }
994 : }
995 : }
996 17768 : S->sub = sub;
997 17768 : M = 1000;
998 17768 : L = cgetg(M+2, t_VECSMALL);
999 17768 : a = S->k1 + A;
1000 :
1001 17768 : B0 = 5 + E - S->sub + logC + S->k1*logN2; /* 5 extra bits */
1002 17768 : B1 = hd * (S->MAXs - S->k1);
1003 17768 : Lestimate = dblcoro526(a + S->MAXs - 2./d, d/2.,
1004 17768 : E - S->sub + logC - log(2*M_PI*hd) + S->MAXs*logN2);
1005 17768 : Mestimate = ((Lestimate > 0? log(Lestimate): 0) + logN2) / hd;
1006 17768 : nmax = 0;
1007 17768 : flag = 0;
1008 17768 : for (m = 0;; m++)
1009 1657510 : {
1010 1675278 : double x, H = logN2 - m*hd, B = B0 + m*B1;
1011 : long n;
1012 1675278 : x = dblcoro526(a, d/2., B);
1013 1675278 : n = floor(x*exp(H));
1014 1675278 : if (n > nmax) nmax = n;
1015 1675278 : if (m > M) { M *= 2; L = vecsmall_lengthen(L,M+2); }
1016 1675278 : L[m+1] = n;
1017 1675278 : if (n == 0) { if (++flag > 2 && m > Mestimate) break; } else flag = 0;
1018 : }
1019 18608 : m -= 2; while (m > 0 && !L[m]) m--;
1020 17768 : if (m == 0) { nmax = 1; L[1] = 1; m = 1; } /* can happen for tiny bitprec */
1021 17768 : setlg(L, m+1); S->M = m-1;
1022 17768 : S->L = L;
1023 17768 : S->nmax = nmax;
1024 :
1025 17768 : S->Dmax = S->D + (long)ceil((S->M * hd * S->MAXs - S->sub) / M_LN2);
1026 17768 : if (S->Dmax < S->D) S->Dmax = S->D;
1027 17768 : S->precmax = nbits2prec(S->Dmax);
1028 17768 : if (DEBUGLEVEL > 1)
1029 0 : err_printf("Dmax=%ld, D=%ld, M = %ld, nmax = %ld, m0 = %ld\n",
1030 : S->Dmax,S->D,S->M,S->nmax, S->m0);
1031 17768 : }
1032 :
1033 : static GEN
1034 5810 : lfuninit_pol(GEN v, GEN poqk, long prec)
1035 : {
1036 5810 : long m, M = lg(v) - 2;
1037 5810 : GEN pol = cgetg(M+3, t_POL);
1038 5810 : pol[1] = evalsigne(1) | evalvarn(0);
1039 5810 : gel(pol, 2) = gprec_w(gmul2n(gel(v,1), -1), prec);
1040 5810 : if (poqk)
1041 301158 : for (m = 2; m <= M+1; m++)
1042 295404 : gel(pol, m+1) = gprec_w(gmul(gel(poqk,m), gel(v,m)), prec);
1043 : else
1044 2324 : for (m = 2; m <= M+1; m++)
1045 2268 : gel(pol, m+1) = gprec_w(gel(v,m), prec);
1046 5810 : return RgX_renormalize_lg(pol, M+3);
1047 : }
1048 :
1049 : static void
1050 51804 : worker_init(long q, GEN *an, GEN *bn, GEN *AB, GEN *A, GEN *B)
1051 : {
1052 51804 : if (typ(*bn) == t_INT) *bn = NULL;
1053 51804 : if (*bn)
1054 : {
1055 1016 : *AB = cgetg(3, t_VEC);
1056 1016 : gel(*AB,1) = *A = cgetg(q+1, t_VEC);
1057 1016 : gel(*AB,2) = *B = cgetg(q+1, t_VEC);
1058 1016 : if (typ(an) == t_VEC) *an = RgV_kill0(*an);
1059 1016 : if (typ(bn) == t_VEC) *bn = RgV_kill0(*bn);
1060 : }
1061 : else
1062 : {
1063 50788 : *B = NULL;
1064 50788 : *AB = *A = cgetg(q+1, t_VEC);
1065 50787 : if (typ(*an) == t_VEC) *an = RgV_kill0(*an);
1066 : }
1067 51802 : }
1068 : GEN
1069 14508 : lfuninit_theta2_worker(long r, GEN L, GEN qk, GEN a, GEN di, GEN an, GEN bn)
1070 : {
1071 14508 : long q, m, prec = di[1], M = di[2], m0 = di[3], L0 = lg(an)-1;
1072 : GEN AB, A, B;
1073 14508 : worker_init((M - r) / m0 + 1, &an, &bn, &AB, &A, &B);
1074 155887 : for (q = 0, m = r; m <= M; m += m0, q++)
1075 : {
1076 141381 : GEN t = gel(qk, m+1);
1077 141381 : long N = minss(L[m+1],L0);
1078 141382 : gel(A, q+1) = theta2(an, N, t, a, prec); /* theta(exp(mh)) */
1079 141380 : if (bn) gel(B, q+1) = theta2(bn, N, t, a, prec);
1080 : }
1081 14506 : return AB;
1082 : }
1083 :
1084 : /* theta(exp(mh)) ~ sum_{n <= N} a(n) k[m,n] */
1085 : static GEN
1086 155090 : an_msum(GEN an, long N, GEN vKm)
1087 : {
1088 155090 : pari_sp av = avma;
1089 155090 : GEN s = gen_0;
1090 : long n;
1091 10047214 : for (n = 1; n <= N; n++)
1092 9892639 : if (gel(vKm,n))
1093 : {
1094 4423496 : GEN c = mul_an(an, n, gel(vKm,n));
1095 4423876 : if (c) s = gadd(s, c);
1096 : }
1097 154575 : return gerepileupto(av, s);
1098 : }
1099 :
1100 : GEN
1101 37294 : lfuninit_worker(long r, GEN K, GEN L, GEN peh2d, GEN vroots, GEN dr, GEN di,
1102 : GEN an, GEN bn)
1103 : {
1104 37294 : pari_sp av0 = avma;
1105 37294 : long m, n, q, L0 = lg(an)-1;
1106 37294 : double sig0 = rtodbl(gel(dr,1)), sub2 = rtodbl(gel(dr,2));
1107 37296 : double k1 = rtodbl(gel(dr,3)), MAXs = rtodbl(gel(dr,4));
1108 37296 : long D = di[1], M = di[2], m0 = di[3];
1109 37296 : double M0 = sig0? sub2 / sig0: 1./0.;
1110 37296 : GEN AB, A, B, vK = cgetg(M/m0 + 2, t_VEC);
1111 :
1112 191405 : for (q = 0, m = r; m <= M; m += m0, q++)
1113 154108 : gel(vK, q+1) = const_vec(L[m+1], NULL);
1114 37297 : worker_init(q, &an, &bn, &AB, &A, &B);
1115 190463 : for (m -= m0, q--; m >= 0; m -= m0, q--)
1116 : {
1117 154110 : double c1 = D + ((m > M0)? m * sig0 - sub2 : 0);
1118 154110 : GEN vKm = gel(vK,q+1); /* conceptually K(m,n) */
1119 10049684 : for (n = 1; n <= L[m+1]; n++)
1120 : {
1121 : GEN t2d, kmn;
1122 9896517 : long nn, mm, qq, p = 0;
1123 : double c, c2;
1124 : pari_sp av;
1125 :
1126 9896517 : if (gel(vKm, n)) continue; /* done already */
1127 7650392 : c = c1 + k1 * log2(n);
1128 : /* n *= 2; m -= m0 => c += c2, provided m >= M0. Else c += k1 */
1129 7650392 : c2 = k1 - MAXs;
1130 : /* p = largest (absolute) accuracy to which we need K(m,n) */
1131 11669178 : for (mm=m,nn=n; mm >= M0;)
1132 : {
1133 8126512 : if (nn <= L[mm+1] && (gel(an, nn) || (bn && gel(bn, nn))))
1134 1489010 : if (c > 0) p = maxuu(p, (ulong)c);
1135 8126637 : nn <<= 1;
1136 8126637 : mm -= m0; if (mm >= M0) c += c2; else { c += k1; break; }
1137 : }
1138 : /* mm < M0 || nn > L[mm+1] */
1139 14742123 : for ( ; mm >= 0; nn<<=1,mm-=m0,c+=k1)
1140 7091653 : if (nn <= L[mm+1] && (gel(an, nn) || (bn && gel(bn, nn))))
1141 1709173 : if (c > 0) p = maxuu(p, (ulong)c);
1142 7650470 : if (!p) continue; /* a_{n 2^v} = 0 for all v in range */
1143 2178922 : av = avma;
1144 2178922 : t2d = mpmul(gel(vroots,n), gel(peh2d,m+1));/*(n exp(mh)/sqrt(N))^(2/d)*/
1145 2179940 : kmn = gerepileupto(av, gammamellininvrt(K, t2d, p));
1146 6732068 : for (qq=q,mm=m,nn=n; mm >= 0; nn<<=1,mm-=m0,qq--)
1147 4554167 : if (nn <= L[mm+1]) gmael(vK, qq+1, nn) = kmn;
1148 : }
1149 : }
1150 190433 : for (q = 0, m = r; m <= M; m += m0, q++)
1151 : {
1152 154076 : long N = minss(L0, L[m+1]);
1153 154077 : gel(A, q+1) = an_msum(an, N, gel(vK,q+1));
1154 154079 : if (bn) gel(B, q+1) = an_msum(bn, N, gel(vK,q+1));
1155 : }
1156 36357 : return gerepileupto(av0, AB);
1157 : }
1158 : /* return A = [\theta(exp(mh)), m=0..M], theta(t) = sum a(n) K(n/sqrt(N) t),
1159 : * h = log(2)/m0. If bn != NULL, return the pair [A, B] */
1160 : static GEN
1161 5656 : lfuninit_ab(GEN theta, GEN h, struct lfunp *S)
1162 : {
1163 5656 : const long M = S->M, prec = S->precmax;
1164 5656 : GEN tech = linit_get_tech(theta), isqN = theta_get_isqrtN(tech);
1165 5656 : GEN an = S->an, bn = S->bn, va, vb;
1166 : struct pari_mt pt;
1167 : GEN worker;
1168 : long m0, r, pending;
1169 :
1170 5656 : if (S->vgaell)
1171 : { /* d=2 and Vga = [a,a+1] */
1172 2800 : GEN a = vecmin(ldata_get_gammavec(linit_get_ldata(theta)));
1173 2800 : GEN qk = gpowers0(mpexp(h), M, isqN);
1174 2800 : m0 = minss(M+1, mt_nbthreads());
1175 2800 : worker = snm_closure(is_entry("_lfuninit_theta2_worker"),
1176 : mkvecn(6, S->L, qk, a, mkvecsmall3(prec, M, m0),
1177 : an, bn? bn: gen_0));
1178 : }
1179 : else
1180 : {
1181 : GEN vroots, peh2d, d2;
1182 2856 : double sig0 = S->MAXs / S->m0, sub2 = S->sub / M_LN2;
1183 : /* For all 0<= m <= M, and all n <= L[m+1] such that a_n!=0, we compute
1184 : * k[m,n] = K(n exp(mh)/sqrt(N))
1185 : * with ln(absolute error) <= E + max(mh sigma - sub, 0) + k1 * log(n).
1186 : * N.B. we use the 'rt' variant and pass (n exp(mh)/sqrt(N))^(2/d).
1187 : * Speedup: if n' = 2n and m' = m - m0 >= 0; then k[m,n] = k[m',n']. */
1188 2856 : vroots = mkvroots(S->d, S->nmax, prec); /* vroots[n] = n^(2/d) */
1189 2856 : d2 = gdivgu(gen_2, S->d);
1190 2856 : peh2d = gpowers0(gexp(gmul(d2,h), prec), M, gpow(isqN, d2, prec));
1191 2856 : m0 = S->m0; /* peh2d[m+1] = (exp(mh)/sqrt(N))^(2/d) */
1192 2856 : worker = snm_closure(is_entry("_lfuninit_worker"),
1193 : mkvecn(8, theta_get_K(tech), S->L, peh2d, vroots,
1194 : mkvec4(dbltor(sig0), dbltor(sub2),
1195 : dbltor(S->k1), dbltor(S->MAXs)),
1196 : mkvecsmall3(S->D, M, m0),
1197 : an, bn? bn: gen_0));
1198 : /* For each 0 <= m <= M, we will sum for n<=L[m+1] a(n) K(m,n)
1199 : * bit accuracy for K(m,n): D + k1*log2(n) + 1_{m > M0} (m*sig0 - sub2)
1200 : * We restrict m to arithmetic progressions r mod m0 to save memory and
1201 : * allow parallelization */
1202 : }
1203 5656 : va = cgetg(M+2, t_VEC);
1204 5656 : vb = bn? cgetg(M+2, t_VEC): NULL;
1205 5656 : mt_queue_start_lim(&pt, worker, m0);
1206 5656 : pending = 0;
1207 73976 : for (r = 0; r < m0 || pending; r++)
1208 : { /* m = q m0 + r */
1209 : GEN done, A, B;
1210 : long q, m, workid;
1211 68320 : mt_queue_submit(&pt, r, r < m0 ? mkvec(utoi(r)): NULL);
1212 68320 : done = mt_queue_get(&pt, &workid, &pending);
1213 68320 : if (!done) continue;
1214 51807 : if (bn) { A = gel(done,1); B = gel(done,2); } else { A = done; B = NULL; }
1215 347309 : for (q = 0, m = workid; m <= M; m += m0, q++)
1216 : {
1217 295502 : gel(va, m+1) = gel(A, q+1);
1218 295502 : if (bn) gel(vb, m+1) = gel(B, q+1);
1219 : }
1220 : }
1221 5656 : mt_queue_end(&pt);
1222 5656 : return bn? mkvec2(va, vb): va;
1223 : }
1224 :
1225 : static void
1226 93132 : parse_dom(double k, GEN dom, struct lfunp *S)
1227 : {
1228 93132 : long l = lg(dom);
1229 93132 : if (typ(dom)!=t_VEC) pari_err_TYPE("lfuninit [domain]", dom);
1230 93132 : if (l == 2)
1231 : {
1232 26125 : S->dc = k/2.;
1233 26125 : S->dw = 0.;
1234 26125 : S->dh = gtodouble(gel(dom,1));
1235 : }
1236 67007 : else if (l == 3)
1237 : {
1238 301 : S->dc = k/2.;
1239 301 : S->dw = gtodouble(gel(dom,1));
1240 301 : S->dh = gtodouble(gel(dom,2));
1241 : }
1242 66706 : else if (l == 4)
1243 : {
1244 66706 : S->dc = gtodouble(gel(dom,1));
1245 66706 : S->dw = gtodouble(gel(dom,2));
1246 66706 : S->dh = gtodouble(gel(dom,3));
1247 : }
1248 : else
1249 : {
1250 0 : pari_err_TYPE("lfuninit [domain]", dom);
1251 0 : S->dc = S->dw = S->dh = 0; /*-Wall*/
1252 : }
1253 93132 : if (S->dw < 0 || S->dh < 0) pari_err_TYPE("lfuninit [domain]", dom);
1254 93132 : }
1255 :
1256 : /* do we have dom \subset dom0 ? dom = [center, width, height] */
1257 : int
1258 16084 : sdomain_isincl(double k, GEN dom, GEN dom0)
1259 : {
1260 : struct lfunp S0, S;
1261 16084 : parse_dom(k, dom, &S);
1262 16084 : parse_dom(k, dom0, &S0);
1263 16084 : return S0.dc - S0.dw <= S.dc - S.dw
1264 16084 : && S0.dc + S0.dw >= S.dc + S.dw && S0.dh >= S.dh;
1265 : }
1266 :
1267 : static int
1268 16084 : checklfuninit(GEN linit, GEN dom, long der, long bitprec)
1269 : {
1270 16084 : GEN ldata = linit_get_ldata(linit);
1271 16084 : GEN domain = lfun_get_domain(linit_get_tech(linit));
1272 16084 : return domain_get_der(domain) >= der
1273 16084 : && domain_get_bitprec(domain) >= bitprec
1274 32168 : && sdomain_isincl(gtodouble(ldata_get_k(ldata)), dom, domain_get_dom(domain));
1275 : }
1276 :
1277 : static GEN
1278 1043 : ginvsqrtvec(GEN x, long prec)
1279 : {
1280 1043 : if (is_vec_t(typ(x)))
1281 1064 : pari_APPLY_same(ginv(gsqrt(gel(x,i), prec)))
1282 742 : else return ginv(gsqrt(x, prec));
1283 : }
1284 :
1285 : GEN
1286 6489 : lfuninit_make(long t, GEN ldata, GEN tech, GEN domain)
1287 : {
1288 6489 : GEN Vga = ldata_get_gammavec(ldata);
1289 6489 : long d = lg(Vga)-1;
1290 6489 : GEN w2 = gen_1, k2 = gmul2n(ldata_get_k(ldata), -1);
1291 6489 : GEN expot = gdivgu(gadd(gmulsg(d, gsubgs(k2, 1)), sumVga(Vga)), 4);
1292 6489 : if (typ(ldata_get_dual(ldata))==t_INT)
1293 : {
1294 6335 : GEN eno = ldata_get_rootno(ldata);
1295 6335 : long prec = nbits2prec( domain_get_bitprec(domain) );
1296 6335 : if (!isint1(eno)) w2 = ginvsqrtvec(eno, prec);
1297 : }
1298 6489 : tech = mkvec3(domain, tech, mkvec4(k2, w2, expot, gammafactor(Vga)));
1299 6489 : return mkvec3(mkvecsmall(t), ldata, tech);
1300 : }
1301 :
1302 : static void
1303 2856 : lfunparams2(struct lfunp *S)
1304 : {
1305 2856 : GEN L = S->L, an = S->an, bn = S->bn;
1306 : double pmax;
1307 2856 : long m, nan, nmax, neval, M = S->M;
1308 :
1309 2856 : S->vgaell = 0;
1310 : /* try to reduce parameters now we know the a_n (some may be 0) */
1311 2856 : if (typ(an) == t_VEC) an = RgV_kill0(an);
1312 2856 : nan = S->nmax; /* lg(an)-1 may be large than this */
1313 2856 : nmax = neval = 0;
1314 2856 : if (!bn)
1315 155936 : for (m = 0; m <= M; m++)
1316 : {
1317 153101 : long n = minss(nan, L[m+1]);
1318 231025 : while (n > 0 && !gel(an,n)) n--;
1319 153101 : if (n > nmax) nmax = n;
1320 153101 : neval += n;
1321 153101 : L[m+1] = n; /* reduce S->L[m+1] */
1322 : }
1323 : else
1324 : {
1325 21 : if (typ(bn) == t_VEC) bn = RgV_kill0(bn);
1326 1036 : for (m = 0; m <= M; m++)
1327 : {
1328 1015 : long n = minss(nan, L[m+1]);
1329 1057 : while (n > 0 && !gel(an,n) && !gel(bn,n)) n--;
1330 1015 : if (n > nmax) nmax = n;
1331 1015 : neval += n;
1332 1015 : L[m+1] = n; /* reduce S->L[m+1] */
1333 : }
1334 : }
1335 2856 : if (DEBUGLEVEL >= 1) err_printf("expected evaluations: %ld\n", neval);
1336 2856 : for (; M > 0; M--)
1337 2856 : if (L[M+1]) break;
1338 2856 : setlg(L, M+2);
1339 2856 : S->M = M;
1340 2856 : S->nmax = nmax;
1341 :
1342 : /* need K(n*exp(mh)/sqrt(N)) to absolute accuracy
1343 : * D + k1*log(n) + max(m * sig0 - sub2, 0) */
1344 2856 : pmax = S->D + S->k1 * log2(L[1]);
1345 2856 : if (S->MAXs)
1346 : {
1347 2856 : double sig0 = S->MAXs/S->m0, sub2 = S->sub / M_LN2;
1348 130428 : for (m = ceil(sub2 / sig0); m <= S->M; m++)
1349 : {
1350 127572 : double c = S->D + m*sig0 - sub2;
1351 127572 : if (S->k1 > 0) c += S->k1 * log2(L[m+1]);
1352 127572 : pmax = maxdd(pmax, c);
1353 : }
1354 : }
1355 2856 : S->Dmax = pmax;
1356 2856 : S->precmax = nbits2prec(pmax);
1357 2856 : }
1358 :
1359 : static GEN
1360 5670 : lfun_init_theta(GEN ldata, GEN eno, struct lfunp *S)
1361 : {
1362 5670 : GEN an2, dual, tdom = NULL, Vga = ldata_get_gammavec(ldata);
1363 5670 : long L, prec = S->precmax;
1364 5670 : if (eno)
1365 4004 : L = S->nmax;
1366 : else
1367 : {
1368 1666 : tdom = dbltor(sqrt(0.5));
1369 1666 : L = maxss(S->nmax, lfunthetacost(ldata, tdom, 0, S->D));
1370 : }
1371 5670 : dual = ldata_get_dual(ldata);
1372 5670 : S->an = ldata_vecan(ldata_get_an(ldata), L, prec);
1373 5656 : S->bn = typ(dual)==t_INT? NULL: ldata_vecan(dual, S->nmax, prec);
1374 5656 : if (!vgaell(Vga)) lfunparams2(S);
1375 : else
1376 : {
1377 2800 : S->an = antwist(S->an, Vga, prec);
1378 2800 : if (S->bn) S->bn = antwist(S->bn, Vga, prec);
1379 2800 : S->vgaell = 1;
1380 : }
1381 5656 : an2 = lg(Vga)-1 == 1? antwist(S->an, Vga, prec): S->an;
1382 5656 : return lfunthetainit0(ldata, tdom, an2, 0, S->Dmax, 0);
1383 : }
1384 :
1385 : GEN
1386 12098 : lfuncost(GEN L, GEN dom, long der, long bit)
1387 : {
1388 12098 : pari_sp av = avma;
1389 12098 : GEN ldata = lfunmisc_to_ldata_shallow(L);
1390 12098 : GEN w, k = ldata_get_k(ldata);
1391 : struct lfunp S;
1392 :
1393 12098 : parse_dom(gtodouble(k), dom, &S);
1394 12098 : lfunp_set(ldata, der, bit, &S);
1395 12098 : w = ldata_get_rootno(ldata);
1396 12098 : if (isintzero(w)) /* for lfunrootres */
1397 7 : S.nmax = maxss(S.nmax, lfunthetacost(ldata, dbltor(sqrt(0.5)), 0, bit+1));
1398 12098 : set_avma(av); return mkvecsmall2(S.nmax, S.Dmax);
1399 : }
1400 : GEN
1401 49 : lfuncost0(GEN L, GEN dom, long der, long bitprec)
1402 : {
1403 49 : pari_sp av = avma;
1404 : GEN C;
1405 :
1406 49 : if (is_linit(L))
1407 : {
1408 28 : GEN tech = linit_get_tech(L);
1409 28 : GEN domain = lfun_get_domain(tech);
1410 28 : dom = domain_get_dom(domain);
1411 28 : der = domain_get_der(domain);
1412 28 : bitprec = domain_get_bitprec(domain);
1413 28 : if (linit_get_type(L) == t_LDESC_PRODUCT)
1414 : {
1415 21 : GEN v = lfunprod_get_fact(linit_get_tech(L)), F = gel(v,1);
1416 21 : long i, l = lg(F);
1417 21 : C = cgetg(l, t_VEC);
1418 70 : for (i = 1; i < l; ++i)
1419 49 : gel(C, i) = zv_to_ZV( lfuncost(gel(F,i), dom, der, bitprec) );
1420 21 : return gerepileupto(av, C);
1421 : }
1422 : }
1423 28 : if (!dom) pari_err_TYPE("lfuncost [missing s domain]", L);
1424 28 : C = lfuncost(L,dom,der,bitprec);
1425 28 : return gerepileupto(av, zv_to_ZV(C));
1426 : }
1427 :
1428 : GEN
1429 22412 : lfuninit(GEN lmisc, GEN dom, long der, long bitprec)
1430 : {
1431 22412 : pari_sp av = avma;
1432 : GEN poqk, AB, R, h, theta, ldata, eno, r, domain, tech, k;
1433 : struct lfunp S;
1434 :
1435 22412 : if (is_linit(lmisc))
1436 : {
1437 16126 : long t = linit_get_type(lmisc);
1438 16126 : if (t==t_LDESC_INIT || t==t_LDESC_PRODUCT)
1439 : {
1440 16084 : if (checklfuninit(lmisc, dom, der, bitprec)) return lmisc;
1441 21 : pari_warn(warner,"lfuninit: insufficient initialization");
1442 : }
1443 : }
1444 6349 : ldata = lfunmisc_to_ldata_shallow(lmisc);
1445 :
1446 6349 : switch (ldata_get_type(ldata))
1447 : {
1448 595 : case t_LFUN_NF:
1449 : {
1450 595 : GEN T = gel(ldata_get_an(ldata), 2);
1451 595 : return gerepilecopy(av, lfunzetakinit(T, dom, der, bitprec));
1452 : }
1453 84 : case t_LFUN_ABELREL:
1454 : {
1455 84 : GEN T = gel(ldata_get_an(ldata), 2);
1456 84 : return gerepilecopy(av, lfunabelianrelinit(gel(T,1), gel(T,2), dom, der, bitprec));
1457 : }
1458 : }
1459 5670 : k = ldata_get_k(ldata);
1460 5670 : parse_dom(gtodouble(k), dom, &S);
1461 5670 : lfunp_set(ldata, der, bitprec, &S);
1462 5670 : ldata = ldata_newprec(ldata, nbits2prec(S.Dmax));
1463 5670 : r = ldata_get_residue(ldata);
1464 : /* Note: all guesses should already have been performed (thetainit more
1465 : * expensive than needed: should be either tdom = 1 or bitprec = S.D).
1466 : * BUT if the root number / polar part do not have an algebraic
1467 : * expression, there is no way to do this until we know the
1468 : * precision, i.e. now. So we can't remove guessing code from here and
1469 : * lfun_init_theta */
1470 5670 : if (r && isintzero(r)) eno = NULL;
1471 : else
1472 : {
1473 5670 : eno = ldata_get_rootno(ldata);
1474 5670 : if (isintzero(eno)) eno = NULL;
1475 : }
1476 5670 : theta = lfun_init_theta(ldata, eno, &S);
1477 5656 : if (eno && !r)
1478 2338 : R = gen_0;
1479 : else
1480 : {
1481 3318 : GEN v = lfunrootres(theta, S.D);
1482 3318 : ldata = shallowcopy(ldata);
1483 3318 : gel(ldata, 6) = gel(v,3);
1484 3318 : r = gel(v,1);
1485 3318 : R = gel(v,2);
1486 3318 : if (isintzero(r)) setlg(ldata,7); else gel(ldata, 7) = r;
1487 : }
1488 5656 : h = divru(mplog2(S.precmax), S.m0);
1489 : /* exp(kh/2 . [0..M]) */
1490 5656 : poqk = gequal0(k) ? NULL
1491 5656 : : gpowers(gprec_w(mpexp(gmul2n(gmul(k,h), -1)), S.precmax), S.M);
1492 5656 : AB = lfuninit_ab(theta, h, &S);
1493 5656 : if (S.bn)
1494 : {
1495 154 : GEN A = gel(AB,1), B = gel(AB,2);
1496 154 : A = lfuninit_pol(A, poqk, S.precmax);
1497 154 : B = lfuninit_pol(B, poqk, S.precmax);
1498 154 : AB = mkvec2(A, B);
1499 : }
1500 : else
1501 5502 : AB = lfuninit_pol(AB, poqk, S.precmax);
1502 5656 : tech = mkvec3(h, AB, R);
1503 5656 : domain = mkvec2(dom, mkvecsmall2(der, bitprec));
1504 5656 : return gerepilecopy(av, lfuninit_make(t_LDESC_INIT, ldata, tech, domain));
1505 : }
1506 :
1507 : GEN
1508 483 : lfuninit0(GEN lmisc, GEN dom, long der, long bitprec)
1509 : {
1510 483 : GEN z = lfuninit(lmisc, dom, der, bitprec);
1511 483 : return z == lmisc? gcopy(z): z;
1512 : }
1513 :
1514 : /* If s is a pole of Lambda, return polar part at s; else return NULL */
1515 : static GEN
1516 5066 : lfunpoleresidue(GEN R, GEN s)
1517 : {
1518 : long j;
1519 14428 : for (j = 1; j < lg(R); j++)
1520 : {
1521 9915 : GEN Rj = gel(R, j), be = gel(Rj, 1);
1522 9915 : if (gequal(s, be)) return gel(Rj, 2);
1523 : }
1524 4513 : return NULL;
1525 : }
1526 :
1527 : /* Compute contribution of polar part at s when not a pole. */
1528 : static GEN
1529 8221 : veccothderivn(GEN a, long n)
1530 : {
1531 : long i;
1532 8221 : pari_sp av = avma;
1533 8221 : GEN c = pol_x(0), cp = mkpoln(3, gen_m1, gen_0, gen_1);
1534 8221 : GEN v = cgetg(n+2, t_VEC);
1535 8221 : gel(v, 1) = poleval(c, a);
1536 24782 : for(i = 2; i <= n+1; i++)
1537 : {
1538 16561 : c = ZX_mul(ZX_deriv(c), cp);
1539 16561 : gel(v, i) = gdiv(poleval(c, a), mpfact(i-1));
1540 : }
1541 8221 : return gerepilecopy(av, v);
1542 : }
1543 :
1544 : static GEN
1545 8340 : polepart(long n, GEN h, GEN C)
1546 : {
1547 8340 : GEN h2n = gpowgs(gdiv(h, gen_2), n-1);
1548 8340 : GEN res = gmul(h2n, gel(C,n));
1549 8340 : return odd(n)? res : gneg(res);
1550 : }
1551 :
1552 : static GEN
1553 4079 : lfunsumcoth(GEN R, GEN s, GEN h, long prec)
1554 : {
1555 : long i,j;
1556 4079 : GEN S = gen_0;
1557 12300 : for (j = 1; j < lg(R); ++j)
1558 : {
1559 8221 : GEN r = gel(R,j), be = gel(r,1), Rj = gel(r, 2);
1560 8221 : long e = valser(Rj);
1561 8221 : GEN z1 = gexpm1(gmul(h, gsub(s,be)), prec); /* exp(h(s-beta))-1 */
1562 8221 : GEN c1 = gaddgs(gdivsg(2, z1), 1); /* coth((h/2)(s-beta)) */
1563 8221 : GEN C1 = veccothderivn(c1, 1-e);
1564 16561 : for (i = e; i < 0; i++)
1565 : {
1566 8340 : GEN Rbe = mysercoeff(Rj, i);
1567 8340 : GEN p1 = polepart(-i, h, C1);
1568 8340 : S = gadd(S, gmul(Rbe, p1));
1569 : }
1570 : }
1571 4079 : return gmul2n(S, -1);
1572 : }
1573 :
1574 : static GEN lfunlambda_OK(GEN linit, GEN s, GEN sdom, long bitprec);
1575 : /* L is a t_LDESC_PRODUCT Linit */
1576 : static GEN
1577 1426 : lfunlambda_product(GEN L, GEN s, GEN sdom, long bitprec)
1578 : {
1579 1426 : GEN ldata = linit_get_ldata(L), v = lfunprod_get_fact(linit_get_tech(L));
1580 1426 : GEN r = gen_1, F = gel(v,1), E = gel(v,2), C = gel(v,3), cs = conj_i(s);
1581 1426 : long i, l = lg(F), isreal = gequal(imag_i(s), imag_i(cs));
1582 4677 : for (i = 1; i < l; ++i)
1583 : {
1584 3251 : GEN f = lfunlambda_OK(gel(F, i), s, sdom, bitprec);
1585 3251 : if( DEBUGLEVEL>=2) err_printf("lfunlambda(%ld): %Ps\n",i,f);
1586 3251 : if (typ(f)==t_VEC) f = RgV_prod(f);
1587 3251 : if (E[i]) r = gmul(r, gpowgs(f, E[i]));
1588 3251 : if (C[i])
1589 : {
1590 0 : GEN fc = isreal? f: conj_i(lfunlambda_OK(gel(F, i), cs, sdom, bitprec));
1591 0 : r = gmul(r, gpowgs(fc, C[i]));
1592 : }
1593 : }
1594 1426 : return (ldata_isreal(ldata) && gequal0(imag_i(s)))? real_i(r): r;
1595 : }
1596 :
1597 : /* s a t_SER */
1598 : static long
1599 1401 : der_level(GEN s)
1600 1401 : { return signe(s)? lg(s)-3: valser(s)-1; }
1601 :
1602 : /* s a t_SER; return coeff(s, X^0) */
1603 : static GEN
1604 385 : ser_coeff0(GEN s) { return simplify_shallow(polcoef_i(s, 0, -1)); }
1605 :
1606 : static GEN
1607 9338 : get_domain(GEN s, GEN *dom, long *der)
1608 : {
1609 9338 : GEN sa = s;
1610 9338 : *der = 0;
1611 9338 : switch(typ(s))
1612 : {
1613 7 : case t_POL:
1614 7 : case t_RFRAC: s = toser_i(s);
1615 385 : case t_SER:
1616 385 : *der = der_level(s);
1617 385 : sa = ser_coeff0(s);
1618 : }
1619 9338 : *dom = mkvec3(real_i(sa), gen_0, gabs(imag_i(sa),DEFAULTPREC));
1620 9338 : return s;
1621 : }
1622 : /* assume s went through get_domain and s/bitprec belong to domain */
1623 : static GEN
1624 23024 : lfunlambda_OK(GEN linit, GEN s, GEN sdom, long bitprec)
1625 : {
1626 : GEN eno, ldata, tech, h, pol;
1627 23024 : GEN S, S0 = NULL, k2, cost;
1628 : long prec, prec0;
1629 : struct lfunp D, D0;
1630 :
1631 23024 : if (linit_get_type(linit) == t_LDESC_PRODUCT)
1632 1426 : return lfunlambda_product(linit, s, sdom, bitprec);
1633 21598 : ldata = linit_get_ldata(linit);
1634 21598 : eno = ldata_get_rootno(ldata);
1635 21598 : tech = linit_get_tech(linit);
1636 21598 : h = lfun_get_step(tech); prec = realprec(h);
1637 : /* try to reduce accuracy */
1638 21598 : parse_dom(0, sdom, &D0);
1639 21598 : parse_dom(0, domain_get_dom(lfun_get_domain(tech)), &D);
1640 21598 : if (0.8 * D.dh > D0.dh)
1641 : {
1642 12021 : cost = lfuncost(linit, sdom, typ(s)==t_SER? der_level(s): 0, bitprec);
1643 12021 : prec0 = nbits2prec(cost[2]);
1644 12021 : if (prec0 < prec) { prec = prec0; h = gprec_w(h, prec); }
1645 : }
1646 21598 : pol = lfun_get_pol(tech);
1647 21598 : s = gprec_w(s, prec);
1648 21598 : if (ldata_get_residue(ldata))
1649 : {
1650 4478 : GEN R = lfun_get_Residue(tech);
1651 4478 : GEN Ra = lfunpoleresidue(R, s);
1652 4478 : if (Ra) return gprec_w(Ra, nbits2prec(bitprec));
1653 4079 : S0 = lfunsumcoth(R, s, h, prec);
1654 : }
1655 21199 : k2 = lfun_get_k2(tech);
1656 21199 : if (typ(pol)==t_POL && typ(s) != t_SER && gequal(real_i(s), k2))
1657 15102 : { /* on critical line: shortcut */
1658 15102 : GEN polz, b = imag_i(s);
1659 15102 : polz = gequal0(b)? poleval(pol,gen_1): poleval(pol, expIr(gmul(h,b)));
1660 15102 : S = gadd(polz, gmulvec(eno, conj_i(polz)));
1661 : }
1662 : else
1663 : {
1664 6097 : GEN z = gexp(gmul(h, gsub(s, k2)), prec);
1665 6097 : GEN zi = ginv(z), zc = conj_i(zi);
1666 6097 : if (typ(pol)==t_POL)
1667 5901 : S = gadd(poleval(pol, z), gmulvec(eno, conj_i(poleval(pol, zc))));
1668 : else
1669 196 : S = gadd(poleval(gel(pol,1), z), gmulvec(eno, poleval(gel(pol,2), zi)));
1670 : }
1671 21199 : if (S0) S = gadd(S,S0);
1672 21199 : return gprec_w(gmul(S,h), nbits2prec(bitprec));
1673 : }
1674 : GEN
1675 1008 : lfunlambda(GEN lmisc, GEN s, long bitprec)
1676 : {
1677 1008 : pari_sp av = avma;
1678 : GEN linit, dom, z;
1679 : long der;
1680 1008 : s = get_domain(s, &dom, &der);
1681 1008 : linit = lfuninit(lmisc, dom, der, bitprec);
1682 1008 : z = lfunlambda_OK(linit,s, dom, bitprec);
1683 1008 : return gerepilecopy(av, z);
1684 : }
1685 :
1686 : static long
1687 8820 : is_ser(GEN x)
1688 : {
1689 8820 : long t = typ(x);
1690 8820 : if (t == t_SER) return 1;
1691 6874 : if (!is_vec_t(t) || lg(x)==1) return 0;
1692 350 : if (typ(gel(x,1))==t_SER) return 1;
1693 252 : return 0;
1694 : }
1695 :
1696 : static GEN
1697 357 : lfunser(GEN L)
1698 : {
1699 357 : long v = valser(L);
1700 357 : if (v > 0) return gen_0;
1701 315 : if (v == 0) L = gel(L, 2);
1702 : else
1703 189 : setlg(L, minss(lg(L), 2-v));
1704 315 : return L;
1705 : }
1706 :
1707 : static GEN
1708 357 : lfunservec(GEN x)
1709 : {
1710 357 : if (typ(x)==t_SER) return lfunser(x);
1711 0 : pari_APPLY_same(lfunser(gel(x,i)))
1712 : }
1713 : static GEN
1714 105 : lfununext(GEN L)
1715 : {
1716 105 : setlg(L, maxss(lg(L)-1, valser(L)? 2: 3));
1717 105 : return normalizeser(L);
1718 : }
1719 : static GEN
1720 105 : lfununextvec(GEN x)
1721 : {
1722 105 : if (typ(x)==t_SER) return lfununext(x);
1723 0 : pari_APPLY_same(lfununext(gel(x,i)));
1724 : }
1725 :
1726 : /* assume lmisc is an linit, s went through get_domain and s/bitprec belong
1727 : * to domain */
1728 : static GEN
1729 5033 : lfun_OK(GEN linit, GEN s, GEN sdom, long bitprec)
1730 : {
1731 5033 : GEN N, gas, S, FVga, res, ss = s;
1732 5033 : long prec = nbits2prec(bitprec), ext;
1733 :
1734 5033 : FVga = lfun_get_factgammavec(linit_get_tech(linit));
1735 5033 : S = lfunlambda_OK(linit, s, sdom, bitprec);
1736 5033 : if (is_ser(S))
1737 : {
1738 1624 : GEN r = typ(S)==t_SER ? S : gel(S,1);
1739 1624 : long d = lg(r) - 2 + fracgammadegree(FVga);
1740 1624 : if (typ(s) == t_SER)
1741 1309 : ss = sertoser(s, d);
1742 : else
1743 315 : ss = deg1ser_shallow(gen_1, s, varn(r), d);
1744 : }
1745 5033 : gas = gammafactproduct(FVga, ss, &ext, prec);
1746 5033 : N = ldata_get_conductor(linit_get_ldata(linit));
1747 5033 : res = gdiv(S, gmul(gpow(N, gdivgu(ss, 2), prec), gas));
1748 5033 : if (typ(s) != t_SER && is_ser(res)) res = lfunservec(res);
1749 4676 : else if (ext) res = lfununextvec(res);
1750 5033 : return gprec_w(res, prec);
1751 : }
1752 :
1753 : GEN
1754 6762 : lfun(GEN lmisc, GEN s, long bitprec)
1755 : {
1756 6762 : pari_sp av = avma;
1757 : GEN linit, dom, z;
1758 : long der;
1759 6762 : s = get_domain(s, &dom, &der);
1760 6762 : if (!der && typ(s) == t_INT && !is_bigint(s))
1761 : { /* special value ? */
1762 : GEN ldata;
1763 5327 : long t, ss = itos(s);
1764 5327 : if (is_linit(lmisc))
1765 616 : ldata = linit_get_ldata(lmisc);
1766 : else
1767 4711 : lmisc = ldata = lfunmisc_to_ldata_shallow(lmisc);
1768 5320 : t = ldata_get_type(ldata);
1769 5320 : if (t == t_LFUN_KRONECKER || t == t_LFUN_ZETA)
1770 : {
1771 2744 : long D = itos_or_0(gel(ldata_get_an(ldata), 2));
1772 2744 : if (D)
1773 : {
1774 2744 : if (ss <= 0) return lfunquadneg(D, ss);
1775 : /* ss > 0 */
1776 238 : if ((!odd(ss) && D > 0) || (odd(ss) && D < 0))
1777 : {
1778 168 : long prec = nbits2prec(bitprec), q = labs(D);
1779 168 : ss = 1 - ss; /* <= 0 */
1780 168 : z = powrs(divrs(mppi(prec + EXTRAPREC64), q), 1-ss);
1781 168 : z = mulrr(shiftr(z, -ss), sqrtr_abs(utor(q, prec)));
1782 168 : z = gdiv(z, mpfactr(-ss, prec));
1783 168 : if (smodss(ss, 4) > 1) togglesign(z);
1784 168 : return gmul(z, lfunquadneg(D, ss));
1785 : }
1786 : }
1787 : }
1788 : }
1789 4081 : linit = lfuninit(lmisc, dom, der, bitprec);
1790 4067 : z = lfun_OK(linit, s, dom, bitprec);
1791 4067 : return gerepilecopy(av, z);
1792 : }
1793 :
1794 : /* given a t_SER a+x*s(x), return x*s(x), shallow */
1795 : static GEN
1796 42 : sersplit1(GEN s, GEN *head)
1797 : {
1798 42 : long i, l = lg(s);
1799 : GEN y;
1800 42 : *head = simplify_shallow(mysercoeff(s, 0));
1801 42 : if (valser(s) > 0) return s;
1802 28 : y = cgetg(l-1, t_SER); y[1] = s[1];
1803 28 : setvalser(y, 1);
1804 140 : for (i=3; i < l; i++) gel(y,i-1) = gel(s,i);
1805 28 : return normalizeser(y);
1806 : }
1807 :
1808 : /* order of pole of Lambda at s (0 if regular point) */
1809 : static long
1810 2198 : lfunlambdaord(GEN linit, GEN s)
1811 : {
1812 2198 : GEN tech = linit_get_tech(linit);
1813 2198 : if (linit_get_type(linit)==t_LDESC_PRODUCT)
1814 : {
1815 280 : GEN v = lfunprod_get_fact(linit_get_tech(linit));
1816 280 : GEN F = gel(v, 1), E = gel(v, 2), C = gel(v, 3);
1817 280 : long i, ex = 0, l = lg(F);
1818 952 : for (i = 1; i < l; i++)
1819 672 : ex += lfunlambdaord(gel(F,i), s) * (E[i]+C[i]);
1820 280 : return ex;
1821 : }
1822 1918 : if (ldata_get_residue(linit_get_ldata(linit)))
1823 : {
1824 588 : GEN r = lfunpoleresidue(lfun_get_Residue(tech), s);
1825 588 : if (r) return lg(r)-2;
1826 : }
1827 1764 : return 0;
1828 : }
1829 :
1830 : static GEN
1831 126 : derser(GEN res, long m)
1832 : {
1833 126 : long v = valser(res);
1834 126 : if (v > m) return gen_0;
1835 126 : if (v >= 0)
1836 126 : return gmul(mysercoeff(res, m), mpfact(m));
1837 : else
1838 0 : return derivn(res, m, -1);
1839 : }
1840 :
1841 : static GEN
1842 189 : derservec(GEN x, long m) { pari_APPLY_same(derser(gel(x,i),m)) }
1843 :
1844 : /* derivative of order m > 0 of L (flag = 0) or Lambda (flag = 1) */
1845 : static GEN
1846 1575 : lfunderiv(GEN lmisc, long m, GEN s, long flag, long bitprec)
1847 : {
1848 1575 : pari_sp ltop = avma;
1849 1575 : GEN res, S = NULL, linit, dom;
1850 1575 : long der, prec = nbits2prec(bitprec);
1851 1575 : if (m <= 0) pari_err_DOMAIN("lfun", "D", "<=", gen_0, stoi(m));
1852 1568 : s = get_domain(s, &dom, &der);
1853 1568 : linit = lfuninit(lmisc, dom, der + m, bitprec);
1854 1568 : if (typ(s) == t_SER)
1855 : {
1856 42 : long v, l = lg(s)-1;
1857 : GEN sh;
1858 42 : if (valser(s) < 0) pari_err_DOMAIN("lfun","valuation", "<", gen_0, s);
1859 42 : S = sersplit1(s, &sh);
1860 42 : v = valser(S);
1861 42 : s = deg1ser_shallow(gen_1, sh, varn(S), m + (l+v-1)/v);
1862 : }
1863 : else
1864 : {
1865 1526 : long ex = lfunlambdaord(linit, s);
1866 : /* HACK: pretend lfuninit was done to right accuracy */
1867 1526 : if (gequal0(s)) s = gen_0;
1868 1526 : s = deg1ser_shallow(gen_1, s, 0, m+1+ex);
1869 : }
1870 1568 : res = flag ? lfunlambda_OK(linit, s, dom, bitprec):
1871 966 : lfun_OK(linit, s, dom, bitprec);
1872 1568 : if (S)
1873 42 : res = gsubst(derivn(res, m, -1), varn(S), S);
1874 1526 : else if (typ(res)==t_SER)
1875 : {
1876 1463 : long v = valser(res);
1877 1463 : if (v > m) { set_avma(ltop); return gen_0; }
1878 1456 : if (v >= 0)
1879 1330 : res = gmul(mysercoeff(res, m), mpfact(m));
1880 : else
1881 126 : res = derivn(res, m, -1);
1882 : }
1883 63 : else if (is_ser(res))
1884 63 : res = derservec(res, m);
1885 1561 : return gerepilecopy(ltop, gprec_w(res, prec));
1886 : }
1887 :
1888 : GEN
1889 1463 : lfunlambda0(GEN lmisc, GEN s, long der, long bitprec)
1890 : {
1891 602 : return der? lfunderiv(lmisc, der, s, 1, bitprec)
1892 2065 : : lfunlambda(lmisc, s, bitprec);
1893 : }
1894 :
1895 : GEN
1896 6545 : lfun0(GEN lmisc, GEN s, long der, long bitprec)
1897 : {
1898 973 : return der? lfunderiv(lmisc, der, s, 0, bitprec)
1899 7511 : : lfun(lmisc, s, bitprec);
1900 : }
1901 :
1902 : GEN
1903 13137 : lfunhardy(GEN lmisc, GEN t, long bitprec)
1904 : {
1905 13137 : pari_sp ltop = avma;
1906 13137 : long prec = nbits2prec(bitprec), d;
1907 : GEN argz, z, linit, ldata, tech, dom, w2, k2, E, h, a, k;
1908 :
1909 13137 : switch(typ(t))
1910 : {
1911 13130 : case t_INT: case t_FRAC: case t_REAL: break;
1912 7 : default: pari_err_TYPE("lfunhardy",t);
1913 : }
1914 :
1915 13130 : ldata = lfunmisc_to_ldata_shallow(lmisc);
1916 13130 : if (!is_linit(lmisc)) lmisc = ldata;
1917 13130 : k = ldata_get_k(ldata);
1918 13130 : d = ldata_get_degree(ldata);
1919 13130 : dom = mkvec3(gmul2n(k, -1), gen_0, gabs(t,LOWDEFAULTPREC));
1920 13130 : linit = lfuninit(lmisc, dom, 0, bitprec);
1921 13130 : tech = linit_get_tech(linit);
1922 13130 : w2 = lfun_get_w2(tech);
1923 13130 : k2 = lfun_get_k2(tech);
1924 13130 : E = lfun_get_expot(tech); /* 4E = d(k2 - 1) + real(vecsum(Vga)) */
1925 13130 : z = mkcomplex(k2, t);
1926 : /* more accurate than garg: k/2 in Q */
1927 13130 : argz = gequal0(k2)? Pi2n(-1, prec): gatan(gdiv(t, k2), prec);
1928 13130 : prec = precision(argz);
1929 : /* prec may have increased: don't lose accuracy if |z|^2 is exact */
1930 13130 : a = gsub(gmulsg(d, gmul(t, gmul2n(argz,-1))),
1931 : gmul(E, glog(gnorm(z),prec)));
1932 13130 : h = lfunlambda_OK(linit, z, dom, bitprec);
1933 13130 : if (!isint1(w2) && typ(ldata_get_dual(ldata))==t_INT)
1934 10773 : h = mulrealvec(h, w2);
1935 13130 : if (typ(h) == t_COMPLEX && gexpo(imag_i(h)) < -(bitprec >> 1))
1936 2308 : h = real_i(h);
1937 13130 : return gerepileupto(ltop, gmul(h, gexp(a, prec)));
1938 : }
1939 :
1940 : /* L = log(t); return \sum_{i = 0}^{v-1} R[-i-1] L^i/i! */
1941 : static GEN
1942 1806 : theta_pole_contrib(GEN R, long v, GEN L)
1943 : {
1944 1806 : GEN s = mysercoeff(R,-v);
1945 : long i;
1946 1911 : for (i = v-1; i >= 1; i--)
1947 105 : s = gadd(mysercoeff(R,-i), gdivgu(gmul(s,L), i));
1948 1806 : return s;
1949 : }
1950 : /* subtract successively rather than adding everything then subtracting.
1951 : * The polar part is "large" and suffers from cancellation: a little stabler
1952 : * this way */
1953 : static GEN
1954 4235 : theta_add_polar_part(GEN S, GEN R, GEN t, long prec)
1955 : {
1956 4235 : GEN logt = NULL;
1957 4235 : long j, l = lg(R);
1958 6041 : for (j = 1; j < l; j++)
1959 : {
1960 1806 : GEN Rj = gel(R,j), b = gel(Rj,1), Rb = gel(Rj,2);
1961 1806 : long v = -valser(Rb);
1962 1806 : if (v > 1 && !logt) logt = glog(t, prec);
1963 1806 : S = gsub(S, gmul(theta_pole_contrib(Rb,v,logt), gpow(t,b,prec)));
1964 : }
1965 4235 : return S;
1966 : }
1967 :
1968 : static long
1969 3528 : lfuncheckfeq_i(GEN theta, GEN thetad, GEN t0, GEN t0i, long bitprec)
1970 : {
1971 3528 : GEN ldata = linit_get_ldata(theta);
1972 : GEN S0, S0i, w, eno;
1973 3528 : long prec = nbits2prec(bitprec);
1974 3528 : if (thetad)
1975 70 : S0 = lfuntheta(thetad, t0, 0, bitprec);
1976 : else
1977 3458 : S0 = conj_i(lfuntheta(theta, conj_i(t0), 0, bitprec));
1978 3528 : S0i = lfuntheta(theta, t0i, 0, bitprec);
1979 :
1980 3528 : eno = ldata_get_rootno(ldata);
1981 3528 : if (ldata_get_residue(ldata))
1982 : {
1983 903 : GEN R = theta_get_R(linit_get_tech(theta));
1984 903 : if (gequal0(R))
1985 : {
1986 : GEN v, r;
1987 105 : long t = ldata_get_type(ldata);
1988 105 : if (t == t_LFUN_NF || t == t_LFUN_ABELREL)
1989 : { /* inefficient since theta not needed; no need to optimize for this
1990 : (artificial) query [e.g. lfuncheckfeq(t_POL)] */
1991 42 : GEN L = lfuninit(ldata,zerovec(3),0,bitprec);
1992 42 : return lfuncheckfeq(L,t0,bitprec);
1993 : }
1994 63 : v = lfunrootres(theta, bitprec);
1995 63 : r = gel(v,1);
1996 63 : if (gequal0(eno)) eno = gel(v,3);
1997 63 : R = lfunrtoR_i(ldata, r, eno, nbits2prec(bitprec));
1998 : }
1999 861 : S0i = theta_add_polar_part(S0i, R, t0, prec);
2000 : }
2001 3486 : if (gequal0(S0i) || gequal0(S0)) pari_err_PREC("lfuncheckfeq");
2002 :
2003 3486 : w = gdivvec(S0i, gmul(S0, gpow(t0, ldata_get_k(ldata), prec)));
2004 : /* missing rootno: guess it */
2005 3486 : if (gequal0(eno)) eno = lfunrootno(theta, bitprec);
2006 3486 : w = gsubvec(w, eno);
2007 3486 : if (thetad) w = gdivvec(w, eno); /* |eno| may be large in non-dual case */
2008 3486 : return gexpo(w);
2009 : }
2010 :
2011 : /* Check whether the coefficients, conductor, weight, polar part and root
2012 : * number are compatible with the functional equation at t0 and 1/t0.
2013 : * Different from lfunrootres. */
2014 : long
2015 3647 : lfuncheckfeq(GEN lmisc, GEN t0, long bitprec)
2016 : {
2017 : GEN ldata, theta, thetad, t0i;
2018 : pari_sp av;
2019 :
2020 3647 : if (is_linit(lmisc) && linit_get_type(lmisc)==t_LDESC_PRODUCT)
2021 : {
2022 154 : GEN v = lfunprod_get_fact(linit_get_tech(lmisc)), F = gel(v,1);
2023 154 : long i, b = -bitprec, l = lg(F);
2024 518 : for (i = 1; i < l; i++) b = maxss(b, lfuncheckfeq(gel(F,i), t0, bitprec));
2025 154 : return b;
2026 : }
2027 3493 : av = avma;
2028 3493 : if (!t0)
2029 : { /* ~Pi/3 + I/7, some random complex number */
2030 3318 : t0 = mkcomplex(uutoQ(355,339), uutoQ(1,7));
2031 3318 : t0i = ginv(t0);
2032 : }
2033 175 : else if (gcmpgs(gnorm(t0), 1) < 0) { t0i = t0; t0 = ginv(t0); }
2034 119 : else t0i = ginv(t0);
2035 : /* |t0| >= 1 */
2036 3493 : theta = lfunthetacheckinit(lmisc, t0i, 0, bitprec);
2037 3486 : ldata = linit_get_ldata(theta);
2038 3486 : thetad = theta_dual(theta, ldata_get_dual(ldata));
2039 3486 : return gc_long(av, lfuncheckfeq_i(theta, thetad, t0, t0i, bitprec));
2040 : }
2041 :
2042 : /*******************************************************************/
2043 : /* Compute root number and residues */
2044 : /*******************************************************************/
2045 : /* round root number to \pm 1 if close to integer. */
2046 : static GEN
2047 3724 : ropm1(GEN w, long prec)
2048 : {
2049 : long e;
2050 : GEN r;
2051 3724 : if (typ(w) == t_INT) return w;
2052 3332 : r = grndtoi(w, &e);
2053 3332 : return (e < -prec2nbits(prec)/2)? r: w;
2054 : }
2055 :
2056 : /* theta for t=1/sqrt(2) and t2==2t simultaneously, saving 25% of the work.
2057 : * Assume correct initialization (no thetacheck) */
2058 : static void
2059 406 : lfunthetaspec(GEN linit, long bitprec, GEN *pv, GEN *pv2)
2060 : {
2061 406 : pari_sp av = avma, av2;
2062 : GEN t, Vga, an, K, ldata, thetainit, v, v2, vroots;
2063 : long L, prec, n, d;
2064 :
2065 406 : ldata = linit_get_ldata(linit);
2066 406 : thetainit = linit_get_tech(linit);
2067 406 : prec = nbits2prec(bitprec);
2068 406 : Vga = ldata_get_gammavec(ldata); d = lg(Vga)-1;
2069 406 : if (Vgaeasytheta(Vga))
2070 : {
2071 196 : GEN v2 = sqrtr(real2n(1, nbits2prec(bitprec)));
2072 196 : GEN v = shiftr(v2,-1);
2073 196 : *pv = lfuntheta(linit, v, 0, bitprec);
2074 196 : *pv2= lfuntheta(linit, v2, 0, bitprec);
2075 196 : return;
2076 : }
2077 210 : an = RgV_kill0( theta_get_an(thetainit) );
2078 210 : L = lg(an)-1;
2079 : /* to compute theta(1/sqrt(2)) */
2080 210 : t = ginv(gsqrt(gmul2n(ldata_get_conductor(ldata), 1), prec));
2081 : /* t = 1/sqrt(2N) */
2082 :
2083 : /* From then on, the code is generic and could be used to compute
2084 : * theta(t) / theta(2t) without assuming t = 1/sqrt(2) */
2085 210 : K = theta_get_K(thetainit);
2086 210 : vroots = mkvroots(d, L, prec);
2087 210 : t = gpow(t, gdivgu(gen_2, d), prec); /* rt variant: t->t^(2/d) */
2088 : /* v = \sum_{n <= L, n odd} a_n K(nt) */
2089 1815212 : for (v = gen_0, n = 1; n <= L; n+=2)
2090 : {
2091 1815002 : GEN tn, Kn, a = gel(an, n);
2092 :
2093 1815002 : if (!a) continue;
2094 113729 : av2 = avma;
2095 113729 : tn = gmul(t, gel(vroots,n));
2096 113729 : Kn = gammamellininvrt(K, tn, bitprec);
2097 113729 : v = gerepileupto(av2, gadd(v, gmul(a,Kn)));
2098 : }
2099 : /* v += \sum_{n <= L, n even} a_n K(nt), v2 = \sum_{n <= L/2} a_n K(2n t) */
2100 1815114 : for (v2 = gen_0, n = 1; n <= L/2; n++)
2101 : {
2102 1814904 : GEN t2n, K2n, a = gel(an, n), a2 = gel(an,2*n);
2103 :
2104 1814904 : if (!a && !a2) continue;
2105 120484 : av2 = avma;
2106 120484 : t2n = gmul(t, gel(vroots,2*n));
2107 120484 : K2n = gerepileupto(av2, gammamellininvrt(K, t2n, bitprec));
2108 120484 : if (a) v2 = gadd(v2, gmul(a, K2n));
2109 120484 : if (a2) v = gadd(v, gmul(a2,K2n));
2110 : }
2111 210 : *pv = v;
2112 210 : *pv2 = v2;
2113 210 : gerepileall(av, 2, pv,pv2);
2114 : }
2115 :
2116 : static GEN
2117 399 : Rtor(GEN a, GEN R, GEN ldata, long prec)
2118 : {
2119 399 : GEN FVga = gammafactor(ldata_get_gammavec(ldata));
2120 399 : GEN Na = gpow(ldata_get_conductor(ldata), gdivgu(a,2), prec);
2121 : long ext;
2122 399 : return gdiv(R, gmul(Na, gammafactproduct(FVga, a, &ext, prec)));
2123 : }
2124 :
2125 : /* v = theta~(t), vi = theta(1/t) */
2126 : static GEN
2127 3374 : get_eno(GEN R, GEN k, GEN t, GEN v, GEN vi, long vx, long bitprec, long force)
2128 : {
2129 3374 : long prec = nbits2prec(bitprec);
2130 3374 : GEN a0, a1, S = deg1pol(gmul(gpow(t,k,prec), gneg(v)), vi, vx);
2131 :
2132 3374 : S = theta_add_polar_part(S, R, t, prec);
2133 3374 : if (typ(S) != t_POL || degpol(S) != 1) return NULL;
2134 3374 : a1 = gel(S,3); if (!force && gexpo(a1) < -bitprec/4) return NULL;
2135 3325 : a0 = gel(S,2);
2136 3325 : return gdivvec(a0, gneg(a1));
2137 :
2138 : }
2139 : /* Return w using theta(1/t) - w t^k \bar{theta}(t) = polar_part(t,w).
2140 : * The full Taylor development of L must be known */
2141 : GEN
2142 3325 : lfunrootno(GEN linit, long bitprec)
2143 : {
2144 : GEN ldata, t, eno, v, vi, R, thetad;
2145 3325 : long c = 0, prec = nbits2prec(bitprec), vx = fetch_var();
2146 : GEN k;
2147 : pari_sp av;
2148 :
2149 : /* initialize for t > 1/sqrt(2) */
2150 3325 : linit = lfunthetacheckinit(linit, dbltor(sqrt(0.5)), 0, bitprec);
2151 3325 : ldata = linit_get_ldata(linit);
2152 3325 : k = ldata_get_k(ldata);
2153 3339 : R = ldata_get_residue(ldata)? lfunrtoR_eno(ldata, pol_x(vx), prec)
2154 3325 : : cgetg(1, t_VEC);
2155 3325 : t = gen_1;
2156 3325 : v = lfuntheta(linit, t, 0, bitprec);
2157 3325 : thetad = theta_dual(linit, ldata_get_dual(ldata));
2158 3325 : vi = !thetad ? conj_i(v): lfuntheta(thetad, t, 0, bitprec);
2159 3325 : eno = get_eno(R,k,t,vi,v, vx, bitprec, 0);
2160 3325 : if (!eno && !thetad)
2161 : { /* t = sqrt(2), vi = theta(1/t), v = theta(t) */
2162 7 : lfunthetaspec(linit, bitprec, &vi, &v);
2163 7 : t = sqrtr(utor(2, prec));
2164 7 : eno = get_eno(R,k,t,conj_i(v),vi, vx, bitprec, 0);
2165 : }
2166 3325 : av = avma;
2167 3367 : while (!eno)
2168 : {
2169 42 : t = addsr(1, shiftr(utor(pari_rand(), prec), -2-BITS_IN_LONG));
2170 : /* t in [1,1.25[ */
2171 0 : v = thetad? lfuntheta(thetad, t, 0, bitprec)
2172 42 : : conj_i(lfuntheta(linit, t, 0, bitprec));
2173 42 : vi = lfuntheta(linit, ginv(t), 0, bitprec);
2174 42 : eno = get_eno(R,k,t,v,vi, vx, bitprec, c++ == 5);
2175 42 : set_avma(av);
2176 : }
2177 3325 : delete_var(); return ropm1(eno,prec);
2178 : }
2179 :
2180 : /* Find root number and/or residues when L-function coefficients and
2181 : conductor are known. For the moment at most a single residue allowed. */
2182 : GEN
2183 3745 : lfunrootres(GEN data, long bitprec)
2184 : {
2185 3745 : pari_sp ltop = avma;
2186 : GEN k, w, r, R, a, b, e, v, v2, be, ldata, linit;
2187 : long prec;
2188 :
2189 3745 : ldata = lfunmisc_to_ldata_shallow(data);
2190 3745 : r = ldata_get_residue(ldata);
2191 3745 : k = ldata_get_k(ldata);
2192 3745 : w = ldata_get_rootno(ldata);
2193 3745 : if (r) r = normalize_simple_pole(r, k);
2194 3745 : if (!r || residues_known(r))
2195 : {
2196 3346 : if (isintzero(w)) w = lfunrootno(data, bitprec);
2197 3346 : if (!r)
2198 1708 : r = R = gen_0;
2199 : else
2200 1638 : R = lfunrtoR_eno(ldata, w, nbits2prec(bitprec));
2201 3346 : return gerepilecopy(ltop, mkvec3(r, R, w));
2202 : }
2203 399 : linit = lfunthetacheckinit(data, dbltor(sqrt(0.5)), 0, bitprec);
2204 399 : prec = nbits2prec(bitprec);
2205 399 : if (lg(r) > 2) pari_err_IMPL("multiple poles in lfunrootres");
2206 : /* Now residue unknown, and r = [[be,0]]. */
2207 399 : be = gmael(r, 1, 1);
2208 399 : if (ldata_isreal(ldata) && gequalm1(w))
2209 0 : R = lfuntheta(linit, gen_1, 0, bitprec);
2210 : else
2211 : {
2212 399 : GEN p2k = gpow(gen_2,k,prec);
2213 399 : lfunthetaspec(linit, bitprec, &v2, &v);
2214 399 : if (gequal(gmulsg(2, be), k)) pari_err_IMPL("pole at k/2 in lfunrootres");
2215 399 : if (gequal(be, k))
2216 : {
2217 133 : a = conj_i(gsub(gmul(p2k, v), v2));
2218 133 : b = subiu(p2k, 1);
2219 133 : e = gmul(gsqrt(p2k, prec), gsub(v2, v));
2220 : }
2221 : else
2222 : {
2223 266 : GEN tk2 = gsqrt(p2k, prec);
2224 266 : GEN tbe = gpow(gen_2, be, prec);
2225 266 : GEN tkbe = gpow(gen_2, gdivgu(gsub(k, be), 2), prec);
2226 266 : a = conj_i(gsub(gmul(tbe, v), v2));
2227 266 : b = gsub(gdiv(tbe, tkbe), tkbe);
2228 266 : e = gsub(gmul(gdiv(tbe, tk2), v2), gmul(tk2, v));
2229 : }
2230 399 : if (isintzero(w))
2231 : { /* Now residue unknown, r = [[be,0]], and w unknown. */
2232 7 : GEN t0 = mkfrac(utoi(11),utoi(10));
2233 7 : GEN th1 = lfuntheta(linit, t0, 0, bitprec);
2234 7 : GEN th2 = lfuntheta(linit, ginv(t0), 0, bitprec);
2235 7 : GEN tbe = gpow(t0, gmulsg(2, be), prec);
2236 7 : GEN tkbe = gpow(t0, gsub(k, be), prec);
2237 7 : GEN tk2 = gpow(t0, k, prec);
2238 7 : GEN c = conj_i(gsub(gmul(tbe, th1), th2));
2239 7 : GEN d = gsub(gdiv(tbe, tkbe), tkbe);
2240 7 : GEN f = gsub(gmul(gdiv(tbe, tk2), th2), gmul(tk2, th1));
2241 7 : GEN D = gsub(gmul(a, d), gmul(b, c));
2242 7 : w = gdiv(gsub(gmul(d, e), gmul(b, f)), D);
2243 : }
2244 399 : w = ropm1(w, prec);
2245 399 : R = gdiv(gsub(e, gmul(a, w)), b);
2246 : }
2247 399 : r = normalize_simple_pole(Rtor(be, R, ldata, prec), be);
2248 399 : R = lfunrtoR_i(ldata, r, w, prec);
2249 399 : return gerepilecopy(ltop, mkvec3(r, R, w));
2250 : }
2251 :
2252 : /*******************************************************************/
2253 : /* Zeros */
2254 : /*******************************************************************/
2255 : struct lhardyz_t {
2256 : long bitprec, prec;
2257 : GEN linit;
2258 : };
2259 :
2260 : static GEN
2261 12507 : lfunhardyzeros(void *E, GEN t)
2262 : {
2263 12507 : struct lhardyz_t *S = (struct lhardyz_t*)E;
2264 12507 : GEN z = gprec_wensure(lfunhardy(S->linit, t, S->bitprec), S->prec);
2265 12507 : return typ(z) == t_VEC ? RgV_prod(z): z;
2266 : }
2267 :
2268 : /* initialize for computation on critical line up to height h, zero
2269 : * of order <= m */
2270 : static GEN
2271 525 : lfuncenterinit(GEN lmisc, double h, long m, long bitprec)
2272 : {
2273 525 : if (m < 0)
2274 : { /* choose a sensible default */
2275 525 : if (!is_linit(lmisc) || linit_get_type(lmisc) != t_LDESC_INIT) m = 4;
2276 : else
2277 : {
2278 448 : GEN domain = lfun_get_domain(linit_get_tech(lmisc));
2279 448 : m = domain_get_der(domain);
2280 : }
2281 : }
2282 525 : return lfuninit(lmisc, mkvec(dbltor(h)), m, bitprec);
2283 : }
2284 :
2285 : long
2286 525 : lfunorderzero(GEN lmisc, long m, long bitprec)
2287 : {
2288 525 : pari_sp ltop = avma;
2289 : GEN eno, ldata, linit, k2;
2290 : long G, c0, c, st;
2291 :
2292 525 : if (is_linit(lmisc) && linit_get_type(lmisc) == t_LDESC_PRODUCT)
2293 : {
2294 77 : GEN M = gmael(linit_get_tech(lmisc), 2,1);
2295 77 : long i, l = lg(M);
2296 259 : for (c=0, i=1; i < l; i++) c += lfunorderzero(gel(M,i), m, bitprec);
2297 77 : return c;
2298 : }
2299 448 : linit = lfuncenterinit(lmisc, 0, m, bitprec);
2300 448 : ldata = linit_get_ldata(linit);
2301 448 : eno = ldata_get_rootno(ldata);
2302 448 : k2 = gmul2n(ldata_get_k(ldata), -1);
2303 448 : G = -bitprec/2;
2304 448 : c0 = 0; st = 1;
2305 448 : if (typ(eno) == t_VEC)
2306 : {
2307 42 : long i, l = lg(eno), cnt = l-1, s = 0;
2308 42 : GEN v = zero_zv(l-1);
2309 42 : if (ldata_isreal(ldata)) st = 2;
2310 84 : for (c = c0; cnt; c += st)
2311 : {
2312 42 : GEN L = lfun0(linit, k2, c, bitprec);
2313 154 : for (i = 1; i < l; i++)
2314 : {
2315 112 : if (v[i]==0 && gexpo(gel(L,i)) > G)
2316 : {
2317 112 : v[i] = c; cnt--; s += c;
2318 : }
2319 : }
2320 : }
2321 42 : return gc_long(ltop,s);
2322 : }
2323 : else
2324 : {
2325 406 : if (ldata_isreal(ldata)) { st = 2; if (!gequal1(eno)) c0 = 1; }
2326 434 : for (c = c0;; c += st)
2327 434 : if (gexpo(lfun0(linit, k2, c, bitprec)) > G) return gc_long(ltop, c);
2328 : }
2329 : }
2330 :
2331 : /* assume T1 * T2 > 0, T1 <= T2 */
2332 : static void
2333 84 : lfunzeros_i(struct lhardyz_t *S, GEN *pw, long *ct, GEN T1, GEN T2, long d,
2334 : GEN cN, GEN pi2, GEN pi2div, long precinit, long prec)
2335 : {
2336 84 : GEN T = T1, w = *pw;
2337 84 : long W = lg(w)-1, s = gsigne(lfunhardyzeros(S, T1));
2338 : for(;;)
2339 392 : {
2340 476 : pari_sp av = avma;
2341 : GEN D, T0, z;
2342 476 : D = gcmp(T, pi2) < 0? cN
2343 476 : : gadd(cN, gmulsg(d, glog(gdiv(T, pi2), prec)));
2344 476 : D = gdiv(pi2div, D);
2345 : for(;;)
2346 7721 : {
2347 : long s0;
2348 8197 : T0 = T; T = gadd(T, D);
2349 8197 : if (gcmp(T, T2) >= 0) T = T2;
2350 8197 : s0 = gsigne(lfunhardyzeros(S, T));
2351 8197 : if (s0 != s) { s = s0; break; }
2352 7805 : if (T == T2) { setlg(w, *ct); *pw = w; return; }
2353 : }
2354 392 : z = zbrent(S, lfunhardyzeros, T0, T, prec); /* T <= T2 */
2355 392 : gerepileall(av, 2, &T, &z);
2356 392 : if (*ct > W) { W *= 2; w = vec_lengthen(w, W); }
2357 392 : if (typ(z) == t_REAL) z = rtor(z, precinit);
2358 392 : gel(w, (*ct)++) = z;
2359 : }
2360 : setlg(w, *ct); *pw = w;
2361 : }
2362 : GEN
2363 84 : lfunzeros(GEN ldata, GEN lim, long divz, long bitprec)
2364 : {
2365 84 : pari_sp ltop = avma;
2366 : GEN linit, pi2, pi2div, cN, w, T, h1, h2;
2367 84 : long i, d, NEWD, c, ct, s1, s2, prec, prec0 = nbits2prec(bitprec);
2368 : double maxt;
2369 : struct lhardyz_t S;
2370 :
2371 84 : if (is_linit(ldata) && linit_get_type(ldata) == t_LDESC_PRODUCT)
2372 : {
2373 0 : GEN M = gmael(linit_get_tech(ldata), 2,1);
2374 0 : long l = lg(M);
2375 0 : w = cgetg(l, t_VEC);
2376 0 : for (i = 1; i < l; i++) gel(w,i) = lfunzeros(gel(M,i), lim, divz, bitprec);
2377 0 : return gerepileupto(ltop, vecsort0(shallowconcat1(w), NULL, 0));
2378 : }
2379 84 : if (typ(lim) == t_VEC)
2380 : {
2381 49 : if (lg(lim) != 3 || gcmp(gel(lim,1),gel(lim,2)) >= 0)
2382 7 : pari_err_TYPE("lfunzeros",lim);
2383 42 : h1 = gel(lim,1);
2384 42 : h2 = gel(lim,2);
2385 42 : maxt = maxdd(fabs(gtodouble(h1)), fabs(gtodouble(h2)));
2386 : }
2387 : else
2388 : {
2389 35 : if (gcmp(lim,gen_0) <= 0) pari_err_TYPE("lfunzeros",lim);
2390 35 : h1 = gen_0;
2391 35 : h2 = lim;
2392 35 : maxt = gtodouble(h2);
2393 : }
2394 77 : S.linit = linit = lfuncenterinit(ldata, maxt, -1, bitprec);
2395 77 : S.bitprec = bitprec;
2396 77 : S.prec = prec0;
2397 77 : ldata = linit_get_ldata(linit);
2398 77 : d = ldata_get_degree(ldata);
2399 :
2400 77 : NEWD = minss((long) ceil(bitprec + (M_PI/(4*M_LN2)) * d * maxt),
2401 : lfun_get_bitprec(linit_get_tech(linit)));
2402 77 : prec = nbits2prec(NEWD);
2403 77 : cN = gdiv(ldata_get_conductor(ldata), gpowgs(Pi2n(-1, prec), d));
2404 77 : cN = gexpo(cN) >= 0? gaddsg(d, gmulsg(2, glog(cN, prec))): utoi(d);
2405 77 : pi2 = Pi2n(1, prec);
2406 77 : pi2div = gdivgu(pi2, labs(divz));
2407 77 : s1 = gsigne(h1);
2408 77 : s2 = gsigne(h2);
2409 77 : w = cgetg(100+1, t_VEC); c = 1; ct = 0; T = NULL;
2410 77 : if (s1 <= 0 && s2 >= 0)
2411 : {
2412 56 : GEN r = ldata_get_residue(ldata);
2413 56 : if (!r || gequal0(r))
2414 : {
2415 35 : ct = lfunorderzero(linit, -1, bitprec);
2416 35 : if (ct) T = real2n(-prec2nbits(prec) / (2*ct), prec);
2417 : }
2418 : }
2419 77 : if (s1 <= 0)
2420 : {
2421 63 : if (s1 < 0)
2422 21 : lfunzeros_i(&S, &w, &c, h1, T? negr(T): h2,
2423 : d, cN, pi2, pi2div, prec0, prec);
2424 63 : if (ct)
2425 : {
2426 21 : long n = lg(w)-1;
2427 21 : if (c + ct >= n) w = vec_lengthen(w, n + ct);
2428 84 : for (i = 1; i <= ct; i++) gel(w,c++) = gen_0;
2429 : }
2430 : }
2431 77 : if (s2 > 0 && (T || s1 >= 0))
2432 63 : lfunzeros_i(&S, &w, &c, T? T: h1, h2, d, cN, pi2, pi2div, prec0, prec);
2433 77 : return gerepilecopy(ltop, w);
2434 : }
2435 :
2436 : /*******************************************************************/
2437 : /* Guess conductor */
2438 : /*******************************************************************/
2439 : struct huntcond_t {
2440 : GEN k;
2441 : GEN theta, thetad;
2442 : GEN *pM, *psqrtM, *pMd, *psqrtMd;
2443 : };
2444 :
2445 : static void
2446 11962 : condset(struct huntcond_t *S, GEN M, long prec)
2447 : {
2448 11962 : *(S->pM) = M;
2449 11962 : *(S->psqrtM) = gsqrt(ginv(M), prec);
2450 11962 : if (S->thetad != S->theta)
2451 : {
2452 0 : *(S->pMd) = *(S->pM);
2453 0 : *(S->psqrtMd) = *(S->psqrtM);
2454 : }
2455 11962 : }
2456 :
2457 : /* M should eventually converge to N, the conductor. L has no pole. */
2458 : static GEN
2459 6888 : wrap1(void *E, GEN M)
2460 : {
2461 6888 : struct huntcond_t *S = (struct huntcond_t*)E;
2462 : GEN thetainit, tk, p1, p1inv;
2463 6888 : GEN t = mkfrac(stoi(11), stoi(10));
2464 : long prec, bitprec;
2465 :
2466 6888 : thetainit = linit_get_tech(S->theta);
2467 6888 : bitprec = theta_get_bitprec(thetainit);
2468 6888 : prec = nbits2prec(bitprec);
2469 6888 : condset(S, M, prec);
2470 6888 : tk = gpow(t, S->k, prec);
2471 6888 : p1 = lfuntheta(S->thetad, t, 0, bitprec);
2472 6888 : p1inv = lfuntheta(S->theta, ginv(t), 0, bitprec);
2473 6888 : return glog(gabs(gmul(tk, gdiv(p1, p1inv)), prec), prec);
2474 : }
2475 :
2476 : /* M should eventually converge to N, the conductor. L has a pole. */
2477 : static GEN
2478 5032 : wrap2(void *E, GEN M)
2479 : {
2480 5032 : struct huntcond_t *S = (struct huntcond_t*)E;
2481 : GEN t1k, t2k, p1, p1inv, p2, p2inv, thetainit, R;
2482 5032 : GEN t1 = mkfrac(stoi(11), stoi(10)), t2 = mkfrac(stoi(13), stoi(11));
2483 : GEN t1be, t2be, t1bemk, t2bemk, t1kmbe, t2kmbe;
2484 : GEN F11, F12, F21, F22, P1, P2, res;
2485 : long prec, bitprec;
2486 5032 : GEN k = S->k;
2487 :
2488 5032 : thetainit = linit_get_tech(S->theta);
2489 5032 : bitprec = theta_get_bitprec(thetainit);
2490 5032 : prec = nbits2prec(bitprec);
2491 5032 : condset(S, M, prec);
2492 :
2493 5032 : p1 = lfuntheta(S->thetad, t1, 0, bitprec);
2494 5032 : p2 = lfuntheta(S->thetad, t2, 0, bitprec);
2495 5032 : p1inv = lfuntheta(S->theta, ginv(t1), 0, bitprec);
2496 5032 : p2inv = lfuntheta(S->theta, ginv(t2), 0, bitprec);
2497 5032 : t1k = gpow(t1, k, prec);
2498 5032 : t2k = gpow(t2, k, prec);
2499 5032 : R = theta_get_R(thetainit);
2500 5032 : if (typ(R) == t_VEC)
2501 : {
2502 0 : GEN be = gmael(R, 1, 1);
2503 0 : t1be = gpow(t1, be, prec); t1bemk = gdiv(gsqr(t1be), t1k);
2504 0 : t2be = gpow(t2, be, prec); t2bemk = gdiv(gsqr(t2be), t2k);
2505 0 : t1kmbe = gdiv(t1k, t1be);
2506 0 : t2kmbe = gdiv(t2k, t2be);
2507 : }
2508 : else
2509 : { /* be = k */
2510 5032 : t1be = t1k; t1bemk = t1k; t1kmbe = gen_1;
2511 5032 : t2be = t2k; t2bemk = t2k; t2kmbe = gen_1;
2512 : }
2513 5032 : F11 = conj_i(gsub(gmul(gsqr(t1be), p1), p1inv));
2514 5032 : F12 = conj_i(gsub(gmul(gsqr(t2be), p2), p2inv));
2515 5032 : F21 = gsub(gmul(t1k, p1), gmul(t1bemk, p1inv));
2516 5032 : F22 = gsub(gmul(t2k, p2), gmul(t2bemk, p2inv));
2517 5032 : P1 = gsub(gmul(t1bemk, t1be), t1kmbe);
2518 5032 : P2 = gsub(gmul(t2bemk, t2be), t2kmbe);
2519 5032 : res = gdiv(gsub(gmul(P2,F21), gmul(P1,F22)),
2520 : gsub(gmul(P2,F11), gmul(P1,F12)));
2521 5032 : return glog(gabs(res, prec), prec);
2522 : }
2523 :
2524 : /* If flag = 0 (default) return all conductors found as integers. If
2525 : flag = 1, return the approximations, not the integers. If flag = 2,
2526 : return all, even nonintegers. */
2527 :
2528 : static GEN
2529 84 : checkconductor(GEN v, long bit, long flag)
2530 : {
2531 : GEN w;
2532 84 : long e, j, k, l = lg(v);
2533 84 : if (flag == 2) return v;
2534 84 : w = cgetg(l, t_VEC);
2535 322 : for (j = k = 1; j < l; j++)
2536 : {
2537 238 : GEN N = grndtoi(gel(v,j), &e);
2538 238 : if (e < -bit) gel(w,k++) = flag ? gel(v,j): N;
2539 : }
2540 84 : if (k == 2) return gel(w,1);
2541 7 : setlg(w,k); return w;
2542 : }
2543 :
2544 : static GEN
2545 98 : parse_maxcond(GEN maxN)
2546 : {
2547 : GEN M;
2548 98 : if (!maxN)
2549 49 : M = utoipos(10000);
2550 49 : else if (typ(maxN) == t_VEC)
2551 : {
2552 14 : if (!RgV_is_ZV(maxN)) pari_err_TYPE("lfunconductor",maxN);
2553 14 : return ZV_sort_shallow(maxN);
2554 : }
2555 : else
2556 35 : M = maxN;
2557 84 : return (typ(M) == t_INT)? addiu(M, 1): gceil(M);
2558 : }
2559 :
2560 : GEN
2561 98 : lfunconductor(GEN data, GEN maxcond, long flag, long bitprec)
2562 : {
2563 : struct huntcond_t S;
2564 98 : pari_sp av = avma;
2565 98 : GEN ldata = lfunmisc_to_ldata_shallow(data);
2566 98 : GEN ld, r, v, theta, thetad, M, tdom, t0 = NULL, t0i = NULL;
2567 : GEN (*eval)(void *, GEN);
2568 : long prec;
2569 98 : M = parse_maxcond(maxcond);
2570 98 : r = ldata_get_residue(ldata);
2571 98 : if (typ(M) == t_VEC) /* select in list */
2572 : {
2573 14 : if (lg(M) == 1) { set_avma(av); return cgetg(1,t_VEC); }
2574 7 : eval = NULL; tdom = dbltor(0.7);
2575 : }
2576 84 : else if (!r) { eval = wrap1; tdom = uutoQ(10,11); }
2577 : else
2578 : {
2579 21 : if (typ(r) == t_VEC && lg(r) > 2)
2580 0 : pari_err_IMPL("multiple poles in lfunconductor");
2581 21 : eval = wrap2; tdom = uutoQ(11,13);
2582 : }
2583 91 : if (eval) bitprec += bitprec/2;
2584 91 : prec = nbits2prec(bitprec);
2585 91 : ld = shallowcopy(ldata);
2586 91 : gel(ld, 5) = eval? M: veclast(M);
2587 91 : theta = lfunthetainit_i(ld, tdom, 0, bitprec);
2588 91 : thetad = theta_dual(theta, ldata_get_dual(ldata));
2589 91 : gel(theta,3) = shallowcopy(linit_get_tech(theta));
2590 91 : S.k = ldata_get_k(ldata);
2591 91 : S.theta = theta;
2592 91 : S.thetad = thetad? thetad: theta;
2593 91 : S.pM = &gel(linit_get_ldata(theta),5);
2594 91 : S.psqrtM = &gel(linit_get_tech(theta),7);
2595 91 : if (thetad)
2596 : {
2597 0 : S.pMd = &gel(linit_get_ldata(thetad),5);
2598 0 : S.psqrtMd = &gel(linit_get_tech(thetad),7);
2599 : }
2600 91 : if (!eval)
2601 : {
2602 7 : long i, besti = 0, beste = -10, l = lg(M);
2603 7 : t0 = uutoQ(11,10); t0i = uutoQ(10,11);
2604 49 : for (i = 1; i < l; i++)
2605 : {
2606 42 : pari_sp av2 = avma;
2607 : long e;
2608 42 : condset(&S, gel(M,i), prec);
2609 42 : e = lfuncheckfeq_i(theta, thetad, t0, t0i, bitprec);
2610 42 : set_avma(av2);
2611 42 : if (e < beste) { beste = e; besti = i; }
2612 35 : else if (e == beste) beste = besti = 0; /* tie: forget */
2613 : }
2614 7 : if (!besti) { set_avma(av); return cgetg(1,t_VEC); }
2615 7 : return gerepilecopy(av, mkvec2(gel(M,besti), stoi(beste)));
2616 : }
2617 84 : v = solvestep((void*)&S, eval, ghalf, M, gen_2, 14, prec);
2618 84 : return gerepilecopy(av, checkconductor(v, bitprec/2, flag));
2619 : }
2620 :
2621 : /* assume chi primitive */
2622 : static GEN
2623 560 : znchargauss_i(GEN G, GEN chi, long bitprec)
2624 : {
2625 560 : GEN z, q, F = znstar_get_N(G);
2626 : long prec;
2627 :
2628 560 : if (equali1(F)) return gen_1;
2629 560 : prec = nbits2prec(bitprec);
2630 560 : q = sqrtr_abs(itor(F, prec));
2631 560 : z = lfuntheta(mkvec2(G,chi), gen_1, 0, bitprec);
2632 560 : if (gexpo(z) < 10 - bitprec)
2633 : {
2634 28 : if (equaliu(F,300))
2635 : {
2636 14 : GEN z = rootsof1u_cx(25, prec);
2637 14 : GEN n = znconreyexp(G, chi);
2638 14 : if (equaliu(n, 131)) return gmul(q, gpowgs(z,14));
2639 7 : if (equaliu(n, 71)) return gmul(q, gpowgs(z,11));
2640 : }
2641 14 : if (equaliu(F,600))
2642 : {
2643 14 : GEN z = rootsof1u_cx(25, prec);
2644 14 : GEN n = znconreyexp(G, chi);
2645 14 : if (equaliu(n, 491)) return gmul(q, gpowgs(z,7));
2646 7 : if (equaliu(n, 11)) return gmul(q, gpowgs(z,18));
2647 : }
2648 0 : pari_err_BUG("znchargauss [ Theta(chi,1) = 0 ]");
2649 : }
2650 532 : z = gmul(gdiv(z, conj_i(z)), q);
2651 532 : if (zncharisodd(G,chi)) z = mulcxI(z);
2652 532 : return z;
2653 : }
2654 : static GEN
2655 560 : Z_radical(GEN N, long *om)
2656 : {
2657 560 : GEN P = gel(Z_factor(N), 1);
2658 560 : *om = lg(P)-1; return ZV_prod(P);
2659 : }
2660 : GEN
2661 1218 : znchargauss(GEN G, GEN chi, GEN a, long bitprec)
2662 : {
2663 : GEN v, T, N, F, b0, b1, b2, bF, a1, aF, A, r, GF, tau, B, faB, u, S;
2664 1218 : long omb0, prec = nbits2prec(bitprec);
2665 1218 : pari_sp av = avma;
2666 :
2667 1218 : if (typ(chi) != t_COL) chi = znconreylog(G,chi);
2668 1218 : T = znchartoprimitive(G, chi);
2669 1218 : GF = gel(T,1);
2670 1218 : chi = gel(T,2); /* now primitive */
2671 1218 : N = znstar_get_N(G);
2672 1218 : F = znstar_get_N(GF);
2673 1218 : if (equalii(N,F)) b1 = bF = gen_1;
2674 : else
2675 : {
2676 245 : v = Z_ppio(diviiexact(N,F), F);
2677 245 : bF = gel(v,2); /* (N/F, F^oo) */
2678 245 : b1 = gel(v,3); /* cofactor */
2679 : }
2680 1218 : if (!a) a = a1 = aF = gen_1;
2681 : else
2682 : {
2683 1169 : if (typ(a) != t_INT) pari_err_TYPE("znchargauss",a);
2684 1169 : a = modii(a, N);
2685 1169 : if (!signe(a)) { set_avma(av); return is_pm1(F)? eulerphi(N): gen_0; }
2686 728 : v = Z_ppio(a, F);
2687 728 : aF = gel(v,2);
2688 728 : a1 = gel(v,3);
2689 : }
2690 777 : if (!equalii(aF, bF)) { set_avma(av); return gen_0; }
2691 560 : b0 = Z_radical(b1, &omb0);
2692 560 : b2 = diviiexact(b1, b0);
2693 560 : A = dvmdii(a1, b2, &r);
2694 560 : if (r != gen_0) { set_avma(av); return gen_0; }
2695 560 : B = gcdii(A,b0); faB = Z_factor(B); /* squarefree */
2696 560 : S = eulerphi(mkvec2(B,faB));
2697 560 : if (odd(omb0 + lg(gel(faB,1))-1)) S = negi(S); /* moebius(b0/B) * phi(B) */
2698 560 : S = mulii(S, mulii(aF,b2));
2699 560 : tau = znchargauss_i(GF, chi, bitprec);
2700 560 : u = Fp_div(b0, A, F);
2701 560 : if (!equali1(u))
2702 : {
2703 7 : GEN ord = zncharorder(GF, chi), z = rootsof1_cx(ord, prec);
2704 7 : tau = gmul(tau, znchareval(GF, chi, u, mkvec2(z,ord)));
2705 : }
2706 560 : return gerepileupto(av, gmul(tau, S));
2707 : }
|