Code coverage tests

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

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

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

LCOV - code coverage report
Current view: top level - basemath - lfun.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.11.0 lcov report (development 22860-5579deb0b) Lines: 1341 1408 95.2 %
Date: 2018-07-19 05:36:41 Functions: 136 137 99.3 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2015  The PARI group.
       2             : 
       3             : This file is part of the PARI/GP package.
       4             : 
       5             : PARI/GP is free software; you can redistribute it and/or modify it under the
       6             : terms of the GNU General Public License as published by the Free Software
       7             : Foundation. It is distributed in the hope that it will be useful, but WITHOUT
       8             : ANY WARRANTY WHATSOEVER.
       9             : 
      10             : Check the License for details. You should have received a copy of it, along
      11             : with the package; see the file 'COPYING'. If not, write to the Free Software
      12             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      13             : 
      14             : /********************************************************************/
      15             : /**                                                                **/
      16             : /**                       L-functions                              **/
      17             : /**                                                                **/
      18             : /********************************************************************/
      19             : 
      20             : #include "pari.h"
      21             : #include "paripriv.h"
      22             : 
      23             : /*******************************************************************/
      24             : /*  Accessors                                                      */
      25             : /*******************************************************************/
      26             : 
      27             : static GEN
      28       11721 : mysercoeff(GEN x, long n)
      29             : {
      30       11721 :   long N = n - valp(x);
      31       11721 :   return (N < 0)? gen_0: gel(x, N+2);
      32             : }
      33             : 
      34             : long
      35        5523 : ldata_get_type(GEN ldata) { return mael3(ldata, 1, 1, 1); }
      36             : 
      37             : GEN
      38       13636 : ldata_get_an(GEN ldata) { return gel(ldata, 1); }
      39             : 
      40             : GEN
      41       28614 : ldata_get_dual(GEN ldata) { return gel(ldata, 2); }
      42             : 
      43             : long
      44        1706 : ldata_isreal(GEN ldata) { return isintzero(gel(ldata, 2)); }
      45             : 
      46             : GEN
      47      167968 : ldata_get_gammavec(GEN ldata) { return gel(ldata, 3); }
      48             : 
      49             : long
      50       11079 : ldata_get_degree(GEN ldata) { return lg(gel(ldata, 3))-1; }
      51             : 
      52             : long
      53       89389 : ldata_get_k(GEN ldata)
      54             : {
      55       89389 :   GEN w = gel(ldata,4);
      56       89389 :   if (typ(w) == t_VEC) w = gel(w,1);
      57       89389 :   return itos(w);
      58             : }
      59             : /* a_n = O(n^{k1 + epsilon}) */
      60             : static double
      61       53850 : ldata_get_k1(GEN ldata)
      62             : {
      63       53850 :   GEN w = gel(ldata,4);
      64             :   long k;
      65       53850 :   if (typ(w) == t_VEC) return gtodouble(gel(w,2));
      66             :   /* by default, assume that k1 = k-1 and even (k-1)/2 for entire functions */
      67       53521 :   k = itos(w);
      68       53521 :   return ldata_get_residue(ldata)? k-1: (k-1)/2.;
      69             : }
      70             : 
      71             : GEN
      72      121256 : ldata_get_conductor(GEN ldata) { return gel(ldata, 5); }
      73             : 
      74             : GEN
      75       44264 : ldata_get_rootno(GEN ldata) { return gel(ldata, 6); }
      76             : 
      77             : GEN
      78      102916 : ldata_get_residue(GEN ldata) { return lg(ldata) == 7 ? NULL: gel(ldata, 7); }
      79             : 
      80             : long
      81       72009 : linit_get_type(GEN linit) { return mael(linit, 1, 1); }
      82             : 
      83             : GEN
      84      106235 : linit_get_ldata(GEN linit) { return gel(linit, 2); }
      85             : 
      86             : GEN
      87      131075 : linit_get_tech(GEN linit) { return gel(linit, 3); }
      88             : 
      89             : long
      90      109165 : is_linit(GEN data)
      91             : {
      92      196969 :   return lg(data) == 4 && typ(data) == t_VEC
      93      196966 :                        && typ(gel(data, 1)) == t_VECSMALL;
      94             : }
      95             : 
      96             : GEN
      97       18070 : lfun_get_step(GEN tech) { return gmael(tech, 2, 1);}
      98             : 
      99             : GEN
     100       18070 : lfun_get_pol(GEN tech) { return gmael(tech, 2, 2);}
     101             : 
     102             : GEN
     103        4205 : lfun_get_Residue(GEN tech) { return gmael(tech, 2, 3);}
     104             : 
     105             : GEN
     106       28484 : lfun_get_k2(GEN tech) { return gmael(tech, 3, 1);}
     107             : 
     108             : GEN
     109       10834 : lfun_get_w2(GEN tech) { return gmael(tech, 3, 2);}
     110             : 
     111             : GEN
     112       10834 : lfun_get_expot(GEN tech) { return gmael(tech, 3, 3);}
     113             : 
     114             : GEN
     115        3955 : lfun_get_factgammavec(GEN tech) { return gmael(tech, 3, 4); }
     116             : 
     117             : static long
     118       13531 : vgaell(GEN Vga)
     119             : {
     120             :   GEN c;
     121       13531 :   long d = lg(Vga)-1;
     122       13531 :   if (d != 2) return 0;
     123        8733 :   c = gsub(gel(Vga,1), gel(Vga,2));
     124        8733 :   return gequal1(c) || gequalm1(c);
     125             : }
     126             : static long
     127       38956 : vgaeasytheta(GEN Vga) { return lg(Vga)-1 == 1 || vgaell(Vga); }
     128             : /* return b(n) := a(n) * n^c, when vgaeasytheta(Vga) is set */
     129             : static GEN
     130        9758 : antwist(GEN an, GEN Vga, long prec)
     131             : {
     132             :   long l, i;
     133        9758 :   GEN b, c = vecmin(Vga);
     134        9758 :   if (gequal0(c)) return an;
     135        1043 :   l = lg(an); b = cgetg(l, t_VEC);
     136        1043 :   if (gequal1(c))
     137             :   {
     138         672 :     if (typ(an) == t_VECSMALL)
     139         231 :       for (i = 1; i < l; i++) gel(b,i) = mulss(an[i], i);
     140             :     else
     141         441 :       for (i = 1; i < l; i++) gel(b,i) = gmulgs(gel(an,i), i);
     142             :   }
     143             :   else
     144             :   {
     145         371 :     GEN v = vecpowug(l-1, c, prec);
     146         371 :     if (typ(an) == t_VECSMALL)
     147           0 :       for (i = 1; i < l; i++) gel(b,i) = gmulsg(an[i], gel(v,i));
     148             :     else
     149         371 :       for (i = 1; i < l; i++) gel(b,i) = gmul(gel(an,i), gel(v,i));
     150             :   }
     151        1043 :   return b;
     152             : }
     153             : 
     154             : static GEN
     155        6321 : theta_dual(GEN theta, GEN bn)
     156             : {
     157        6321 :   if (typ(bn)==t_INT) return NULL;
     158             :   else
     159             :   {
     160         105 :     GEN thetad = shallowcopy(theta), ldata = linit_get_ldata(theta);
     161         105 :     GEN Vga = ldata_get_gammavec(ldata);
     162         105 :     GEN tech = shallowcopy(linit_get_tech(theta));
     163         105 :     GEN an = theta_get_an(tech);
     164         105 :     long prec = nbits2prec(theta_get_bitprec(tech));
     165         105 :     GEN vb = ldata_vecan(bn, lg(an)-1, prec);
     166         105 :     if (!theta_get_m(tech) && vgaeasytheta(Vga)) vb = antwist(vb, Vga, prec);
     167         105 :     gel(tech,1) = vb;
     168         105 :     gel(thetad,3) = tech; return thetad;
     169             :   }
     170             : }
     171             : 
     172             : static GEN
     173       31431 : domain_get_dom(GEN domain)  { return gel(domain,1); }
     174             : static long
     175       13732 : domain_get_der(GEN domain)  { return mael2(domain, 2, 1); }
     176             : static long
     177       18506 : domain_get_bitprec(GEN domain)  { return mael2(domain, 2, 2); }
     178             : GEN
     179       31844 : lfun_get_domain(GEN tech) { return gel(tech,1); }
     180             : long
     181          42 : lfun_get_bitprec(GEN tech){ return domain_get_bitprec(lfun_get_domain(tech)); }
     182             : GEN
     183           0 : lfun_get_dom(GEN tech) { return domain_get_dom(lfun_get_domain(tech)); }
     184             : 
     185             : GEN
     186        1797 : lfunprod_get_fact(GEN tech)  { return gel(tech, 2); }
     187             : 
     188             : GEN
     189       32250 : theta_get_an(GEN tdata)        { return gel(tdata, 1);}
     190             : GEN
     191        5204 : theta_get_K(GEN tdata)         { return gel(tdata, 2);}
     192             : GEN
     193        2663 : theta_get_R(GEN tdata)         { return gel(tdata, 3);}
     194             : long
     195       43479 : theta_get_bitprec(GEN tdata)   { return itos(gel(tdata, 4));}
     196             : long
     197       64217 : theta_get_m(GEN tdata)         { return itos(gel(tdata, 5));}
     198             : GEN
     199       34707 : theta_get_tdom(GEN tdata)      { return gel(tdata, 6);}
     200             : GEN
     201       36534 : theta_get_sqrtN(GEN tdata)     { return gel(tdata, 7);}
     202             : 
     203             : /*******************************************************************/
     204             : /*  Helper functions related to Gamma products                     */
     205             : /*******************************************************************/
     206             : 
     207             : /* return -itos(s) >= 0 if s is (approximately) equal to a non-positive
     208             :  * integer, and -1 otherwise */
     209             : static long
     210       15246 : isnegint(GEN s)
     211             : {
     212       15246 :   GEN r = ground(real_i(s));
     213       15246 :   if (signe(r) <= 0 && gequal(s, r)) return -itos(r);
     214       15204 :   return -1;
     215             : }
     216             : 
     217             : /* pi^(-s/2) Gamma(s/2) */
     218             : static GEN
     219       10801 : gamma_R(GEN s, long prec)
     220             : {
     221       10801 :   GEN s2 = gdivgs(s, 2), pi = mppi(prec);
     222       10801 :   long ms = isnegint(s2);
     223       10801 :   if (ms >= 0)
     224             :   {
     225          42 :     GEN pr = gmul(powru(pi, ms), gdivsg(odd(ms)? -2: 2, mpfact(ms)));
     226          42 :     GEN S = scalarser(pr, 0, 1);
     227          42 :     setvalp(S,-1); return S;
     228             :   }
     229       10759 :   return gdiv(ggamma(s2,prec), gpow(pi,s2,prec));
     230             : }
     231             : 
     232             : /* gamma_R(s)gamma_R(s+1) = 2 (2pi)^(-s) Gamma(s) */
     233             : static GEN
     234        4445 : gamma_C(GEN s, long prec)
     235             : {
     236        4445 :   GEN pi2 = Pi2n(1,prec);
     237        4445 :   long ms = isnegint(s);
     238        4445 :   if (ms >= 0)
     239             :   {
     240           0 :     GEN pr = gmul(powrs(pi2, ms), gdivsg(odd(ms)? -2: 2, mpfact(ms)));
     241           0 :     GEN S = scalarser(pr, 0, 1);
     242           0 :     setvalp(S,-1); return S;
     243             :   }
     244        4445 :   return gmul2n(gdiv(ggamma(s,prec), gpow(pi2,s,prec)), 1);
     245             : }
     246             : 
     247             : static GEN
     248        1099 : gammafrac(GEN r, long d)
     249             : {
     250        1099 :   GEN pr, a = gmul2n(r, -1);
     251        1099 :   GEN polj = cgetg(labs(d)+1, t_COL);
     252        1099 :   long i, v=0;
     253        1099 :   if (d > 0)
     254           0 :     for (i = 1; i <= d; ++i)
     255           0 :       gel(polj, i) = deg1pol_shallow(ghalf, gaddgs(a, i-1), v);
     256             :   else
     257        2198 :     for (i = 1; i <= -d; ++i)
     258        1099 :       gel(polj, i) = deg1pol_shallow(ghalf, gsubgs(a, i), v);
     259        1099 :   pr = RgV_prod(polj);
     260        1099 :   return d < 0 ? ginv(pr): pr;
     261             : }
     262             : 
     263             : static GEN
     264       15477 : gammafactor(GEN Vga)
     265             : {
     266       15477 :   pari_sp av = avma;
     267       15477 :   long i, m, d = lg(Vga)-1, dr, dc;
     268       15477 :   GEN pol = pol_1(0), pi = gen_0, R = cgetg(d+1,t_VEC);
     269             :   GEN P, F, FR, FC, E, ER, EC;
     270       39620 :   for (i = 1; i <= d; ++i)
     271             :   {
     272       24143 :     GEN a = gel(Vga,i), qr = gdiventres(real_i(a), gen_2);
     273       24143 :     long q = itos(gel(qr,1));
     274       24143 :     gel(R, i) = gadd(gel(qr,2), imag_i(a));
     275       24143 :     if (q)
     276             :     {
     277        1099 :       pol = gmul(pol, gammafrac(gel(R,i), q));
     278        1099 :       pi  = addis(pi, q);
     279             :     }
     280             :   }
     281       15477 :   gen_sort_inplace(R, (void*)cmp_universal, cmp_nodata, &P);
     282       15477 :   F = cgetg(d+1, t_VEC); E = cgetg(d+1, t_VECSMALL);
     283       51709 :   for (i = 1, m = 0; i <= d;)
     284             :   {
     285             :     long k;
     286       20755 :     GEN u = gel(R, i);
     287       24143 :     for(k = i + 1; k <= d; ++k)
     288        8666 :       if (cmp_universal(gel(R, k), u)) break;
     289       20755 :     m++;
     290       20755 :     E[m] = k - i;
     291       20755 :     gel(F, m) = u;
     292       20755 :     i = k;
     293             :   }
     294       15477 :   setlg(F, m+1); setlg(E, m+1);
     295       15477 :   R = cgetg(m+1, t_VEC);
     296       36232 :   for (i = 1; i <= m; i++)
     297             :   {
     298       20755 :     GEN qr = gdiventres(gel(F,i), gen_1);
     299       20755 :     gel(R, i) = mkvec2(gel(qr,2), stoi(E[i]));
     300             :   }
     301       15477 :   gen_sort_inplace(R, (void*)cmp_universal, cmp_nodata, &P);
     302       15477 :   FR = cgetg(m+1, t_VEC); ER = cgetg(m+1, t_VECSMALL);
     303       15477 :   FC = cgetg(m+1, t_VEC); EC = cgetg(m+1, t_VECSMALL);
     304       47152 :   for (i = 1, dr = 1, dc = 1; i <= m;)
     305             :   {
     306       16198 :     if (i==m || cmp_universal(gel(R,i), gel(R,i+1)))
     307             :     {
     308       11641 :       gel(FR, dr) = gel(F, P[i]);
     309       11641 :       ER[dr] = E[P[i]];
     310       11641 :       dr++; i++;
     311             :     } else
     312             :     {
     313        4557 :       if (gequal(gaddgs(gmael(R,i,1), 1), gmael(R,i+1,1)))
     314           0 :         gel(FC, dc) = gel(F, P[i+1]);
     315             :       else
     316        4557 :         gel(FC, dc) = gel(F, P[i]);
     317        4557 :       EC[dc] = E[P[i]];
     318        4557 :       dc++; i+=2;
     319             :     }
     320             :   }
     321       15477 :   setlg(FR, dr); setlg(ER, dr);
     322       15477 :   setlg(FC, dc); setlg(EC, dc);
     323       15477 :   return gerepilecopy(av, mkvec4(pol, pi, mkvec2(FR,ER), mkvec2(FC,EC)));
     324             : }
     325             : 
     326             : static GEN
     327       13489 : deg1ser_shallow(GEN a1, GEN a0, long v, long e)
     328             : {
     329       13489 :   return RgX_to_ser(deg1pol_shallow(a1, a0, v), e+2);
     330             : }
     331             : /*
     332             : To test:
     333             : GR(s)=Pi^-(s/2)*gamma(s/2);
     334             : GC(s)=2*(2*Pi)^-s*gamma(s)
     335             : gam_direct(F,s)=prod(i=1,#F,GR(s+F[i]))
     336             : gam_fact(F,s)=my([P,p,R,C]=gammafactor(F));subst(P,x,s)*Pi^-p*prod(i=1,#R[1],GR(s+R[1][i])^R[2][i])*prod(i=1,#C[1],GC(s+C[1][i])^C[2][i])
     337             : */
     338             : 
     339             : static GEN
     340       15673 : polgammaeval(GEN F, GEN s)
     341             : {
     342       15673 :   GEN r = poleval(F, s);
     343       15673 :   if (typ(s)!=t_SER && gequal0(r))
     344             :   {
     345           0 :     long e = gvaluation(F, deg1pol(gen_1, gneg(s), varn(F)));
     346           0 :     r = poleval(F, deg1ser_shallow(gen_1, s, 0, e+1));
     347             :   }
     348       15673 :   return r;
     349             : }
     350             : 
     351             : static GEN
     352       14483 : fracgammaeval(GEN F, GEN s)
     353             : {
     354       14483 :   if (typ(F)==t_POL)
     355       13293 :     return polgammaeval(F, s);
     356        1190 :   else if (typ(F)==t_RFRAC)
     357        1190 :     return gdiv(polgammaeval(gel(F,1), s), polgammaeval(gel(F,2), s));
     358           0 :   return F;
     359             : }
     360             : 
     361             : static GEN
     362       14483 : gammafactproduct(GEN F, GEN s, long prec)
     363             : {
     364       14483 :   pari_sp av = avma;
     365       14483 :   GEN P = fracgammaeval(gel(F,1), s);
     366       14483 :   GEN p = gpow(mppi(prec),gneg(gel(F,2)), prec), z = gmul(P, p);
     367       14483 :   GEN R = gel(F,3), Rw = gel(R,1), Re=gel(R,2);
     368       14483 :   GEN C = gel(F,4), Cw = gel(C,1), Ce=gel(C,2);
     369       14483 :   long i, lR = lg(Rw), lC = lg(Cw);
     370       25284 :   for (i=1; i< lR; i++)
     371       10801 :     z = gmul(z, gpowgs(gamma_R(gadd(s,gel(Rw, i)), prec), Re[i]));
     372       18928 :   for (i=1; i< lC; i++)
     373        4445 :     z = gmul(z, gpowgs(gamma_C(gadd(s,gel(Cw, i)), prec), Ce[i]));
     374       14483 :   return gerepileupto(av, z);
     375             : }
     376             : 
     377             : static int
     378        4025 : gammaordinary(GEN Vga, GEN s)
     379             : {
     380        4025 :   long i, d = lg(Vga)-1;
     381       10738 :   for (i = 1; i <= d; i++)
     382             :   {
     383        6797 :     GEN z = gadd(s, gel(Vga,i));
     384             :     long e;
     385        6797 :     if (gsigne(z) <= 0) { (void)grndtoi(z, &e); if (e < -10) return 0; }
     386             :   }
     387        3941 :   return 1;
     388             : }
     389             : 
     390             : /* Exponent A of t in asymptotic expansion; K(t) ~ C t^A exp(-pi d t^(2/d)).
     391             :  * suma = vecsum(Vga)*/
     392             : static double
     393       53843 : gammavec_expo(long d, double suma) { return (1 - d + suma) / d; }
     394             : 
     395             : /*******************************************************************/
     396             : /*       First part: computations only involving Theta(t)          */
     397             : /*******************************************************************/
     398             : 
     399             : static void
     400       84919 : get_cone(GEN t, double *r, double *a)
     401             : {
     402       84919 :   const long prec = LOWDEFAULTPREC;
     403       84919 :   if (typ(t) == t_COMPLEX)
     404             :   {
     405       15162 :     t  = gprec_w(t, prec);
     406       15162 :     *r = gtodouble(gabs(t, prec));
     407       15162 :     *a = fabs(gtodouble(garg(t, prec)));
     408             :   }
     409             :   else
     410             :   {
     411       69757 :     *r = fabs(gtodouble(t));
     412       69757 :     *a = 0.;
     413             :   }
     414       84919 :   if (!*r && !*a) pari_err_DOMAIN("lfunthetainit","t","=",gen_0,t);
     415       84912 : }
     416             : /* slightly larger cone than necessary, to avoid round-off problems */
     417             : static void
     418       50212 : get_cone_fuzz(GEN t, double *r, double *a)
     419       50212 : { get_cone(t, r, a); *r -= 1e-10; if (*a) *a += 1e-10; }
     420             : 
     421             : /* Initialization m-th Theta derivative. tdom is either
     422             :  * - [rho,alpha]: assume |t| >= rho and |arg(t)| <= alpha
     423             :  * - a positive real scalar: assume t real, t >= tdom;
     424             :  * - a complex number t: compute at t;
     425             :  * N is the conductor (either the true one from ldata or a guess from
     426             :  * lfunconductor) */
     427             : long
     428       39012 : lfunthetacost(GEN ldata, GEN tdom, long m, long bitprec)
     429             : {
     430       39012 :   pari_sp av = avma;
     431       39012 :   GEN Vga = ldata_get_gammavec(ldata);
     432       39012 :   long d = lg(Vga)-1;
     433       39012 :   long k1 = ldata_get_k1(ldata);
     434       39012 :   double c = d/2., a, A, B, logC, al, rho, T;
     435       39012 :   double N = gtodouble(ldata_get_conductor(ldata));
     436             : 
     437       39012 :   if (!N) pari_err_TYPE("lfunthetaneed [missing conductor]", ldata);
     438       39012 :   if (typ(tdom) == t_VEC && lg(tdom) == 3)
     439             :   {
     440           7 :     rho= gtodouble(gel(tdom,1));
     441           7 :     al = gtodouble(gel(tdom,2));
     442             :   }
     443             :   else
     444       39005 :     get_cone_fuzz(tdom, &rho, &al);
     445       39005 :   A = gammavec_expo(d, gtodouble(vecsum(Vga))); avma = av;
     446       39005 :   a = (A+k1+1) + (m-1)/c;
     447       39005 :   if (fabs(a) < 1e-10) a = 0.;
     448       39005 :   logC = c*M_LN2 - log(c)/2;
     449             :   /* +1: fudge factor */
     450       39005 :   B = M_LN2*bitprec+logC+m*log(2*M_PI) + 1 + (k1+1)*log(N)/2 - (k1+m+1)*log(rho);
     451       39005 :   if (al)
     452             :   { /* t = rho e^(i*al), T^(1/c) = Re(t^(1/c)) > 0, T = rho cos^c(al/c) */
     453        7581 :     double z = cos(al/c);
     454        7581 :     T = (d == 2 && typ(tdom) != t_VEC)? gtodouble(real_i(tdom)): rho*pow(z,c);
     455        7581 :     if (z <= 0)
     456           0 :       pari_err_DOMAIN("lfunthetaneed", "arg t", ">", dbltor(c*M_PI/2), tdom);
     457        7581 :     B -= log(z) * (c * (k1+A+1) + m);
     458             :   }
     459             :   else
     460       31424 :     T = rho;
     461       39005 :   return B <= 0? 0: floor(0.9 + dblcoro526(a,c,B) / T * sqrt(N));
     462             : }
     463             : long
     464          14 : lfunthetacost0(GEN L, GEN tdom, long m, long bitprec)
     465             : {
     466             :   long n;
     467          14 :   if (is_linit(L) && linit_get_type(L)==t_LDESC_THETA)
     468           7 :   {
     469           7 :     GEN tech = linit_get_tech(L);
     470           7 :     n = lg(theta_get_an(tech))-1;
     471             :   }
     472             :   else
     473             :   {
     474           7 :     pari_sp av = avma;
     475           7 :     GEN ldata = lfunmisc_to_ldata_shallow(L);
     476           7 :     n = lfunthetacost(ldata, tdom? tdom: gen_1, m, bitprec);
     477           7 :     avma = av;
     478             :   }
     479          14 :   return n;
     480             : }
     481             : 
     482             : static long
     483        7854 : fracgammadegree(GEN FVga)
     484        7854 : { GEN F = gel(FVga,1); return (typ(F)==t_RFRAC)? degpol(gel(F,2)): 0; }
     485             : 
     486             : /* Poles of a L-function can be represented in the following ways:
     487             :  * 1) Nothing (ldata has only 6 components, ldata_get_residue = NULL).
     488             :  * 2) a complex number (single pole at s = k with given residue, unknown if 0).
     489             :  * 3) A vector (possibly empty) of 2-component vectors [a, ra], where a is the
     490             :  * pole, ra a t_SER: its Taylor expansion at a. A t_VEC encodes the polar
     491             :  * part of L, a t_COL, the polar part of Lambda */
     492             : 
     493             : /* 'a' a complex number (pole), 'r' the polar part of L at 'a';
     494             :  * return 'R' the polar part of Lambda at 'a' */
     495             : static GEN
     496        6531 : rtoR(GEN a, GEN r, GEN FVga, GEN N, long prec)
     497             : {
     498        6531 :   long v = lg(r)-2;
     499        6531 :   GEN as = deg1ser_shallow(gen_1, a, varn(r), v);
     500        6531 :   GEN Na = gpow(N, gdivgs(as, 2), prec);
     501        6531 :   long d = fracgammadegree(FVga);
     502        6531 :   if (d) as = sertoser(as, v+d); /* make up for a possible loss of accuracy */
     503        6531 :   return gmul(gmul(r, Na), gammafactproduct(FVga, as, prec));
     504             : }
     505             : 
     506             : /* assume r in normalized form: t_VEC of pairs [be,re] */
     507             : GEN
     508        6482 : lfunrtopoles(GEN r)
     509             : {
     510        6482 :   long j, l = lg(r);
     511        6482 :   GEN v = cgetg(l, t_VEC);
     512       13118 :   for (j = 1; j < l; j++)
     513             :   {
     514        6636 :     GEN rj = gel(r,j), a = gel(rj,1);
     515        6636 :     gel(v,j) = a;
     516             :   }
     517        6482 :   gen_sort_inplace(v, (void*)&cmp_universal, cmp_nodata, NULL);
     518        6482 :   return v;
     519             : }
     520             : 
     521             : /* r / x + O(1) */
     522             : static GEN
     523        5334 : simple_pole(GEN r)
     524             : {
     525             :   GEN S;
     526        5334 :   if (isintzero(r)) return gen_0;
     527        5327 :   S = deg1ser_shallow(gen_0, r, 0, 1);
     528        5327 :   setvalp(S, -1); return S;
     529             : }
     530             : static GEN
     531        5712 : normalize_simple_pole(GEN r, GEN k)
     532             : {
     533        5712 :   long tx = typ(r);
     534        5712 :   if (is_vec_t(tx)) return r;
     535        5334 :   if (!is_scalar_t(tx)) pari_err_TYPE("lfunrootres [poles]", r);
     536        5334 :   return mkvec(mkvec2(k, simple_pole(r)));
     537             : }
     538             : /* normalize the description of a polar part */
     539             : static GEN
     540        7273 : normalizepoles(GEN r, long k)
     541             : {
     542             :   long iv, j, l;
     543             :   GEN v;
     544        7273 :   if (!is_vec_t(typ(r))) return normalize_simple_pole(r, stoi(k));
     545        3003 :   v = cgetg_copy(r, &l);
     546        7056 :   for (j = iv = 1; j < l; j++)
     547             :   {
     548        4053 :     GEN rj = gel(r,j), a = gel(rj,1), ra = gel(rj,2);
     549        4053 :     if (!is_scalar_t(typ(a)) || typ(ra) != t_SER)
     550           0 :       pari_err_TYPE("lfunrootres [poles]",r);
     551        4053 :     if (valp(ra) >= 0) continue;
     552        4053 :     gel(v,iv++) = rj;
     553             :   }
     554        3003 :   setlg(v, iv); return v;
     555             : }
     556             : static int
     557        8883 : residues_known(GEN r)
     558             : {
     559        8883 :   long i, l = lg(r);
     560        8883 :   if (isintzero(r)) return 0;
     561        8757 :   if (!is_vec_t(typ(r))) return 1;
     562       10066 :   for (i = 1; i < l; i++)
     563             :   {
     564        5733 :     GEN ri = gel(r,i);
     565        5733 :     if (!is_vec_t(typ(ri)) || lg(ri)!=3)
     566           0 :       pari_err_TYPE("lfunrootres [poles]",r);
     567        5733 :     if (isintzero(gel(ri, 2))) return 0;
     568             :   }
     569        4333 :   return 1;
     570             : }
     571             : 
     572             : /* Compute R's from r's (r = Taylor devts of L(s), R of Lambda(s)).
     573             :  * 'r/eno' passed to override the one from ldata  */
     574             : static GEN
     575       16373 : lfunrtoR_i(GEN ldata, GEN r, GEN eno, long prec)
     576             : {
     577       16373 :   GEN Vga = ldata_get_gammavec(ldata), N = ldata_get_conductor(ldata);
     578             :   GEN R, vr, FVga;
     579       16373 :   pari_sp av = avma;
     580       16373 :   long lr, j, jR, k = ldata_get_k(ldata);
     581             : 
     582       16373 :   if (!r || isintzero(eno) || !residues_known(r))
     583        9100 :     return gen_0;
     584        7273 :   r = normalizepoles(r, k);
     585        7273 :   if (typ(r) == t_COL) return gerepilecopy(av, r);
     586        6377 :   if (typ(ldata_get_dual(ldata)) != t_INT)
     587           0 :     pari_err(e_MISC,"please give the Taylor development of Lambda");
     588        6377 :   vr = lfunrtopoles(r); lr = lg(vr);
     589        6377 :   FVga = gammafactor(Vga);
     590        6377 :   R = cgetg(2*lr, t_VEC);
     591       12908 :   for (j = jR = 1; j < lr; j++)
     592             :   {
     593        6531 :     GEN rj = gel(r,j), a = gel(rj,1), ra = gel(rj,2);
     594        6531 :     GEN Ra = rtoR(a, ra, FVga, N, prec);
     595        6531 :     GEN b = gsubsg(k, conj_i(a));
     596        6531 :     if (lg(Ra)-2 < -valp(Ra))
     597           0 :       pari_err(e_MISC,
     598             :         "please give more terms in L function's Taylor development at %Ps", a);
     599        6531 :     gel(R,jR++) = mkvec2(a, Ra);
     600        6531 :     if (!tablesearch(vr, b, (int (*)(GEN,GEN))&cmp_universal))
     601             :     {
     602        6454 :       GEN mX = gneg(pol_x(varn(Ra)));
     603        6454 :       GEN Rb = gmul(eno, gsubst(conj_i(Ra), varn(Ra), mX));
     604        6454 :       gel(R,jR++) = mkvec2(b, Rb);
     605             :     }
     606             :   }
     607        6377 :   setlg(R, jR); return gerepilecopy(av, R);
     608             : }
     609             : static GEN
     610       16296 : lfunrtoR_eno(GEN ldata, GEN eno, long prec)
     611       16296 : { return lfunrtoR_i(ldata, ldata_get_residue(ldata), eno, prec); }
     612             : static GEN
     613       13636 : lfunrtoR(GEN ldata, long prec)
     614       13636 : { return lfunrtoR_eno(ldata, ldata_get_rootno(ldata), prec); }
     615             : 
     616             : /* thetainit using {an: n <= L}; if (m = 0 && easytheta), an2 is an * n^al */
     617             : static GEN
     618       11214 : lfunthetainit0(GEN ldata, GEN tdom, GEN an2, long m,
     619             :     long bitprec, long extrabit)
     620             : {
     621       11214 :   long prec = nbits2prec(bitprec);
     622       11214 :   GEN tech, N = ldata_get_conductor(ldata);
     623       11214 :   GEN Vga = ldata_get_gammavec(ldata);
     624       11214 :   GEN K = gammamellininvinit(Vga, m, bitprec + extrabit);
     625       11214 :   GEN R = lfunrtoR(ldata, prec);
     626       11214 :   if (!tdom) tdom = gen_1;
     627       11214 :   if (typ(tdom) != t_VEC)
     628             :   {
     629             :     double r, a;
     630       11207 :     get_cone_fuzz(tdom, &r, &a);
     631       11207 :     tdom = mkvec2(dbltor(r), a? dbltor(a): gen_0);
     632             :   }
     633       11214 :   tech = mkvecn(7, an2,K,R, stoi(bitprec), stoi(m), tdom, gsqrt(N,prec));
     634       11214 :   return mkvec3(mkvecsmall(t_LDESC_THETA), ldata, tech);
     635             : }
     636             : 
     637             : /* tdom: 1) positive real number r, t real, t >= r; or
     638             :  *       2) [r,a], describing the cone |t| >= r, |arg(t)| <= a */
     639             : static GEN
     640        6804 : lfunthetainit_i(GEN data, GEN tdom, long m, long bitprec)
     641             : {
     642        6804 :   GEN ldata = lfunmisc_to_ldata_shallow(data);
     643        6804 :   long L = lfunthetacost(ldata, tdom, m, bitprec), prec = nbits2prec(bitprec);
     644        6797 :   GEN an = ldata_vecan(ldata_get_an(ldata), L, prec);
     645        6797 :   GEN Vga = ldata_get_gammavec(ldata);
     646        6797 :   if (m == 0 && vgaeasytheta(Vga)) an = antwist(an, Vga, prec);
     647        6797 :   return lfunthetainit0(ldata, tdom, an, m, bitprec, 32);
     648             : }
     649             : 
     650             : GEN
     651         259 : lfunthetainit(GEN ldata, GEN tdom, long m, long bitprec)
     652             : {
     653         259 :   pari_sp av = avma;
     654         259 :   GEN S = lfunthetainit_i(ldata, tdom? tdom: gen_1, m, bitprec);
     655         259 :   return gerepilecopy(av, S);
     656             : }
     657             : 
     658             : GEN
     659         819 : lfunan(GEN ldata, long L, long prec)
     660             : {
     661         819 :   pari_sp av = avma;
     662             :   GEN an ;
     663         819 :   ldata = lfunmisc_to_ldata_shallow(ldata);
     664         819 :   an = gerepilecopy(av, ldata_vecan(ldata_get_an(ldata), L, prec));
     665         812 :   if (typ(an) != t_VEC) an = vecsmall_to_vec_inplace(an);
     666         812 :   return an;
     667             : }
     668             : 
     669             : /* [1^B,...,N^B] */
     670             : GEN
     671         112 : vecpowuu(long N, ulong B)
     672             : {
     673             :   GEN v;
     674             :   long p, i;
     675             :   forprime_t T;
     676             : 
     677         112 :   if (B <= 2)
     678             :   {
     679           0 :     if (!B) return const_vec(N,gen_1);
     680           0 :     v = cgetg(N+1, t_VEC); if (N == 0) return v;
     681           0 :     gel(v,1) = gen_1;
     682           0 :     if (B == 1)
     683           0 :       for (i = 2; i <= N; i++) gel(v,i) = utoipos(i);
     684             :     else
     685           0 :       for (i = 2; i <= N; i++) gel(v,i) = sqru(i);
     686           0 :     return v;
     687             :   }
     688         112 :   v = const_vec(N, NULL);
     689         112 :   u_forprime_init(&T, 3, N);
     690        2023 :   while ((p = u_forprime_next(&T)))
     691             :   {
     692             :     long m, pk, oldpk;
     693        1799 :     gel(v,p) = powuu(p, B);
     694        3997 :     for (pk = p, oldpk = p; pk; oldpk = pk, pk = umuluu_le(pk,p,N))
     695             :     {
     696        2198 :       if (pk != p) gel(v,pk) = mulii(gel(v,oldpk), gel(v,p));
     697        9107 :       for (m = N/pk; m > 1; m--)
     698        6909 :         if (gel(v,m) && m%p) gel(v, m*pk) = mulii(gel(v,m), gel(v,pk));
     699             :     }
     700             :   }
     701         112 :   gel(v,1) = gen_1;
     702        3479 :   for (i = 2; i <= N; i+=2)
     703             :   {
     704        3367 :     long vi = vals(i);
     705        3367 :     gel(v,i) = shifti(gel(v,i >> vi), B * vi);
     706             :   }
     707         112 :   return v;
     708             : }
     709             : /* [1^B,...,N^B] */
     710             : GEN
     711        8342 : vecpowug(long N, GEN B, long prec)
     712             : {
     713        8342 :   GEN v = const_vec(N, NULL);
     714        8342 :   long p, eB = gexpo(B);
     715        8342 :   long prec0 = eB < 5? prec: prec + nbits2extraprec(eB);
     716             :   forprime_t T;
     717        8342 :   u_forprime_init(&T, 2, N);
     718        8342 :   gel(v,1) = gen_1;
     719      339242 :   while ((p = u_forprime_next(&T)))
     720             :   {
     721             :     long m, pk, oldpk;
     722      322558 :     gel(v,p) = gpow(utor(p,prec0), B, prec);
     723      322558 :     if (prec0 != prec) gel(v,p) = gprec_wtrunc(gel(v,p), prec);
     724      714585 :     for (pk = p, oldpk = p; pk; oldpk = pk, pk = umuluu_le(pk,p,N))
     725             :     {
     726      392027 :       if (pk != p) gel(v,pk) = gmul(gel(v,oldpk), gel(v,p));
     727     5968577 :       for (m = N/pk; m > 1; m--)
     728     5576550 :         if (gel(v,m) && m%p) gel(v, m*pk) = gmul(gel(v,m), gel(v,pk));
     729             :     }
     730             :   }
     731        8342 :   return v;
     732             : }
     733             : 
     734             : /* return [1^(2/d), 2^(2/d),...,lim^(2/d)] */
     735             : static GEN
     736        5204 : mkvroots(long d, long lim, long prec)
     737             : {
     738        5204 :   if (d <= 4)
     739             :   {
     740        5050 :     GEN v = cgetg(lim+1,t_VEC);
     741             :     long n;
     742        5050 :     switch(d)
     743             :     {
     744             :       case 1:
     745        1855 :         for (n=1; n <= lim; n++) gel(v,n) = sqru(n);
     746        1855 :         return v;
     747             :       case 2:
     748         924 :         for (n=1; n <= lim; n++) gel(v,n) = utoipos(n);
     749         924 :         return v;
     750             :       case 4:
     751        1370 :         for (n=1; n <= lim; n++) gel(v,n) = sqrtr(utor(n, prec));
     752        1370 :         return v;
     753             :     }
     754             :   }
     755        1055 :   return vecpowug(lim, gdivgs(gen_2,d), prec);
     756             : }
     757             : 
     758             : GEN
     759       38417 : lfunthetacheckinit(GEN data, GEN t, long m, long bitprec)
     760             : {
     761       38417 :   if (is_linit(data) && linit_get_type(data)==t_LDESC_THETA)
     762             :   {
     763       34707 :     GEN tdom, thetainit = linit_get_tech(data);
     764       34707 :     long bitprecnew = theta_get_bitprec(thetainit);
     765       34707 :     long m0 = theta_get_m(thetainit);
     766             :     double r, al, rt, alt;
     767       34707 :     if (m0 != m)
     768           0 :       pari_err_DOMAIN("lfuntheta","derivative order","!=", stoi(m),stoi(m0));
     769       34707 :     if (bitprec > bitprecnew) goto INIT;
     770       34707 :     get_cone(t, &rt, &alt);
     771       34707 :     tdom = theta_get_tdom(thetainit);
     772       34707 :     r = rtodbl(gel(tdom,1));
     773       34707 :     al= rtodbl(gel(tdom,2)); if (rt >= r && alt <= al) return data;
     774             :   }
     775             : INIT:
     776        6461 :   return lfunthetainit_i(data, t, m, bitprec);
     777             : }
     778             : 
     779             : static GEN
     780     4667911 : get_an(GEN an, long n)
     781             : {
     782     4667911 :   if (typ(an) == t_VECSMALL) { long a = an[n]; if (a) return stoi(a); }
     783     4667911 :   else { GEN a = gel(an,n); if (a && !gequal0(a)) return a; }
     784     3385790 :   return NULL;
     785             : }
     786             : /* x * an[n] */
     787             : static GEN
     788    11278649 : mul_an(GEN an, long n, GEN x)
     789             : {
     790    11278649 :   if (typ(an) == t_VECSMALL) { long a = an[n]; if (a) return gmulsg(a,x); }
     791     8476343 :   else { GEN a = gel(an,n); if (a && !gequal0(a)) return gmul(a,x); }
     792     5688089 :   return NULL;
     793             : }
     794             : /* 2*t^a * x **/
     795             : static GEN
     796      139858 : mulT(GEN t, GEN a, GEN x, long prec)
     797             : {
     798      139858 :   if (gequal0(a)) return gmul2n(x,1);
     799       10566 :   return gmul(x, gmul2n(gequal1(a)? t: gpow(t,a,prec), 1));
     800             : }
     801             : 
     802             : static GEN
     803    23317882 : vecan_cmul(void *E, GEN P, long a, GEN x)
     804             : {
     805             :   (void)E;
     806    23317882 :   return (a==0 || !gel(P,a))? NULL: gmul(gel(P,a), x);
     807             : }
     808             : /* d=2, 2 sum_{n <= limt} a(n) (n t)^al q^n, q = exp(-2pi t),
     809             :  * an2[n] = a(n) * n^al */
     810             : static GEN
     811      114510 : theta2(GEN an2, long limt, GEN t, GEN al, long prec)
     812             : {
     813      114510 :   GEN S, q, pi2 = Pi2n(1,prec);
     814      114510 :   const struct bb_algebra *alg = get_Rg_algebra();
     815      114510 :   setsigne(pi2,-1); q = gexp(gmul(pi2, t), prec);
     816             :   /* Brent-Kung in case the a_n are small integers */
     817      114510 :   S = gen_bkeval(an2, limt, q, 1, NULL, alg, vecan_cmul);
     818      114510 :   return mulT(t, al, S, prec);
     819             : }
     820             : 
     821             : /* d=1, 2 sum_{n <= limt} a_n (n t)^al q^(n^2), q = exp(-pi t^2),
     822             :  * an2[n] is a_n n^al */
     823             : static GEN
     824       25348 : theta1(GEN an2, long limt, GEN t, GEN al, long prec)
     825             : {
     826       25348 :   GEN q = gexp(gmul(negr(mppi(prec)), gsqr(t)), prec);
     827       25348 :   GEN vexp = gsqrpowers(q, limt), S = gen_0;
     828       25348 :   pari_sp av = avma;
     829             :   long n;
     830     4097007 :   for (n = 1; n <= limt; n++)
     831             :   {
     832     4071659 :     GEN c = mul_an(an2, n, gel(vexp,n));
     833     4071659 :     if (c)
     834             :     {
     835     3054803 :       S = gadd(S, c);
     836     3054803 :       if (gc_needed(av, 3)) S = gerepileupto(av, S);
     837             :     }
     838             :   }
     839       25348 :   return mulT(t, al, S, prec);
     840             : }
     841             : 
     842             : /* If m > 0, compute m-th derivative of theta(t) = theta0(t/sqrt(N))
     843             :  * with absolute error 2^-bitprec; theta(t)=\sum_{n\ge1}a(n)K(nt/N^(1/2)) */
     844             : GEN
     845       32124 : lfuntheta(GEN data, GEN t, long m, long bitprec)
     846             : {
     847       32124 :   pari_sp ltop = avma;
     848             :   long limt, d;
     849             :   GEN sqN, vecan, Vga, ldata, theta, thetainit, S;
     850       32124 :   long n, prec = nbits2prec(bitprec);
     851       32124 :   t = gprec_w(t, prec);
     852       32124 :   theta = lfunthetacheckinit(data, t, m, bitprec);
     853       32117 :   ldata = linit_get_ldata(theta);
     854       32117 :   thetainit = linit_get_tech(theta);
     855       32117 :   vecan = theta_get_an(thetainit);
     856       32117 :   sqN = theta_get_sqrtN(thetainit);
     857       32117 :   limt = lg(vecan)-1;
     858       32117 :   if (theta == data)
     859       30878 :     limt = minss(limt, lfunthetacost(ldata, t, m, bitprec));
     860       32117 :   if (!limt)
     861             :   {
     862          14 :     avma = ltop; S = real_0_bit(-bitprec);
     863          14 :     if (!is_real_t(typ(t)) || !ldata_isreal(ldata))
     864           7 :       S = gerepilecopy(ltop, mkcomplex(S,S));
     865          14 :     return S;
     866             :   }
     867       32103 :   t = gdiv(t, sqN);
     868       32103 :   Vga = ldata_get_gammavec(ldata);
     869       32103 :   d = lg(Vga)-1;
     870       32103 :   if (m == 0 && vgaeasytheta(Vga))
     871             :   {
     872       29405 :     if (theta_get_m(thetainit) > 0) vecan = antwist(vecan, Vga, prec);
     873       29405 :     if (d == 1) S = theta1(vecan, limt, t, gel(Vga,1), prec);
     874        4057 :     else        S = theta2(vecan, limt, t, vecmin(Vga), prec);
     875       29405 :     return gerepileupto(ltop, S);
     876             :   }
     877             :   else
     878             :   {
     879        2698 :     GEN K = theta_get_K(thetainit);
     880        2698 :     GEN vroots = mkvroots(d, limt, prec);
     881             :     pari_sp av;
     882        2698 :     t = gpow(t, gdivgs(gen_2,d), prec);
     883        2698 :     S = gen_0; av = avma;
     884     4670609 :     for (n = 1; n <= limt; ++n)
     885             :     {
     886     4667911 :       GEN nt, an = get_an(vecan, n);
     887     4667911 :       if (!an) continue;
     888     1282121 :       nt = gmul(gel(vroots,n), t);
     889     1282121 :       if (m) an = gmul(an, powuu(n, m));
     890     1282121 :       S = gadd(S, gmul(an, gammamellininvrt(K, nt, bitprec)));
     891     1282121 :       if ((n & 0x1ff) == 0) S = gerepileupto(av, S);
     892             :     }
     893        2698 :     if (m) S = gdiv(S, gpowgs(sqN, m));
     894        2698 :     return gerepileupto(ltop, S);
     895             :   }
     896             : }
     897             : 
     898             : /*******************************************************************/
     899             : /* Second part: Computation of L-Functions.                        */
     900             : /*******************************************************************/
     901             : 
     902             : struct lfunp {
     903             :   long precmax, Dmax, D, M, m0, nmax, d;
     904             :   double k1, E, logN2, logC, A, hd, dc, dw, dh, MAXs, sub;
     905             :   GEN L, vprec, an, bn;
     906             : };
     907             : 
     908             : static void
     909       14838 : lfunparams(GEN ldata, long der, long bitprec, struct lfunp *S)
     910             : {
     911       14838 :   const long derprec = (der > 1)? dbllog2(mpfact(der)): 0; /* log2(der!) */
     912             :   GEN Vga, N, L;
     913             :   long k, k1, d, m, M, flag, nmax;
     914             :   double a, E, hd, Ep, d2, suma, maxs, mins, sub, B0,B1, Lestimate, Mestimate;
     915             : 
     916       14838 :   Vga = ldata_get_gammavec(ldata);
     917       14838 :   S->d = d = lg(Vga)-1; d2 = d/2.;
     918       14838 :   suma = gtodouble(vecsum(Vga));
     919       14838 :   k = ldata_get_k(ldata);
     920       14838 :   N = ldata_get_conductor(ldata);
     921       14838 :   S->logN2 = log(gtodouble(N)) / 2;
     922       14838 :   maxs = S->dc + S->dw;
     923       14838 :   mins = S->dc - S->dw;
     924       14838 :   S->MAXs = maxss(maxs, k-mins);
     925             : 
     926             :   /* we compute Lambda^(der)(s) / der!; need to compensate for L^(der)(s)
     927             :    * ln |gamma(s)| ~ (pi/4) d |t|; max with 1: fudge factor */
     928       14838 :   S->D = (long)ceil(bitprec + derprec + maxdd((M_PI/(4*M_LN2))*d*S->dh, 1));
     929       14838 :   S->E = E = M_LN2*S->D; /* D:= required absolute bitprec */
     930             : 
     931       14838 :   Ep = E + maxdd(M_PI * S->dh * d2, (d*S->MAXs + suma - 1) * log(E));
     932       14838 :   hd = d2*M_PI*M_PI / Ep;
     933       14838 :   S->m0 = (long)ceil(M_LN2/hd);
     934       14838 :   S->hd = M_LN2/S->m0;
     935             : 
     936       14838 :   S->logC = d2*M_LN2 - log(d2)/2;
     937       14838 :   k1 = ldata_get_k1(ldata);
     938       14838 :   S->k1 = k1; /* assume |a_n| << n^k1 with small implied constant */
     939       14838 :   S->A  = gammavec_expo(d, suma);
     940             : 
     941       14838 :   sub = 0.;
     942       14838 :   if (mins > 1)
     943             :   {
     944        4025 :     GEN sig = dbltor(mins);
     945        4025 :     sub += S->logN2*mins;
     946        4025 :     if (gammaordinary(Vga, sig))
     947             :     {
     948        3941 :       GEN FVga = gammafactor(Vga);
     949        3941 :       GEN gas = gammafactproduct(FVga, sig, LOWDEFAULTPREC);
     950        3941 :       if (typ(gas) != t_SER)
     951             :       {
     952        3941 :         double dg = dbllog2(gas);
     953        3941 :         if (dg > 0) sub += dg * M_LN2;
     954             :       }
     955             :     }
     956             :   }
     957       14838 :   S->sub = sub;
     958       14838 :   M = 1000;
     959       14838 :   L = cgetg(M+2, t_VECSMALL);
     960       14838 :   a = S->k1 + S->A;
     961             : 
     962       14838 :   B0 = 5 + S->E - S->sub + S->logC + S->k1*S->logN2; /* 5 extra bits */
     963       14838 :   B1 = S->hd * (S->MAXs - S->k1);
     964       14838 :   Lestimate = dblcoro526(a + S->MAXs - 2./d, d/2.,
     965       14838 :     S->E - S->sub + S->logC - log(2*M_PI*S->hd) + S->MAXs*S->logN2);
     966       14838 :   Mestimate = ((Lestimate > 0? log(Lestimate): 0) + S->logN2) / S->hd;
     967       14838 :   nmax = 0;
     968       14838 :   flag = 0;
     969     1462853 :   for (m = 0;; m++)
     970     1448015 :   {
     971     1462853 :     double x, H = S->logN2 - m*S->hd, B = B0 + m*B1;
     972             :     long n;
     973     1462853 :     x = dblcoro526(a, d/2., B);
     974     1462853 :     n = floor(x*exp(H));
     975     1462853 :     if (n > nmax) nmax = n;
     976     1462853 :     if (m > M) { M *= 2; L = vecsmall_lengthen(L,M+2); }
     977     1462853 :     L[m+1] = n;
     978     1462853 :     if (n == 0) { if (++flag > 2 && m > Mestimate) break; } else flag = 0;
     979             :   }
     980       14838 :   m -= 2; while (m > 0 && !L[m]) m--;
     981       14838 :   if (m == 0) { nmax = 1; L[1] = 1; m = 1; } /* can happen for tiny bitprec */
     982       14838 :   setlg(L, m+1); S->M = m-1;
     983       14838 :   S->L = L;
     984       14838 :   S->nmax = nmax;
     985             : 
     986       14838 :   S->Dmax = S->D + (long)ceil((S->M * S->hd * S->MAXs - S->sub) / M_LN2);
     987       14838 :   if (S->Dmax < S->D) S->Dmax = S->D;
     988       14838 :   S->precmax = nbits2prec(S->Dmax);
     989       14838 :   if (DEBUGLEVEL > 1)
     990           0 :     err_printf("Dmax=%ld, D=%ld, M = %ld, nmax = %ld, m0 = %ld\n",
     991             :                S->Dmax,S->D,S->M,S->nmax, S->m0);
     992       14838 : }
     993             : 
     994             : /* x0 * [1,x,..., x^n] */
     995             : static GEN
     996        1932 : powersshift(GEN x, long n, GEN x0)
     997             : {
     998        1932 :   long i, l = n+2;
     999        1932 :   GEN V = cgetg(l, t_VEC);
    1000        1932 :   gel(V,1) = x0;
    1001        1932 :   for(i = 2; i < l; i++) gel(V,i) = gmul(gel(V,i-1),x);
    1002        1932 :   return V;
    1003             : }
    1004             : 
    1005             : static GEN
    1006        4557 : lfuninit_pol(GEN vecc, GEN poqk, long M, long prec)
    1007             : {
    1008             :   long m;
    1009        4557 :   GEN pol = cgetg(M+3, t_POL); pol[1] = evalsigne(1) | evalvarn(0);
    1010        4557 :   gel(pol, 2) = gprec_w(gmul2n(gel(vecc,1), -1), prec);
    1011      277032 :   for (m = 2; m <= M+1; m++)
    1012      272475 :     gel(pol, m+1) = gprec_w(gmul(gel(poqk,m), gel(vecc,m)), prec);
    1013        4557 :   return RgX_renormalize_lg(pol, M+3);
    1014             : }
    1015             : 
    1016             : static GEN
    1017        2051 : lfuninit_vecc2_sum(GEN an, GEN qk, GEN a, struct lfunp *Q, GEN poqk)
    1018             : {
    1019        2051 :   const long M = Q->M, prec = Q->precmax;
    1020        2051 :   GEN L = Q->L;
    1021        2051 :   long m, L0 = lg(an)-1;
    1022        2051 :   GEN v = cgetg(M + 2, t_VEC);
    1023        2051 :   if (typ(an) == t_VEC) an = RgV_kill0(an);
    1024      112504 :   for (m = 0; m <= M; m++)
    1025             :   {
    1026      110453 :     pari_sp av = avma;
    1027      110453 :     GEN t = gel(qk, m+1), S = theta2(an, minss(L[m+1],L0), t, a, prec);
    1028      110453 :     gel(v, m+1) = gerepileupto(av, S); /* theta(exp(mh)) */
    1029             :   }
    1030        2051 :   return lfuninit_pol(v, poqk, M, prec);
    1031             : }
    1032             : 
    1033             : /* theta(exp(mh)) ~ sum_{n <= L[m]} a(n) k[m,n] */
    1034             : static GEN
    1035        2506 : lfuninit_vecc_sum(GEN L, long M, GEN vecan, GEN vK, GEN pokq, long prec)
    1036             : {
    1037             :   long m;
    1038        2506 :   GEN vecc = cgetg(M+2, t_VEC);
    1039      169085 :   for (m = 0; m <= M; ++m)
    1040             :   {
    1041      166579 :     pari_sp av = avma;
    1042      166579 :     GEN s = gen_0, vKm = gel(vK,m+1);
    1043             :     long n;
    1044     7373569 :     for (n = 1; n <= L[m+1]; n++)
    1045             :     {
    1046     7206990 :       GEN c = mul_an(vecan, n, gel(vKm,n));
    1047     7206990 :       if (c)
    1048             :       {
    1049     2535757 :         s = gadd(s, c);
    1050     2535757 :         if (gc_needed(av, 3)) s = gerepileupto(av, s);
    1051             :       }
    1052             :     }
    1053      166579 :     gel(vecc,m+1) = gerepileupto(av, s);
    1054             :   }
    1055        2506 :   return lfuninit_pol(vecc, pokq, M, prec);
    1056             : }
    1057             : 
    1058             : /* return [\theta(exp(mh)), m=0..M], theta(t) = sum a(n) K(n/sqrt(N) t),
    1059             :  * h = log(2)/m0 */
    1060             : static GEN
    1061        4417 : lfuninit_vecc(GEN theta, GEN h, struct lfunp *S, GEN poqk)
    1062             : {
    1063        4417 :   const long m0 = S->m0, M = S->M;
    1064        4417 :   GEN tech = linit_get_tech(theta);
    1065             :   GEN va, vK, L, K, d2, vroots, eh2d, peh2d;
    1066        4417 :   GEN sqN = theta_get_sqrtN(tech), an = S->an, bn = S->bn, vprec = S->vprec;
    1067             :   long d, prec, m, n, neval;
    1068             : 
    1069        4417 :   if (!vprec)
    1070             :   { /* d=2 and Vga = [a,a+1] */
    1071        1932 :     GEN ldata = linit_get_ldata(theta);
    1072        1932 :     GEN a = vecmin(ldata_get_gammavec(ldata));
    1073        1932 :     GEN qk = powersshift(mpexp(h), M, ginv(sqN));
    1074        1932 :     va = lfuninit_vecc2_sum(an, qk, a, S, poqk);
    1075        1932 :     return bn? mkvec2(va, lfuninit_vecc2_sum(bn, qk, a, S, poqk)): va;
    1076             :   }
    1077        2485 :   d = S->d;
    1078        2485 :   L = S->L;
    1079        2485 :   prec = S->precmax;
    1080        2485 :   K = theta_get_K(tech);
    1081             : 
    1082             :   /* For all 0<= m <= M, and all n <= L[m+1] such that a_n!=0, we must compute
    1083             :    *   k[m,n] = K(n exp(mh)/sqrt(N))
    1084             :    * with ln(absolute error) <= E + max(mh sigma - sub, 0) + k1 * log(n).
    1085             :    * N.B. we use the 'rt' variant and pass argument (n exp(mh)/sqrt(N))^(2/d).
    1086             :    * Speedup: if n' = 2n and m' = m - m0 >= 0; then k[m,n] = k[m',n']. */
    1087             :   /* vroots[n] = n^(2/d) */
    1088        2485 :   vroots = mkvroots(d, S->nmax, prec);
    1089        2485 :   d2 = gdivgs(gen_2, d);
    1090        2485 :   eh2d = gexp(gmul(d2,h), prec); /* exp(2h/d) */
    1091             :   /* peh2d[m+1] = (exp(mh)/sqrt(N))^(2/d) */
    1092        2485 :   peh2d = gpowers0(eh2d, M, invr(gpow(sqN, d2, prec)));
    1093        2485 :   neval = 0;
    1094             :   /* vK[m+1,n] will contain k[m,n]. For each 0 <= m <= M, sum for n<=L[m+1] */
    1095        2485 :   vK = cgetg(M+2, t_VEC);
    1096      168049 :   for (m = 0; m <= M; m++)
    1097      165564 :     gel(vK,m+1) = const_vec(L[m+1], NULL);
    1098             : 
    1099      168049 :   for (m = M; m >= 0; m--)
    1100     7368459 :     for (n = 1; n <= L[m+1]; n++)
    1101             :     {
    1102     7202895 :       GEN t2d, kmn = gmael(vK,m+1,n);
    1103     7202895 :       long nn, mm, p = 0;
    1104             : 
    1105     7202895 :       if (kmn) continue; /* done already */
    1106             :       /* p = largest (absolute) accuracy to which we need k[m,n] */
    1107    10866898 :       for (mm=m,nn=n; mm>=0 && nn <= L[mm+1]; nn<<=1,mm-=m0)
    1108     7203105 :         if (gel(an, nn) || (bn && gel(bn, nn)))
    1109     7198254 :           p = maxuu(p, umael(vprec,mm+1,nn));
    1110     3663793 :       if (!p) continue; /* a_{n 2^v} = 0 for all v in range */
    1111     3662792 :       t2d = mpmul(gel(vroots, n), gel(peh2d,m+1)); /*(n exp(mh)/sqrt(N))^(2/d)*/
    1112     3662792 :       neval++;
    1113     3662792 :       kmn = gammamellininvrt(K, t2d, p);
    1114    10864686 :       for (mm=m,nn=n; mm>=0 && nn <= L[mm+1]; nn<<=1,mm-=m0)
    1115     7201894 :         gmael(vK,mm+1,nn) = kmn;
    1116             :     }
    1117        2485 :   if (DEBUGLEVEL >= 1) err_printf("true evaluations: %ld\n", neval);
    1118        2485 :   va = lfuninit_vecc_sum(L, M, an, vK, poqk, S->precmax);
    1119        2485 :   return bn? mkvec2(va, lfuninit_vecc_sum(L, M, bn, vK, poqk, S->precmax)): va;
    1120             : }
    1121             : 
    1122             : static void
    1123       77644 : parse_dom(long k, GEN dom, struct lfunp *S)
    1124             : {
    1125       77644 :   long l = lg(dom);
    1126       77644 :   if (typ(dom)!=t_VEC) pari_err_TYPE("lfuninit [domain]", dom);
    1127       77644 :   if (l == 2)
    1128             :   {
    1129       41314 :     S->dc = k/2.;
    1130       41314 :     S->dw = 0.;
    1131       41314 :     S->dh = gtodouble(gel(dom,1));
    1132             :   }
    1133       36330 :   else if (l == 3)
    1134             :   {
    1135         301 :     S->dc = k/2.;
    1136         301 :     S->dw = gtodouble(gel(dom,1));
    1137         301 :     S->dh = gtodouble(gel(dom,2));
    1138             :   }
    1139       36029 :   else if (l == 4)
    1140             :   {
    1141       36029 :     S->dc = gtodouble(gel(dom,1));
    1142       36029 :     S->dw = gtodouble(gel(dom,2));
    1143       36029 :     S->dh = gtodouble(gel(dom,3));
    1144             :   }
    1145             :   else
    1146             :   {
    1147           0 :     pari_err_TYPE("lfuninit [domain]", dom);
    1148           0 :     S->dc = S->dw = S->dh = 0; /*-Wall*/
    1149             :   }
    1150       77644 :   if (S->dw < 0 || S->dh < 0) pari_err_TYPE("lfuninit [domain]", dom);
    1151       77644 : }
    1152             : 
    1153             : /* do we have dom \subset dom0 ? dom = [center, width, height] */
    1154             : int
    1155       13333 : sdomain_isincl(long k, GEN dom, GEN dom0)
    1156             : {
    1157             :   struct lfunp S0, S;
    1158       13333 :   parse_dom(k, dom, &S);
    1159       13333 :   parse_dom(k, dom0, &S0);
    1160       13333 :   return S0.dc - S0.dw <= S.dc - S.dw
    1161       13333 :       && S0.dc + S0.dw >= S.dc + S.dw && S0.dh >= S.dh;
    1162             : }
    1163             : 
    1164             : static int
    1165       13333 : checklfuninit(GEN linit, GEN dom, long der, long bitprec)
    1166             : {
    1167       13333 :   GEN ldata = linit_get_ldata(linit);
    1168       13333 :   GEN domain = lfun_get_domain(linit_get_tech(linit));
    1169       13333 :   return domain_get_der(domain) >= der
    1170       13333 :     && domain_get_bitprec(domain) >= bitprec
    1171       26666 :     && sdomain_isincl(ldata_get_k(ldata), dom, domain_get_dom(domain));
    1172             : }
    1173             : 
    1174             : GEN
    1175        5103 : lfuninit_make(long t, GEN ldata, GEN molin, GEN domain)
    1176             : {
    1177        5103 :   GEN Vga = ldata_get_gammavec(ldata);
    1178        5103 :   long d = lg(Vga)-1;
    1179        5103 :   long k = ldata_get_k(ldata);
    1180        5103 :   GEN k2 = gdivgs(stoi(k), 2);
    1181        5103 :   GEN expot = gdivgs(gadd(gmulsg(d, gsubgs(k2, 1)), vecsum(Vga)), 4);
    1182        5103 :   GEN eno = ldata_get_rootno(ldata);
    1183        5103 :   long prec = nbits2prec( domain_get_bitprec(domain) );
    1184        5103 :   GEN w2 = ginv(gsqrt(eno, prec));
    1185        5103 :   GEN hardy = mkvec4(k2, w2, expot, gammafactor(Vga));
    1186        5103 :   return mkvec3(mkvecsmall(t),ldata, mkvec3(domain, molin, hardy));
    1187             : }
    1188             : 
    1189             : static void
    1190        2485 : lfunparams2(struct lfunp *S)
    1191             : {
    1192        2485 :   GEN vprec, L = S->L, an = S->an, bn = S->bn;
    1193             :   double sig0, pmax, sub2;
    1194        2485 :   long m, nan, nmax, neval, M = S->M;
    1195             : 
    1196             :   /* try to reduce parameters now we know the a_n (some may be 0) */
    1197        2485 :   if (typ(an) == t_VEC) an = RgV_kill0(an);
    1198        2485 :   nan = S->nmax; /* lg(an)-1 may be large than this */
    1199        2485 :   nmax = neval = 0;
    1200        2485 :   if (!bn)
    1201      167013 :     for (m = 0; m <= M; m++)
    1202             :     {
    1203      164549 :       long n = minss(nan, L[m+1]);
    1204      164549 :       while (n > 0 && !gel(an,n)) n--;
    1205      164549 :       if (n > nmax) nmax = n;
    1206      164549 :       neval += n;
    1207      164549 :       L[m+1] = n; /* reduce S->L[m+1] */
    1208             :     }
    1209             :   else
    1210             :   {
    1211          21 :     if (typ(bn) == t_VEC) bn = RgV_kill0(bn);
    1212        1036 :     for (m = 0; m <= M; m++)
    1213             :     {
    1214        1015 :       long n = minss(nan, L[m+1]);
    1215        1015 :       while (n > 0 && !gel(an,n) && !gel(bn,n)) n--;
    1216        1015 :       if (n > nmax) nmax = n;
    1217        1015 :       neval += n;
    1218        1015 :       L[m+1] = n; /* reduce S->L[m+1] */
    1219             :     }
    1220             :   }
    1221        2485 :   if (DEBUGLEVEL >= 1) err_printf("expected evaluations: %ld\n", neval);
    1222        2485 :   for (; M > 0; M--)
    1223        2485 :     if (L[M+1]) break;
    1224        2485 :   setlg(L, M+2);
    1225        2485 :   S->M = M;
    1226        2485 :   S->nmax = nmax;
    1227             : 
    1228        2485 :   pmax = 0;
    1229        2485 :   sig0 = S->MAXs/S->m0;
    1230        2485 :   sub2 = S->sub / M_LN2;
    1231        2485 :   vprec = cgetg(S->M+2, t_VEC);
    1232             :   /* compute accuracy to which we will need k[m,n] = K(n*exp(mh)/sqrt(N))
    1233             :    * vprec[m+1,n] = absolute accuracy to which we need k[m,n] */
    1234      168049 :   for (m = 0; m <= S->M; m++)
    1235             :   {
    1236      165564 :     double c = S->D + maxdd(m*sig0 - sub2, 0);
    1237             :     GEN t;
    1238      165564 :     if (!S->k1)
    1239             :     {
    1240      152922 :       t = const_vecsmall(L[m+1]+1, c);
    1241      152922 :       pmax = maxdd(pmax,c);
    1242             :     }
    1243             :     else
    1244             :     {
    1245             :       long n;
    1246       12642 :       t = cgetg(L[m+1]+1, t_VECSMALL);
    1247     2255022 :       for (n = 1; n <= L[m+1]; n++)
    1248             :       {
    1249     2242380 :         t[n] = c + S->k1 * log2(n);
    1250     2242380 :         pmax = maxdd(pmax, t[n]);
    1251             :       }
    1252             :     }
    1253      165564 :     gel(vprec,m+1) = t;
    1254             :   }
    1255        2485 :   S->vprec = vprec;
    1256        2485 :   S->Dmax = pmax;
    1257        2485 :   S->precmax = nbits2prec(pmax);
    1258        2485 : }
    1259             : 
    1260             : static GEN
    1261        4424 : lfun_init_theta(GEN ldata, GEN eno, struct lfunp *S)
    1262             : {
    1263        4424 :   GEN an2, dual, tdom = NULL, Vga = ldata_get_gammavec(ldata);
    1264             :   long L;
    1265        4424 :   if (eno)
    1266        3346 :     L = S->nmax;
    1267             :   else
    1268             :   {
    1269        1078 :     tdom = dbltor(sqrt(0.5));
    1270        1078 :     L = maxss(S->nmax, lfunthetacost(ldata, tdom, 0, S->D));
    1271             :   }
    1272        4424 :   dual = ldata_get_dual(ldata);
    1273        4424 :   S->an = ldata_vecan(ldata_get_an(ldata), L, S->precmax);
    1274        4417 :   S->bn = typ(dual)==t_INT? NULL: ldata_vecan(dual, S->nmax, S->precmax);
    1275        4417 :   if (!vgaell(Vga)) lfunparams2(S);
    1276             :   else
    1277             :   {
    1278        1932 :     S->an = antwist(S->an, Vga, S->precmax);
    1279        1932 :     if (S->bn) S->bn = antwist(S->bn, Vga, S->precmax);
    1280        1932 :     S->vprec = NULL;
    1281             :   }
    1282        4417 :   an2 = lg(Vga)-1 == 1? antwist(S->an, Vga, S->precmax): S->an;
    1283        4417 :   return lfunthetainit0(ldata, tdom, an2, 0, S->Dmax, 0);
    1284             : }
    1285             : 
    1286             : GEN
    1287       10414 : lfuncost(GEN L, GEN dom, long der, long bitprec)
    1288             : {
    1289       10414 :   pari_sp av = avma;
    1290       10414 :   GEN ldata = lfunmisc_to_ldata_shallow(L);
    1291       10414 :   long k = ldata_get_k(ldata);
    1292             :   struct lfunp S;
    1293             : 
    1294       10414 :   parse_dom(k, dom, &S);
    1295       10414 :   lfunparams(ldata, der, bitprec, &S);
    1296       10414 :   avma = av; return mkvecsmall2(S.nmax, S.Dmax);
    1297             : }
    1298             : GEN
    1299          42 : lfuncost0(GEN L, GEN dom, long der, long bitprec)
    1300             : {
    1301          42 :   pari_sp av = avma;
    1302             :   GEN C;
    1303             : 
    1304          42 :   if (is_linit(L))
    1305             :   {
    1306          28 :     GEN tech = linit_get_tech(L);
    1307          28 :     GEN domain = lfun_get_domain(tech);
    1308          28 :     dom = domain_get_dom(domain);
    1309          28 :     der = domain_get_der(domain);
    1310          28 :     bitprec = domain_get_bitprec(domain);
    1311          28 :     if (linit_get_type(L) == t_LDESC_PRODUCT)
    1312             :     {
    1313          21 :       GEN v = lfunprod_get_fact(linit_get_tech(L)), F = gel(v,1);
    1314          21 :       long i, l = lg(F);
    1315          21 :       C = cgetg(l, t_VEC);
    1316          77 :       for (i = 1; i < l; ++i)
    1317          56 :         gel(C, i) = zv_to_ZV( lfuncost(gel(F,i), dom, der, bitprec) );
    1318          21 :       return gerepileupto(av, C);
    1319             :     }
    1320             :   }
    1321          21 :   if (!dom) pari_err_TYPE("lfuncost [missing s domain]", L);
    1322          21 :   C = lfuncost(L,dom,der,bitprec);
    1323          21 :   return gerepileupto(av, zv_to_ZV(C));
    1324             : }
    1325             : 
    1326             : GEN
    1327       18233 : lfuninit(GEN lmisc, GEN dom, long der, long bitprec)
    1328             : {
    1329       18233 :   pari_sp ltop = avma;
    1330             :   GEN R, h, theta, ldata, qk, poqk, pol, eno, r, domain, molin;
    1331             :   long k;
    1332             :   struct lfunp S;
    1333             : 
    1334       18233 :   if (is_linit(lmisc))
    1335             :   {
    1336       13382 :     long t = linit_get_type(lmisc);
    1337       13382 :     if (t==t_LDESC_INIT || t==t_LDESC_PRODUCT)
    1338             :     {
    1339       13333 :       if (checklfuninit(lmisc, dom, der, bitprec)) return lmisc;
    1340          21 :       pari_warn(warner,"lfuninit: insufficient initialization");
    1341             :     }
    1342             :   }
    1343        4921 :   ldata = lfunmisc_to_ldata_shallow(lmisc);
    1344             : 
    1345        4921 :   if (ldata_get_type(ldata)==t_LFUN_NF)
    1346             :   {
    1347         497 :     GEN T = gel(ldata_get_an(ldata), 2);
    1348         497 :     return lfunzetakinit(T, dom, der, 0, bitprec);
    1349             :   }
    1350        4424 :   k = ldata_get_k(ldata);
    1351        4424 :   parse_dom(k, dom, &S);
    1352        4424 :   lfunparams(ldata, der, bitprec, &S);
    1353        4424 :   r = ldata_get_residue(ldata);
    1354             :   /* Note: all guesses should already have been performed (thetainit more
    1355             :    * expensive than needed: should be either tdom = 1 or bitprec = S.D).
    1356             :    * BUT if the root number / polar part do not have an algebraic
    1357             :    * expression, there is no way to do this until we know the
    1358             :    * precision, i.e. now. So we can't remove guessing code from here and
    1359             :    * lfun_init_theta */
    1360        4424 :   if (r && isintzero(r)) eno = NULL;
    1361             :   else
    1362             :   {
    1363        4424 :     eno = ldata_get_rootno(ldata);
    1364        4424 :     if (isintzero(eno)) eno = NULL;
    1365             :   }
    1366        4424 :   theta = lfun_init_theta(ldata, eno, &S);
    1367        4417 :   if (eno && lg(ldata)==7)
    1368        1995 :     R = gen_0;
    1369             :   else
    1370             :   {
    1371        2422 :     GEN v = lfunrootres(theta, S.D);
    1372        2422 :     ldata = shallowcopy(ldata);
    1373        2422 :     gel(ldata, 6) = gel(v,3);
    1374        2422 :     r = gel(v,1);
    1375        2422 :     if (isintzero(r))
    1376        1064 :       setlg(ldata,7); /* no pole */
    1377             :     else
    1378        1358 :       gel(ldata, 7) = r;
    1379        2422 :     R = lfunrtoR(ldata, nbits2prec(S.D));
    1380             :   }
    1381        4417 :   h = divru(mplog2(S.precmax), S.m0);
    1382        4417 :   k = ldata_get_k(ldata);
    1383        4417 :   qk = gprec_w(mpexp(gmul2n(gmulsg(k,h), -1)), S.precmax); /* exp(kh/2) */
    1384        4417 :   poqk = gpowers(qk, S.M);
    1385        4417 :   pol = lfuninit_vecc(theta, h, &S, poqk);
    1386        4417 :   molin = mkvec3(h, pol, R);
    1387        4417 :   domain = mkvec2(dom, mkvecsmall2(der, bitprec));
    1388        4417 :   return gerepilecopy(ltop, lfuninit_make(t_LDESC_INIT, ldata, molin, domain));
    1389             : }
    1390             : 
    1391             : GEN
    1392         399 : lfuninit0(GEN lmisc, GEN dom, long der, long bitprec)
    1393             : {
    1394         399 :   GEN z = lfuninit(lmisc, dom, der, bitprec);
    1395         399 :   return z == lmisc? gcopy(z): z;
    1396             : }
    1397             : 
    1398             : /* If s is a pole of Lambda, return polar part at s; else return NULL */
    1399             : static GEN
    1400        4205 : lfunpoleresidue(GEN R, GEN s)
    1401             : {
    1402             :   long j;
    1403       11782 :   for (j = 1; j < lg(R); j++)
    1404             :   {
    1405        8088 :     GEN Rj = gel(R, j), be = gel(Rj, 1);
    1406        8088 :     if (gequal(s, be)) return gel(Rj, 2);
    1407             :   }
    1408        3694 :   return NULL;
    1409             : }
    1410             : 
    1411             : /* Compute contribution of polar part at s when not a pole. */
    1412             : static GEN
    1413        6639 : veccothderivn(GEN a, long n)
    1414             : {
    1415             :   long i;
    1416        6639 :   pari_sp av = avma;
    1417        6639 :   GEN c = pol_x(0), cp = mkpoln(3, gen_m1, gen_0, gen_1);
    1418        6639 :   GEN v = cgetg(n+2, t_VEC);
    1419        6639 :   gel(v, 1) = poleval(c, a);
    1420       19980 :   for(i = 2; i <= n+1; i++)
    1421             :   {
    1422       13341 :     c = ZX_mul(ZX_deriv(c), cp);
    1423       13341 :     gel(v, i) = gdiv(poleval(c, a), mpfact(i-1));
    1424             :   }
    1425        6639 :   return gerepilecopy(av, v);
    1426             : }
    1427             : 
    1428             : static GEN
    1429        6702 : polepart(long n, GEN h, GEN C)
    1430             : {
    1431        6702 :   GEN h2n = gpowgs(gdiv(h, gen_2), n-1);
    1432        6702 :   GEN res = gmul(h2n, gel(C,n));
    1433        6702 :   return odd(n)? res : gneg(res);
    1434             : }
    1435             : 
    1436             : static GEN
    1437        3316 : lfunsumcoth(GEN R, GEN s, GEN h, long prec)
    1438             : {
    1439             :   long i,j;
    1440        3316 :   GEN S = gen_0;
    1441        9955 :   for (j = 1; j < lg(R); ++j)
    1442             :   {
    1443        6639 :     GEN r = gel(R,j), be = gel(r,1), Rj = gel(r, 2);
    1444        6639 :     long e = valp(Rj);
    1445        6639 :     GEN z1 = gexpm1(gmul(h, gsub(s,be)), prec); /* exp(h(s-beta))-1 */
    1446        6639 :     GEN c1 = gaddgs(gdivsg(2, z1), 1); /* coth((h/2)(s-beta)) */
    1447        6639 :     GEN C1 = veccothderivn(c1, 1-e);
    1448       13341 :     for (i = e; i < 0; i++)
    1449             :     {
    1450        6702 :       GEN Rbe = mysercoeff(Rj, i);
    1451        6702 :       GEN p1 = polepart(-i, h, C1);
    1452        6702 :       S = gadd(S, gmul(Rbe, p1));
    1453             :     }
    1454             :   }
    1455        3316 :   return gmul2n(S, -1);
    1456             : }
    1457             : 
    1458             : static GEN lfunlambda_OK(GEN linit, GEN s, GEN sdom, long bitprec);
    1459             : /* L is a t_LDESC_PRODUCT Linit */
    1460             : static GEN
    1461        1279 : lfunlambda_product(GEN L, GEN s, GEN sdom, long bitprec)
    1462             : {
    1463        1279 :   GEN ldata = linit_get_ldata(L), v = lfunprod_get_fact(linit_get_tech(L));
    1464        1279 :   GEN r = gen_1, F = gel(v,1), E = gel(v,2), C = gel(v,3), cs = conj_i(s);
    1465        1279 :   long i, l = lg(F), isreal = gequal(imag_i(s), imag_i(cs));
    1466        4439 :   for (i = 1; i < l; ++i)
    1467             :   {
    1468        3160 :     GEN f = lfunlambda_OK(gel(F, i), s, sdom, bitprec);
    1469        3160 :     if (E[i])
    1470        3160 :       r = gmul(r, gpowgs(f, E[i]));
    1471        3160 :     if (C[i])
    1472             :     {
    1473         378 :       GEN fc = isreal? f: lfunlambda_OK(gel(F, i), cs, sdom, bitprec);
    1474         378 :       r = gmul(r, gpowgs(conj_i(fc), C[i]));
    1475             :     }
    1476             :   }
    1477        1279 :   return (ldata_isreal(ldata) && gequal0(imag_i(s)))? real_i(r): r;
    1478             : }
    1479             : 
    1480             : /* s a t_SER */
    1481             : static long
    1482        1086 : der_level(GEN s)
    1483        1086 : { return signe(s)? lg(s)-3: valp(s)-1; }
    1484             : 
    1485             : /* s a t_SER; return coeff(s, X^0) */
    1486             : static GEN
    1487         210 : ser_coeff0(GEN s) { return simplify_shallow(polcoef_i(s, 0, -1)); }
    1488             : 
    1489             : static GEN
    1490        5250 : get_domain(GEN s, GEN *dom, long *der)
    1491             : {
    1492        5250 :   GEN sa = s;
    1493        5250 :   *der = 0;
    1494        5250 :   switch(typ(s))
    1495             :   {
    1496             :     case t_POL:
    1497           7 :     case t_RFRAC: s = toser_i(s);
    1498             :     case t_SER:
    1499         210 :       *der = der_level(s);
    1500         210 :       sa = ser_coeff0(s);
    1501             :   }
    1502        5250 :   *dom = mkvec3(real_i(sa), gen_0, gabs(imag_i(sa),DEFAULTPREC));
    1503        5250 :   return s;
    1504             : }
    1505             : 
    1506             : /* assume lmisc is an linit, s went through get_domain and s/bitprec belong
    1507             :  * to domain */
    1508             : static GEN
    1509       19349 : lfunlambda_OK(GEN linit, GEN s, GEN sdom, long bitprec)
    1510             : {
    1511             :   GEN eno, ldata, tech, h, pol;
    1512       19349 :   GEN S, S0 = NULL, k2, cost;
    1513             :   long prec, prec0;
    1514             :   struct lfunp D, D0;
    1515             : 
    1516       19349 :   if (linit_get_type(linit) == t_LDESC_PRODUCT)
    1517        1279 :     return lfunlambda_product(linit, s, sdom, bitprec);
    1518       18070 :   ldata = linit_get_ldata(linit);
    1519       18070 :   eno = ldata_get_rootno(ldata);
    1520       18070 :   tech = linit_get_tech(linit);
    1521       18070 :   h = lfun_get_step(tech); prec = realprec(h);
    1522             :   /* try to reduce accuracy */
    1523       18070 :   parse_dom(0, sdom, &D0);
    1524       18070 :   parse_dom(0, domain_get_dom(lfun_get_domain(tech)), &D);
    1525       18070 :   if (0.8 * D.dh > D0.dh)
    1526             :   {
    1527       10337 :     cost = lfuncost(linit, sdom, typ(s)==t_SER? der_level(s): 0, bitprec);
    1528       10337 :     prec0 = nbits2prec(cost[2]);
    1529       10337 :     if (prec0 < prec) { prec = prec0; h = gprec_w(h, prec); }
    1530             :   }
    1531       18070 :   pol = lfun_get_pol(tech);
    1532       18070 :   s = gprec_w(s, prec);
    1533       18070 :   if (ldata_get_residue(ldata))
    1534             :   {
    1535        3736 :     GEN R = lfun_get_Residue(tech);
    1536        3736 :     GEN Ra = lfunpoleresidue(R, s);
    1537        3736 :     if (Ra) return gprec_w(Ra, nbits2prec(bitprec));
    1538        3316 :     S0 = lfunsumcoth(R, s, h, prec);
    1539             :   }
    1540       17650 :   k2 = lfun_get_k2(tech);
    1541       17650 :   if (typ(pol)==t_POL && gequal(real_i(s), k2))
    1542       12386 :   { /* on critical line: shortcut */
    1543       12386 :     GEN polz, b = imag_i(s);
    1544       12386 :     polz = gequal0(b)? poleval(pol,gen_1): poleval(pol, expIr(gmul(h,b)));
    1545       12386 :     S = gadd(polz, gmul(eno, conj_i(polz)));
    1546             :   }
    1547             :   else
    1548             :   {
    1549        5264 :     GEN z = gexp(gmul(h, gsub(s, k2)), prec);
    1550        5264 :     GEN zi = ginv(z), zc = conj_i(zi);
    1551        5264 :     if (typ(pol)==t_POL)
    1552        5075 :       S = gadd(poleval(pol, z), gmul(eno, conj_i(poleval(pol, zc))));
    1553             :     else
    1554         189 :       S = gadd(poleval(gel(pol,1), z), gmul(eno, poleval(gel(pol,2), zi)));
    1555             :   }
    1556       17650 :   if (S0) S = gadd(S,S0);
    1557       17650 :   return gprec_w(gmul(S,h), nbits2prec(bitprec));
    1558             : }
    1559             : GEN
    1560         812 : lfunlambda(GEN lmisc, GEN s, long bitprec)
    1561             : {
    1562         812 :   pari_sp av = avma;
    1563             :   GEN linit, dom, z;
    1564             :   long der;
    1565         812 :   s = get_domain(s, &dom, &der);
    1566         812 :   linit = lfuninit(lmisc, dom, der, bitprec);
    1567         812 :   z = lfunlambda_OK(linit,s, dom, bitprec);
    1568         812 :   return gerepilecopy(av, z);
    1569             : }
    1570             : 
    1571             : /* assume lmisc is an linit, s went through get_domain and s/bitprec belong
    1572             :  * to domain */
    1573             : static GEN
    1574        3955 : lfun_OK(GEN linit, GEN s, GEN sdom, long bitprec)
    1575             : {
    1576        3955 :   GEN N, gas, S, FVga, res, ss = s;
    1577        3955 :   long prec = nbits2prec(bitprec);
    1578             : 
    1579        3955 :   FVga = lfun_get_factgammavec(linit_get_tech(linit));
    1580        3955 :   S = lfunlambda_OK(linit, s, sdom, bitprec);
    1581        3955 :   if (typ(S)==t_SER)
    1582             :   {
    1583        1323 :     long d = lg(S) - 2 + fracgammadegree(FVga);
    1584        1323 :     if (typ(s) == t_SER)
    1585         952 :       ss = sertoser(s, d);
    1586             :     else
    1587         371 :       ss = deg1ser_shallow(gen_1, s, varn(S), d);
    1588             :   }
    1589        3955 :   gas = gammafactproduct(FVga, ss, prec);
    1590        3955 :   N = ldata_get_conductor(linit_get_ldata(linit));
    1591        3955 :   res = gdiv(S, gmul(gpow(N, gdivgs(ss, 2), prec), gas));
    1592        3955 :   if (typ(s)!=t_SER && typ(res)==t_SER)
    1593             :   {
    1594         406 :     long v = valp(res);
    1595         406 :     if (v > 0) return gen_0;
    1596         371 :     if (v == 0) res = gel(res, 2);
    1597             :     else
    1598         259 :       setlg(res, minss(lg(res), 2-v));
    1599             :   }
    1600        3920 :   return gprec_w(res, prec);
    1601             : }
    1602             : 
    1603             : GEN
    1604        3178 : lfun(GEN lmisc, GEN s, long bitprec)
    1605             : {
    1606        3178 :   pari_sp av = avma;
    1607             :   GEN linit, dom, z;
    1608             :   long der;
    1609        3178 :   s = get_domain(s, &dom, &der);
    1610        3178 :   linit = lfuninit(lmisc, dom, der, bitprec);
    1611        3171 :   z = lfun_OK(linit, s, dom, bitprec);
    1612        3171 :   return gerepilecopy(av, z);
    1613             : }
    1614             : 
    1615             : /* given a t_SER a+x*s(x), return x*s(x), shallow */
    1616             : static GEN
    1617          42 : sersplit1(GEN s, GEN *head)
    1618             : {
    1619          42 :   long i, l = lg(s);
    1620             :   GEN y;
    1621          42 :   *head = simplify_shallow(mysercoeff(s, 0));
    1622          42 :   if (valp(s) > 0) return s;
    1623          28 :   y = cgetg(l-1, t_SER); y[1] = s[1];
    1624          28 :   setvalp(y, 1);
    1625          28 :   for (i=3; i < l; i++) gel(y,i-1) = gel(s,i);
    1626          28 :   return normalize(y);
    1627             : }
    1628             : 
    1629             : /* n-th derivative of t_SER x, n > 0 */
    1630             : static GEN
    1631         112 : derivnser(GEN x, long n)
    1632             : {
    1633         112 :   long i, vx = varn(x), e = valp(x), lx = lg(x);
    1634             :   GEN y;
    1635         112 :   if (ser_isexactzero(x))
    1636             :   {
    1637           0 :     x = gcopy(x);
    1638           0 :     if (e) setvalp(x,e-n);
    1639           0 :     return x;
    1640             :   }
    1641         210 :   if (e < 0 || e >= n)
    1642             :   {
    1643          98 :     y = cgetg(lx,t_SER);
    1644          98 :     y[1] = evalsigne(1)| evalvalp(e-n) | evalvarn(vx);
    1645         490 :     for (i=0; i<lx-2; i++)
    1646         392 :       gel(y,i+2) = gmul(muls_interval(i+e-n+1,i+e), gel(x,i+2));
    1647             :   } else {
    1648          14 :     if (lx <= n+2) return zeroser(vx, 0);
    1649          14 :     lx -= n;
    1650          14 :     y = cgetg(lx,t_SER);
    1651          14 :     y[1] = evalsigne(1)|_evalvalp(0) | evalvarn(vx);
    1652          70 :     for (i=0; i<lx-2; i++)
    1653          56 :       gel(y,i+2) = gmul(muls_interval(i+1,i+n),gel(x,i+2+n-e));
    1654             :   }
    1655         112 :   return normalize(y);
    1656             : }
    1657             : 
    1658             : /* order of pole of Lambda at s (0 if regular point) */
    1659             : static long
    1660        1834 : lfunlambdaord(GEN linit, GEN s)
    1661             : {
    1662        1834 :   GEN tech = linit_get_tech(linit);
    1663        1834 :   if (linit_get_type(linit)==t_LDESC_PRODUCT)
    1664             :   {
    1665         224 :     GEN v = lfunprod_get_fact(linit_get_tech(linit));
    1666         224 :     GEN F = gel(v, 1), E = gel(v, 2), C = gel(v, 3);
    1667         224 :     long i, ex = 0, l = lg(F);
    1668         840 :     for (i = 1; i < l; i++)
    1669         616 :       ex += lfunlambdaord(gel(F,i), s) * (E[i]+C[i]);
    1670         224 :     return ex;
    1671             :   }
    1672        1610 :   if (ldata_get_residue(linit_get_ldata(linit)))
    1673             :   {
    1674         469 :     GEN r = lfunpoleresidue(lfun_get_Residue(tech), s);
    1675         469 :     if (r) return lg(r)-2;
    1676             :   }
    1677        1519 :   return 0;
    1678             : }
    1679             : 
    1680             : /* derivative of order m > 0 of L (flag = 0) or Lambda (flag = 1) */
    1681             : static GEN
    1682        1267 : lfunderiv(GEN lmisc, long m, GEN s, long flag, long bitprec)
    1683             : {
    1684        1267 :   pari_sp ltop = avma;
    1685        1267 :   GEN res, S = NULL, linit, dom;
    1686        1267 :   long der, prec = nbits2prec(bitprec);
    1687        1267 :   if (m <= 0) pari_err_DOMAIN("lfun", "D", "<=", gen_0, stoi(m));
    1688        1260 :   s = get_domain(s, &dom, &der);
    1689        1260 :   linit = lfuninit(lmisc, dom, der + m, bitprec);
    1690        1260 :   if (typ(s) == t_SER)
    1691             :   {
    1692          42 :     long v, l = lg(s)-1;
    1693             :     GEN sh;
    1694          42 :     if (valp(s) < 0) pari_err_DOMAIN("lfun","valuation", "<", gen_0, s);
    1695          42 :     S = sersplit1(s, &sh);
    1696          42 :     v = valp(S);
    1697          42 :     s = deg1ser_shallow(gen_1, sh, varn(S), m + (l+v-1)/v);
    1698             :   }
    1699             :   else
    1700             :   {
    1701        1218 :     long ex = lfunlambdaord(linit, s);
    1702             :     /* HACK: pretend lfuninit was done to right accuracy */
    1703        1218 :     s = deg1ser_shallow(gen_1, s, 0, m+1+ex);
    1704             :   }
    1705        2044 :   res = flag ? lfunlambda_OK(linit, s, dom, bitprec):
    1706         784 :                lfun_OK(linit, s, dom, bitprec);
    1707        1260 :   if (S)
    1708          42 :     res = gsubst(derivnser(res, m), varn(S), S);
    1709        1218 :   else if (typ(res)==t_SER)
    1710             :   {
    1711        1218 :     long v = valp(res);
    1712        1218 :     if (v > m) { avma = ltop; return gen_0; }
    1713        1211 :     if (v >= 0)
    1714        1141 :       res = gmul(mysercoeff(res, m), mpfact(m));
    1715             :     else
    1716          70 :       res = derivnser(res, m);
    1717             :   }
    1718        1253 :   return gerepilecopy(ltop, gprec_w(res, prec));
    1719             : }
    1720             : 
    1721             : GEN
    1722        1211 : lfunlambda0(GEN lmisc, GEN s, long der, long bitprec)
    1723             : {
    1724        1211 :   return der ? lfunderiv(lmisc, der, s, 1, bitprec):
    1725             :                lfunlambda(lmisc, s, bitprec);
    1726             : }
    1727             : 
    1728             : GEN
    1729        3199 : lfun0(GEN lmisc, GEN s, long der, long bitprec)
    1730             : {
    1731        3199 :   return der ? lfunderiv(lmisc, der, s, 0, bitprec):
    1732             :                lfun(lmisc, s, bitprec);
    1733             : }
    1734             : 
    1735             : GEN
    1736       10841 : lfunhardy(GEN lmisc, GEN t, long bitprec)
    1737             : {
    1738       10841 :   pari_sp ltop = avma;
    1739       10841 :   long prec = nbits2prec(bitprec), k, d;
    1740             :   GEN argz, z, linit, ldata, tech, dom, w2, k2, expot, h, a;
    1741             : 
    1742       10841 :   switch(typ(t))
    1743             :   {
    1744       10834 :     case t_INT: case t_FRAC: case t_REAL: break;
    1745           7 :     default: pari_err_TYPE("lfunhardy",t);
    1746             :   }
    1747             : 
    1748       10834 :   ldata = lfunmisc_to_ldata_shallow(lmisc);
    1749       10834 :   if (!is_linit(lmisc)) lmisc = ldata;
    1750       10834 :   k = ldata_get_k(ldata);
    1751       10834 :   d = ldata_get_degree(ldata);
    1752       10834 :   dom = mkvec3(dbltor(k/2.), gen_0, gabs(t,LOWDEFAULTPREC));
    1753       10834 :   linit = lfuninit(lmisc, dom, 0, bitprec);
    1754       10834 :   tech = linit_get_tech(linit);
    1755       10834 :   w2 = lfun_get_w2(tech);
    1756       10834 :   k2 = lfun_get_k2(tech);
    1757       10834 :   expot = lfun_get_expot(tech);
    1758       10834 :   z = mkcomplex(k2, t);
    1759       10834 :   argz = gatan(gdiv(t, k2), prec); /* more accurate than garg since k/2 \in Q */
    1760             :   /* prec may have increased: don't lose accuracy if |z|^2 is exact */
    1761       10834 :   prec = precision(argz);
    1762       10834 :   a = gsub(gmulsg(d, gmul(t, gmul2n(argz,-1))),
    1763             :            gmul(expot,glog(gnorm(z),prec)));
    1764       10834 :   h = lfunlambda_OK(linit, z, mkvec(t), bitprec);
    1765       10834 :   if (typ(ldata_get_dual(ldata))==t_INT)
    1766       10806 :     h = mulreal(h, w2);
    1767             :   else
    1768          28 :     h = gmul(h, w2);
    1769       10834 :   if (typ(h) == t_COMPLEX && gexpo(imag_i(h)) < -(bitprec >> 1))
    1770           0 :     h = real_i(h);
    1771       10834 :   return gerepileupto(ltop, gmul(h, gexp(a, prec)));
    1772             : }
    1773             : 
    1774             : /* L = log(t); return  \sum_{i = 0}^{v-1}  R[-i-1] L^i/i! */
    1775             : static GEN
    1776        3675 : theta_pole_contrib(GEN R, long v, GEN L)
    1777             : {
    1778        3675 :   GEN s = mysercoeff(R,-v);
    1779             :   long i;
    1780        3836 :   for (i = v-1; i >= 1; i--)
    1781         161 :     s = gadd(mysercoeff(R,-i), gdivgs(gmul(s,L), i));
    1782        3675 :   return s;
    1783             : }
    1784             : /* subtract successively rather than adding everything then subtracting.
    1785             :  * The polar part is "large" and suffers from cancellation: a little stabler
    1786             :  * this way */
    1787             : static GEN
    1788        4270 : theta_add_polar_part(GEN S, GEN R, GEN t, long prec)
    1789             : {
    1790        4270 :   GEN logt = NULL;
    1791        4270 :   long j, l = lg(R);
    1792        7945 :   for (j = 1; j < l; j++)
    1793             :   {
    1794        3675 :     GEN Rj = gel(R,j), b = gel(Rj,1), Rb = gel(Rj,2);
    1795        3675 :     long v = -valp(Rb);
    1796        3675 :     if (v > 1 && !logt) logt = glog(t, prec);
    1797        3675 :     S = gsub(S, gmul(theta_pole_contrib(Rb,v,logt), gpow(t,b,prec)));
    1798             :   }
    1799        4270 :   return S;
    1800             : }
    1801             : 
    1802             : /* Check whether the coefficients, conductor, weight, polar part and root
    1803             :  * number are compatible with the functional equation at t0 and 1/t0.
    1804             :  * Different from lfunrootres. */
    1805             : long
    1806        2562 : lfuncheckfeq(GEN lmisc, GEN t0, long bitprec)
    1807             : {
    1808             :   GEN ldata, theta, thetad, t0i, S0, S0i, w, eno;
    1809             :   long e, prec;
    1810             :   pari_sp av;
    1811             : 
    1812        2562 :   if (is_linit(lmisc) && linit_get_type(lmisc)==t_LDESC_PRODUCT)
    1813             :   {
    1814         112 :     GEN v = lfunprod_get_fact(linit_get_tech(lmisc)), F = gel(v,1);
    1815         112 :     long i, b = -bitprec, l = lg(F);
    1816         413 :     for (i = 1; i < l; i++)
    1817         301 :       b = maxss(b, lfuncheckfeq(gel(F,i), t0, bitprec));
    1818         112 :     return b;
    1819             :   }
    1820        2450 :   av = avma;
    1821        2450 :   prec = nbits2prec(bitprec);
    1822        2450 :   if (!t0)
    1823             :   {
    1824        2380 :     t0 = gadd(gdivgs(mppi(prec), 3), gdivgs(gen_I(), 7));
    1825        2380 :     t0i = ginv(t0);
    1826             :   }
    1827          70 :   else if (gcmpgs(gnorm(t0), 1) < 0)
    1828             :   {
    1829           7 :     t0i = t0;
    1830           7 :     t0 = ginv(t0);
    1831             :   }
    1832             :   else
    1833          63 :     t0i = ginv(t0);
    1834             :   /* |t0| >= 1 */
    1835        2450 :   theta = lfunthetacheckinit(lmisc, t0i, 0, bitprec);
    1836        2450 :   ldata = linit_get_ldata(theta);
    1837        2450 :   thetad = theta_dual(theta, ldata_get_dual(ldata));
    1838        2450 :   if (thetad)
    1839          35 :     S0 = lfuntheta(thetad, t0, 0, bitprec);
    1840             :   else
    1841        2415 :     S0 = conj_i(lfuntheta(theta, conj_i(t0), 0, bitprec));
    1842        2450 :   S0i = lfuntheta(theta, t0i, 0, bitprec);
    1843             : 
    1844        2450 :   eno = ldata_get_rootno(ldata);
    1845        2450 :   if (ldata_get_residue(ldata))
    1846             :   {
    1847         469 :     GEN R = theta_get_R(linit_get_tech(theta));
    1848         469 :     if (gequal0(R))
    1849             :     {
    1850             :       GEN v, r;
    1851          42 :       if (ldata_get_type(ldata) == t_LFUN_NF)
    1852             :       { /* inefficient since theta not needed; no need to optimize for this
    1853             :            (artificial) query [e.g. lfuncheckfeq(t_POL)] */
    1854          21 :         GEN T = gel(ldata_get_an(ldata), 2);
    1855          21 :         GEN L = lfunzetakinit(T,zerovec(3),0,0,bitprec);
    1856          21 :         long e = lfuncheckfeq(L,t0,bitprec);
    1857          21 :         avma = av; return e;
    1858             :       }
    1859          21 :       v = lfunrootres(theta, bitprec);
    1860          21 :       r = gel(v,1);
    1861          21 :       if (gequal0(eno)) eno = gel(v,3);
    1862          21 :       R = lfunrtoR_i(ldata, r, eno, nbits2prec(bitprec));
    1863             :     }
    1864         448 :     S0i = theta_add_polar_part(S0i, R, t0, prec);
    1865             :   }
    1866        2429 :   if (gequal0(S0i) || gequal0(S0)) pari_err_PREC("lfuncheckfeq");
    1867        2429 :   w = gdiv(S0i, gmul(S0, gpowgs(t0, ldata_get_k(ldata))));
    1868             :   /* missing rootno: guess it */
    1869        2429 :   if (gequal0(eno)) eno = lfunrootno(theta, bitprec);
    1870        2429 :   w = gsub(w, eno);
    1871        2429 :   if (thetad) w = gdiv(w, eno); /* |eno| may be large in non-dual case */
    1872        2429 :   e = gexpo(w);
    1873        2429 :   avma = av; return e;
    1874             : }
    1875             : 
    1876             : /*******************************************************************/
    1877             : /*       Compute root number and residues                          */
    1878             : /*******************************************************************/
    1879             : /* round root number to \pm 1 if close to integer. */
    1880             : static GEN
    1881        3843 : ropm1(GEN eno, long prec)
    1882             : {
    1883             :   long e;
    1884        3843 :   GEN r = grndtoi(eno, &e);
    1885        3843 :   return (e < -prec2nbits(prec)/2)? r: eno;
    1886             : }
    1887             : 
    1888             : /* theta for t=1/sqrt(2) and t2==2t simultaneously, saving 25% of the work.
    1889             :  * Assume correct initialization (no thetacheck) */
    1890             : static void
    1891          91 : lfunthetaspec(GEN linit, long bitprec, GEN *pv, GEN *pv2)
    1892             : {
    1893          91 :   pari_sp av = avma;
    1894             :   GEN t, Vga, an, K, ldata, thetainit, v, v2, vroots;
    1895             :   long L, prec, n, d;
    1896             : 
    1897          91 :   ldata = linit_get_ldata(linit);
    1898          91 :   thetainit = linit_get_tech(linit);
    1899          91 :   prec = nbits2prec(bitprec);
    1900          91 :   Vga = ldata_get_gammavec(ldata); d = lg(Vga)-1;
    1901          91 :   if (vgaeasytheta(Vga))
    1902             :   {
    1903          70 :     GEN v2 = sqrtr(real2n(1, nbits2prec(bitprec)));
    1904          70 :     GEN v = shiftr(v2,-1);
    1905          70 :     *pv = lfuntheta(linit, v,  0, bitprec);
    1906          70 :     *pv2= lfuntheta(linit, v2, 0, bitprec);
    1907          70 :     return;
    1908             :   }
    1909          21 :   an = RgV_kill0( theta_get_an(thetainit) );
    1910          21 :   L = lg(an)-1;
    1911             :   /* to compute theta(1/sqrt(2)) */
    1912          21 :   t = ginv(gsqrt(gmul2n(ldata_get_conductor(ldata), 1), prec));
    1913             :   /* t = 1/sqrt(2N) */
    1914             : 
    1915             :   /* From then on, the code is generic and could be used to compute
    1916             :    * theta(t) / theta(2t) without assuming t = 1/sqrt(2) */
    1917          21 :   K = theta_get_K(thetainit);
    1918          21 :   vroots = mkvroots(d, L, prec);
    1919          21 :   t = gpow(t, gdivgs(gen_2, d), prec); /* rt variant: t->t^(2/d) */
    1920             :   /* v = \sum_{n <= L, n odd} a_n K(nt) */
    1921       83342 :   for (v = gen_0, n = 1; n <= L; n+=2)
    1922             :   {
    1923       83321 :     GEN tn, Kn, a = gel(an, n);
    1924             : 
    1925       83321 :     if (!a) continue;
    1926       12299 :     tn = gmul(t, gel(vroots,n));
    1927       12299 :     Kn = gammamellininvrt(K, tn, bitprec);
    1928       12299 :     v = gadd(v, gmul(a,Kn));
    1929             :   }
    1930             :   /* v += \sum_{n <= L, n even} a_n K(nt), v2 = \sum_{n <= L/2} a_n K(2n t) */
    1931       83328 :   for (v2 = gen_0, n = 1; n <= L/2; n++)
    1932             :   {
    1933       83307 :     GEN t2n, K2n, a = gel(an, n), a2 = gel(an,2*n);
    1934             : 
    1935       83307 :     if (!a && !a2) continue;
    1936        9198 :     t2n = gmul(t, gel(vroots,2*n));
    1937        9198 :     K2n = gammamellininvrt(K, t2n, bitprec);
    1938        9198 :     if (a) v2 = gadd(v2, gmul(a, K2n));
    1939        9198 :     if (a2) v = gadd(v,  gmul(a2,K2n));
    1940             :   }
    1941          21 :   *pv = v;
    1942          21 :   *pv2 = v2;
    1943          21 :   gerepileall(av, 2, pv,pv2);
    1944             : }
    1945             : 
    1946             : static GEN
    1947          56 : Rtor(GEN a, GEN R, GEN ldata, long prec)
    1948             : {
    1949          56 :   GEN FVga = gammafactor(ldata_get_gammavec(ldata));
    1950          56 :   GEN Na = gpow(ldata_get_conductor(ldata), gdivgs(a,2), prec);
    1951          56 :   return gdiv(R, gmul(Na, gammafactproduct(FVga, a, prec)));
    1952             : }
    1953             : 
    1954             : /* v = theta~(t), vi = theta(1/t) */
    1955             : static GEN
    1956        3822 : get_eno(GEN R, long k, GEN t, GEN v, GEN vi, long vx, long bitprec)
    1957             : {
    1958        3822 :   long prec = nbits2prec(bitprec);
    1959             :   GEN a0, a1;
    1960        3822 :   GEN S = deg1pol(gmul(gpowgs(t,k), gneg(v)), vi, vx);
    1961             : 
    1962        3822 :   S = theta_add_polar_part(S, R, t, prec);
    1963        3822 :   if (typ(S) != t_POL || degpol(S) != 1) return NULL;
    1964        3822 :   a1 = gel(S,3); if (gexpo(a1) < -bitprec/4) return NULL;
    1965        3787 :   a0 = gel(S,2);
    1966        3787 :   return gdiv(a0, gneg(a1));
    1967             : 
    1968             : }
    1969             : /* Return w using theta(1/t) - w t^k \bar{theta}(t) = polar_part(t,w).
    1970             :  * The full Taylor development of L must be known */
    1971             : GEN
    1972        3787 : lfunrootno(GEN linit, long bitprec)
    1973             : {
    1974             :   GEN ldata, t, eno, v, vi, R, thetad;
    1975        3787 :   long k, prec = nbits2prec(bitprec), vx = fetch_var();
    1976             :   pari_sp av;
    1977             : 
    1978             :   /* initialize for t > 1/sqrt(2) */
    1979        3787 :   linit = lfunthetacheckinit(linit, dbltor(sqrt(0.5)), 0, bitprec);
    1980        3787 :   ldata = linit_get_ldata(linit);
    1981        3787 :   k = ldata_get_k(ldata);
    1982        8904 :   R = ldata_get_residue(ldata)? lfunrtoR_eno(ldata, pol_x(vx), prec)
    1983        5117 :                               : cgetg(1, t_VEC);
    1984        3787 :   t = gen_1;
    1985        3787 :   v = lfuntheta(linit, t, 0, bitprec);
    1986        3787 :   thetad = theta_dual(linit, ldata_get_dual(ldata));
    1987        3787 :   vi = !thetad ? conj_i(v): lfuntheta(thetad, t, 0, bitprec);
    1988        3787 :   eno = get_eno(R,k,t,vi,v, vx, bitprec);
    1989        3787 :   if (!eno && !thetad)
    1990             :   { /* t = sqrt(2), vi = theta(1/t), v = theta(t) */
    1991          35 :     lfunthetaspec(linit, bitprec, &vi, &v);
    1992          35 :     t = sqrtr(utor(2, prec));
    1993          35 :     eno = get_eno(R,k,t,conj_i(v),vi, vx, bitprec);
    1994             :   }
    1995        3787 :   av = avma;
    1996        7574 :   while (!eno)
    1997             :   {
    1998           0 :     t = addsr(1, shiftr(utor(pari_rand(), prec), -66)); /* in [1,1.25[ */
    1999           0 :     v = !thetad ? conj_i(lfuntheta(linit, t, 0, bitprec)):
    2000             :                   lfuntheta(thetad, t, 0, bitprec);
    2001           0 :     vi= lfuntheta(linit, ginv(t), 0, bitprec);
    2002           0 :     eno = get_eno(R,k,t,v,vi, vx, bitprec);
    2003           0 :     avma = av;
    2004             :   }
    2005        3787 :   delete_var(); return ropm1(eno,prec);
    2006             : }
    2007             : 
    2008             : /* Find root number and/or residues when L-function coefficients and
    2009             :    conductor are known. For the moment at most a single residue allowed. */
    2010             : GEN
    2011        2471 : lfunrootres(GEN data, long bitprec)
    2012             : {
    2013        2471 :   pari_sp ltop = avma;
    2014             :   GEN w, r, R, a, b, e, v, v2, be, ldata, linit;
    2015             :   long k, prec;
    2016             : 
    2017        2471 :   ldata = lfunmisc_to_ldata_shallow(data);
    2018        2471 :   r = ldata_get_residue(ldata);
    2019        2471 :   k = ldata_get_k(ldata);
    2020        2471 :   if (r) r = normalize_simple_pole(r, stoi(k));
    2021        2471 :   if (!r || residues_known(r))
    2022             :   {
    2023        2415 :     w = lfunrootno(data, bitprec);
    2024        2415 :     if (!r)
    2025        1085 :       r = R = gen_0;
    2026             :     else
    2027        1330 :       R = lfunrtoR_eno(ldata, w, nbits2prec(bitprec));
    2028        2415 :     return gerepilecopy(ltop, mkvec3(r, R, w));
    2029             :   }
    2030          56 :   linit = lfunthetacheckinit(data, dbltor(sqrt(0.5)), 0, bitprec);
    2031          56 :   prec = nbits2prec(bitprec);
    2032          56 :   if (lg(r) > 2) pari_err_IMPL("multiple poles in lfunrootres");
    2033             :   /* Now residue unknown, and r = [[be,0]]. */
    2034          56 :   be = gmael(r, 1, 1);
    2035          56 :   w = ldata_get_rootno(ldata);
    2036          56 :   if (ldata_isreal(ldata) && gequalm1(w))
    2037           0 :     R = lfuntheta(linit, gen_1, 0, bitprec);
    2038             :   else
    2039             :   {
    2040          56 :     lfunthetaspec(linit, bitprec, &v2, &v);
    2041          56 :     if (gequalgs(gmulsg(2, be), k)) pari_err_IMPL("pole at k/2 in lfunrootres");
    2042          56 :     if (gequalgs(be, k))
    2043             :     {
    2044           7 :       GEN p2k = int2n(k);
    2045           7 :       a = conj_i(gsub(gmul(p2k, v), v2));
    2046           7 :       b = subiu(p2k, 1);
    2047           7 :       e = gmul(gsqrt(p2k, prec), gsub(v2, v));
    2048             :     }
    2049             :     else
    2050             :     {
    2051          49 :       GEN tk2 = gsqrt(int2n(k), prec);
    2052          49 :       GEN tbe = gpow(gen_2, be, prec);
    2053          49 :       GEN tkbe = gpow(gen_2, gdivgs(gsubsg(k, be), 2), prec);
    2054          49 :       a = conj_i(gsub(gmul(tbe, v), v2));
    2055          49 :       b = gsub(gdiv(tbe, tkbe), tkbe);
    2056          49 :       e = gsub(gmul(gdiv(tbe, tk2), v2), gmul(tk2, v));
    2057             :     }
    2058          56 :     if (!isintzero(w)) R = gdiv(gsub(e, gmul(a, w)), b);
    2059             :     else
    2060             :     { /* Now residue unknown, r = [[be,0]], and w unknown. */
    2061           0 :       GEN t0  = mkfrac(stoi(11),stoi(10));
    2062           0 :       GEN th1 = lfuntheta(linit, t0,  0, bitprec);
    2063           0 :       GEN th2 = lfuntheta(linit, ginv(t0), 0, bitprec);
    2064           0 :       GEN tbe = gpow(t0, gmulsg(2, be), prec);
    2065           0 :       GEN tkbe = gpow(t0, gsubsg(k, be), prec);
    2066           0 :       GEN tk2 = gpowgs(t0, k);
    2067           0 :       GEN c = conj_i(gsub(gmul(tbe, th1), th2));
    2068           0 :       GEN d = gsub(gdiv(tbe, tkbe), tkbe);
    2069           0 :       GEN f = gsub(gmul(gdiv(tbe, tk2), th2), gmul(tk2, th1));
    2070           0 :       GEN D = gsub(gmul(a, d), gmul(b, c));
    2071           0 :       w = gdiv(gsub(gmul(d, e), gmul(b, f)), D);
    2072           0 :       R = gdiv(gsub(gmul(a, f), gmul(c, e)), D);
    2073             :     }
    2074             :   }
    2075          56 :   r = normalize_simple_pole(Rtor(be, R, ldata, prec), be);
    2076          56 :   R = lfunrtoR_i(ldata, r, w, prec);
    2077          56 :   return gerepilecopy(ltop, mkvec3(r, R, ropm1(w, prec)));
    2078             : }
    2079             : 
    2080             : /*******************************************************************/
    2081             : /*                           Zeros                                 */
    2082             : /*******************************************************************/
    2083             : struct lhardyz_t {
    2084             :   long bitprec, prec;
    2085             :   GEN linit;
    2086             : };
    2087             : 
    2088             : static GEN
    2089       10344 : lfunhardyzeros(void *E, GEN t)
    2090             : {
    2091       10344 :   struct lhardyz_t *S = (struct lhardyz_t*)E;
    2092       10344 :   long prec = S->prec;
    2093       10344 :   GEN h = lfunhardy(S->linit, t, S->bitprec);
    2094       10344 :   if (typ(h) == t_REAL && realprec(h) < prec) h = gprec_w(h, prec);
    2095       10344 :   return h;
    2096             : }
    2097             : 
    2098             : /* initialize for computation on critical line up to height h, zero
    2099             :  * of order <= m */
    2100             : static GEN
    2101         406 : lfuncenterinit(GEN lmisc, double h, long m, long bitprec)
    2102             : {
    2103         406 :   if (m < 0)
    2104             :   { /* choose a sensible default */
    2105         406 :     if (!is_linit(lmisc) || linit_get_type(lmisc) != t_LDESC_INIT) m = 4;
    2106             :     else
    2107             :     {
    2108         371 :       GEN domain = lfun_get_domain(linit_get_tech(lmisc));
    2109         371 :       m = domain_get_der(domain);
    2110             :     }
    2111             :   }
    2112         406 :   return lfuninit(lmisc, mkvec(dbltor(h)), m, bitprec);
    2113             : }
    2114             : 
    2115             : long
    2116         427 : lfunorderzero(GEN lmisc, long m, long bitprec)
    2117             : {
    2118         427 :   pari_sp ltop = avma;
    2119             :   GEN eno, ldata, linit, k2;
    2120             :   long G, c0, c, st, k;
    2121             : 
    2122         427 :   if (is_linit(lmisc) && linit_get_type(lmisc) == t_LDESC_PRODUCT)
    2123             :   {
    2124          63 :     GEN M = gmael(linit_get_tech(lmisc), 2,1);
    2125             :     long i;
    2126          63 :     for (c=0,i=1; i < lg(M); i++) c += lfunorderzero(gel(M,i), m, bitprec);
    2127          63 :     return c;
    2128             :   }
    2129         364 :   linit = lfuncenterinit(lmisc, 0, m, bitprec);
    2130         364 :   ldata = linit_get_ldata(linit);
    2131         364 :   eno = ldata_get_rootno(ldata);
    2132         364 :   G = -bitprec/2;
    2133         364 :   c0 = 0; st = 1;
    2134         364 :   if (ldata_isreal(ldata))
    2135             :   {
    2136         301 :     if (!gequal1(eno)) c0 = 1;
    2137         301 :     st = 2;
    2138             :   }
    2139         364 :   k = ldata_get_k(ldata);
    2140         364 :   k2 = gdivgs(stoi(k), 2);
    2141         378 :   for (c = c0;; c += st)
    2142         392 :     if (gexpo(lfun0(linit, k2, c, bitprec)) > G) break;
    2143         364 :   avma = ltop; return c;
    2144             : }
    2145             : 
    2146             : GEN
    2147          42 : lfunzeros(GEN ldata, GEN lim, long divz, long bitprec)
    2148             : {
    2149          42 :   pari_sp ltop = avma;
    2150             :   GEN ldataf, linit, N, pi2, cN, pi2div, w, T, Vga, h1, h2;
    2151          42 :   long i, d, W, NEWD, precinit, ct, s, prec = nbits2prec(bitprec);
    2152             :   double maxt;
    2153             :   GEN maxtr, maxtr1;
    2154             :   struct lhardyz_t S;
    2155             : 
    2156          42 :   if (typ(lim) == t_VEC)
    2157             :   {
    2158           7 :     if (lg(lim) != 3 || gcmp(gel(lim,1),gel(lim,2)) >= 0
    2159           7 :                      || gcmp(gel(lim,1),gen_0) <= 0)
    2160           0 :       pari_err_TYPE("lfunzeros",lim);
    2161           7 :     h1 = gel(lim,1); h2 = gel(lim,2);
    2162             :   }
    2163             :   else
    2164             :   {
    2165          35 :     if (gcmp(lim,gen_0) <= 0)
    2166           0 :       pari_err_TYPE("lfunzeros",lim);
    2167          35 :     h1 = gen_0; h2 = lim;
    2168             :   }
    2169          42 :   maxt = gtodouble(h2);
    2170             : 
    2171          42 :   if (is_linit(ldata) && linit_get_type(ldata) == t_LDESC_PRODUCT)
    2172             :   {
    2173           0 :     GEN v, M = gmael(linit_get_tech(ldata), 2,1);
    2174           0 :     long l = lg(M);
    2175           0 :     v = cgetg(l, t_VEC);
    2176           0 :     for (i = 1; i < l; i++)
    2177           0 :       gel(v,i) = lfunzeros(gel(M,i), lim, divz, bitprec);
    2178           0 :     return gerepileupto(ltop, vecsort0(shallowconcat1(v), NULL, 0));
    2179             :   }
    2180          42 :   S.linit = linit = lfuncenterinit(ldata, maxt + 1, -1, bitprec);
    2181          42 :   S.bitprec = bitprec;
    2182          42 :   S.prec = prec;
    2183          42 :   ldataf = linit_get_ldata(linit);
    2184          42 :   Vga = ldata_get_gammavec(ldataf); d = lg(Vga) - 1;
    2185          42 :   N = ldata_get_conductor(ldataf);
    2186          42 :   NEWD = minss((long) ceil(bitprec+(M_PI/(4*M_LN2))*d*maxt),
    2187             :                lfun_get_bitprec(linit_get_tech(linit)));
    2188          42 :   precinit = prec; prec = nbits2prec(NEWD);
    2189          42 :   pi2 = Pi2n(1, prec);
    2190          42 :   cN = gdiv(N, gpowgs(Pi2n(-1, prec), d));
    2191          42 :   cN = gexpo(cN) >= 0? gaddsg(d, gmulsg(2, glog(cN, prec))): stoi(d);
    2192          42 :   pi2div = gdivgs(pi2, labs(divz));
    2193          42 :   ct = 0;
    2194          42 :   T = h1;
    2195          42 :   if (gequal0(h1))
    2196             :   {
    2197          35 :     GEN r = ldata_get_residue(ldataf);
    2198          35 :     if (!r || gequal0(r))
    2199             :     {
    2200          21 :       ct = lfunorderzero(linit, -1, bitprec);
    2201          21 :       if (ct) T = real2n(-prec2nbits(prec)/(2*ct), prec);
    2202             :     }
    2203             :   }
    2204             :   /* initialize for 100 further zeros, double later if needed */
    2205          42 :   W = 100 + ct; w = cgetg(W+1,t_VEC);
    2206          42 :   for (i=1; i<=ct; i++) gel(w,i) = gen_0;
    2207          42 :   s = gsigne(lfunhardyzeros(&S, T));
    2208          42 :   maxtr = h2; maxtr1 = gaddsg(1, maxtr);
    2209         427 :   while (gcmp(T, maxtr1) < 0)
    2210             :   {
    2211         385 :     pari_sp av = avma;
    2212         385 :     GEN T0 = T, z;
    2213             :     for(;;)
    2214        5705 :     {
    2215             :       long s0;
    2216             :       GEN L;
    2217        6090 :       if (gcmp(T, pi2) >= 0)
    2218        4221 :         L = gadd(cN, gmulsg(d, glog(gdiv(T, pi2), prec)));
    2219             :       else
    2220        1869 :         L = cN;
    2221        6090 :       T = gadd(T, gdiv(pi2div, L));
    2222        6090 :       if (gcmp(T, maxtr1) > 0) goto END;
    2223        6069 :       s0 = gsigne(lfunhardyzeros(&S, T));
    2224        6069 :       if (s0 != s) { s = s0; break; }
    2225             :     }
    2226         364 :     T = gerepileupto(av, T);
    2227         364 :     z = zbrent(&S, lfunhardyzeros, T0, T, prec);
    2228         364 :     if (gcmp(z, maxtr) > 0) break;
    2229         343 :     if (typ(z) == t_REAL) z  = rtor(z, precinit);
    2230             :     /* room for twice as many zeros */
    2231         343 :     if (ct >= W) { W *= 2; w = vec_lengthen(w, W); }
    2232         343 :     gel(w, ++ct) = z;
    2233             :   }
    2234             : END:
    2235          42 :   setlg(w, ct+1); return gerepilecopy(ltop, w);
    2236             : }
    2237             : 
    2238             : /*******************************************************************/
    2239             : /*       Guess conductor                                           */
    2240             : /*******************************************************************/
    2241             : struct huntcond_t {
    2242             :   long k;
    2243             :   GEN data, thetad;
    2244             :   GEN *pM, *psqrtM, *pMd, *psqrtMd;
    2245             : };
    2246             : 
    2247             : /* M should eventually converge to N, the conductor. L has no pole. */
    2248             : static GEN
    2249        6473 : wrap1(void *E, GEN M)
    2250             : {
    2251        6473 :   struct huntcond_t *S = (struct huntcond_t*)E;
    2252        6473 :   GEN data = S->data, thetainit, tk, p1, p1inv;
    2253        6473 :   GEN t = mkfrac(stoi(11), stoi(10));
    2254             :   long prec, bitprec;
    2255             : 
    2256        6473 :   thetainit = linit_get_tech(data);
    2257        6473 :   bitprec = theta_get_bitprec(thetainit);
    2258        6473 :   prec = nbits2prec(bitprec);
    2259        6473 :   *(S->pM) = M;
    2260        6473 :   *(S->psqrtM) = gsqrt(M, prec);
    2261        6473 :   tk = gpowgs(t, S->k);
    2262        6473 :   if (S->thetad)
    2263             :   {
    2264           0 :     *(S->pMd) = M;
    2265           0 :     *(S->psqrtMd) = *(S->psqrtM);
    2266           0 :     p1 = lfuntheta(S->thetad, t, 0, bitprec);
    2267             :   }
    2268             :   else
    2269        6473 :     p1 = lfuntheta(data, t, 0, bitprec);
    2270        6473 :   p1inv = lfuntheta(data, ginv(t), 0, bitprec);
    2271        6473 :   return glog(gabs(gmul(tk, gdiv(p1, p1inv)), prec), prec);
    2272             : }
    2273             : 
    2274             : /* M should eventually converge to N, the conductor. L has a pole. */
    2275             : static GEN
    2276        2194 : wrap2(void *E, GEN M)
    2277             : {
    2278        2194 :   struct huntcond_t *S = (struct huntcond_t*)E;
    2279        2194 :   GEN data = S->data, t1k, t2k, p1, p1inv, p2, p2inv;
    2280             :   GEN thetainit, R;
    2281        2194 :   GEN t1 = mkfrac(stoi(11), stoi(10)), t2 = mkfrac(stoi(13), stoi(11));
    2282             :   GEN t1be, t2be, t1bemk, t2bemk, t1kmbe, t2kmbe;
    2283             :   GEN F11, F12, F21, F22, P1, P2, res;
    2284        2194 :   long k = S->k, prec, bitprec;
    2285             : 
    2286        2194 :   thetainit = linit_get_tech(data);
    2287        2194 :   bitprec = theta_get_bitprec(thetainit);
    2288        2194 :   prec = nbits2prec(bitprec);
    2289        2194 :   *(S->pM) = M;
    2290        2194 :   *(S->psqrtM) = gsqrt(M, prec);
    2291             : 
    2292        2194 :   p1 = lfuntheta(data, t1, 0, bitprec);
    2293        2194 :   p2 = lfuntheta(data, t2, 0, bitprec);
    2294        2194 :   p1inv = lfuntheta(data, ginv(t1), 0, bitprec);
    2295        2194 :   p2inv = lfuntheta(data, ginv(t2), 0, bitprec);
    2296        2194 :   t1k = gpowgs(t1, k);
    2297        2194 :   t2k = gpowgs(t2, k);
    2298        2194 :   R = theta_get_R(thetainit);
    2299        2194 :   if (typ(R) == t_VEC)
    2300             :   {
    2301        2194 :     GEN be = gmael(R, 1, 1);
    2302        2194 :     t1be = gpow(t1, be, prec); t1bemk = gdiv(gsqr(t1be), t1k);
    2303        2194 :     t2be = gpow(t2, be, prec); t2bemk = gdiv(gsqr(t2be), t2k);
    2304        2194 :     t1kmbe = gdiv(t1k, t1be);
    2305        2194 :     t2kmbe = gdiv(t2k, t2be);
    2306             :   }
    2307             :   else
    2308             :   { /* be = k */
    2309           0 :     t1be = t1k; t1bemk = t1k; t1kmbe = gen_1;
    2310           0 :     t2be = t2k; t2bemk = t2k; t2kmbe = gen_1;
    2311             :   }
    2312        2194 :   F11 = conj_i(gsub(gmul(gsqr(t1be), p1), p1inv));
    2313        2194 :   F12 = conj_i(gsub(gmul(gsqr(t2be), p2), p2inv));
    2314        2194 :   F21 = gsub(gmul(t1k, p1), gmul(t1bemk, p1inv));
    2315        2194 :   F22 = gsub(gmul(t2k, p2), gmul(t2bemk, p2inv));
    2316        2194 :   P1 = gsub(gmul(t1bemk, t1be), t1kmbe);
    2317        2194 :   P2 = gsub(gmul(t2bemk, t2be), t2kmbe);
    2318        2194 :   res = gdiv(gsub(gmul(P2,F21), gmul(P1,F22)),
    2319             :              gsub(gmul(P2,F11), gmul(P1,F12)));
    2320        2194 :   return glog(gabs(res, prec), prec);
    2321             : }
    2322             : 
    2323             : /* If flag = 0 (default) return all conductors found as integers. If
    2324             : flag = 1, return the approximations, not the integers. If flag = 2,
    2325             : return all, even nonintegers. */
    2326             : 
    2327             : static GEN
    2328          84 : checkconductor(GEN v, long bit, long flag)
    2329             : {
    2330             :   GEN w;
    2331          84 :   long e, j, k, l = lg(v);
    2332          84 :   if (flag == 2) return v;
    2333          84 :   w = cgetg(l, t_VEC);
    2334         301 :   for (j = k = 1; j < l; j++)
    2335             :   {
    2336         217 :     GEN N = grndtoi(gel(v,j), &e);
    2337         217 :     if (e < -bit) gel(w,k++) = flag ? gel(v,j): N;
    2338             :   }
    2339          84 :   if (k == 2) return gel(w,1);
    2340           7 :   setlg(w,k); return w;
    2341             : }
    2342             : 
    2343             : static void
    2344          84 : parse_maxcond(GEN maxcond, GEN *pm, GEN *pM)
    2345             : {
    2346          84 :   GEN m = gen_1, M;
    2347          84 :   if (!maxcond)
    2348          49 :     M = utoipos(10000);
    2349          35 :   else if (typ(maxcond) == t_VEC)
    2350             :   {
    2351           7 :     if (lg(maxcond) != 3) pari_err_TYPE("lfunconductor", maxcond);
    2352           7 :     m = gel(maxcond,1);
    2353           7 :     M = gel(maxcond,2);
    2354             :   }
    2355             :   else
    2356          28 :     M = maxcond;
    2357          84 :   m = (typ(m) == t_INT)? gsub(m,ghalf): gfloor(m);
    2358          84 :   if (signe(m) <= 0) m = ghalf;
    2359          84 :   M = (typ(M) == t_INT)? addiu(M, 1): gceil(M);
    2360          84 :   *pm = m;
    2361          84 :   *pM = M;
    2362          84 : }
    2363             : 
    2364             : GEN
    2365          84 : lfunconductor(GEN data, GEN maxcond, long flag, long bitprec)
    2366             : {
    2367             :   struct huntcond_t S;
    2368          84 :   pari_sp ltop = avma;
    2369             :   GEN ld, r, v, ldata, theta, thetad, m, M, tdom;
    2370             :   GEN (*eval)(void *, GEN);
    2371          84 :   bitprec = 3*bitprec/2;
    2372          84 :   ldata = lfunmisc_to_ldata_shallow(data);
    2373          84 :   parse_maxcond(maxcond, &m,&M);
    2374          84 :   r = ldata_get_residue(ldata);
    2375          84 :   if (r && typ(r) == t_VEC)
    2376             :   {
    2377           0 :     if (lg(r) > 2) pari_err_IMPL("multiple poles in lfunconductor");
    2378           0 :     r = gmael(r,1,2);
    2379             :   }
    2380          84 :   if (!r)
    2381          63 :   { eval = wrap1; tdom = mkfrac(stoi(10), stoi(11)); }
    2382             :   else
    2383          21 :   { eval = wrap2; tdom = mkfrac(stoi(11), stoi(13)); }
    2384          84 :   ld = shallowcopy(ldata);
    2385          84 :   gel(ld, 5) = M;
    2386          84 :   theta = lfunthetainit_i(ld, tdom, 0, bitprec);
    2387          84 :   thetad = theta_dual(theta, ldata_get_dual(ldata));
    2388          84 :   gel(theta,3) = shallowcopy(linit_get_tech(theta));
    2389          84 :   S.k = ldata_get_k(ldata);
    2390          84 :   S.data = theta;
    2391          84 :   S.thetad = thetad;
    2392          84 :   S.pM = &gel(linit_get_ldata(theta),5);
    2393          84 :   S.psqrtM = &gel(linit_get_tech(theta),7);
    2394          84 :   if (thetad)
    2395             :   {
    2396           0 :     S.pMd = &gel(linit_get_ldata(thetad),5);
    2397           0 :     S.psqrtMd = &gel(linit_get_tech(thetad),7);
    2398             :   }
    2399          84 :   v = solvestep((void*)&S, eval, m, M, gen_2, 14, nbits2prec(bitprec));
    2400          84 :   return gerepilecopy(ltop, checkconductor(v, bitprec/2, flag));
    2401             : }
    2402             : 
    2403             : /* assume chi primitive */
    2404             : static GEN
    2405         588 : znchargauss_i(GEN G, GEN chi, long bitprec)
    2406             : {
    2407         588 :   GEN z, q, F = znstar_get_N(G);
    2408             :   long prec;
    2409             : 
    2410         588 :   if (equali1(F)) return gen_1;
    2411         343 :   prec = nbits2prec(bitprec);
    2412         343 :   q = sqrtr_abs(itor(F, prec));
    2413         343 :   z = lfuntheta(mkvec2(G,chi), gen_1, 0, bitprec);
    2414         343 :   if (gexpo(z) < 10 - bitprec)
    2415             :   {
    2416          28 :     if (equaliu(F,300))
    2417             :     {
    2418          14 :       GEN z = rootsof1u_cx(25, prec);
    2419          14 :       GEN n = znconreyexp(G, chi);
    2420          14 :       if (equaliu(n, 131)) return gmul(q, gpowgs(z,14));
    2421           7 :       if (equaliu(n, 71)) return gmul(q, gpowgs(z,11));
    2422             :     }
    2423          14 :     if (equaliu(F,600))
    2424             :     {
    2425          14 :       GEN z = rootsof1u_cx(25, prec);
    2426          14 :       GEN n = znconreyexp(G, chi);
    2427          14 :       if (equaliu(n, 491)) return gmul(q, gpowgs(z,7));
    2428           7 :       if (equaliu(n, 11)) return gmul(q, gpowgs(z,18));
    2429             :     }
    2430           0 :     pari_err_BUG("znchargauss [ Theta(chi,1) = 0 ]");
    2431             :   }
    2432         315 :   z = gmul(gdiv(z, conj_i(z)), q);
    2433         315 :   if (zncharisodd(G,chi)) z = mulcxI(z);
    2434         315 :   return z;
    2435             : }
    2436             : static GEN
    2437         588 : Z_radical(GEN N, long *om)
    2438             : {
    2439         588 :   GEN P = gel(Z_factor(N), 1);
    2440         588 :   *om = lg(P)-1; return ZV_prod(P);
    2441             : }
    2442             : GEN
    2443         805 : znchargauss(GEN G, GEN chi, GEN a, long bitprec)
    2444             : {
    2445             :   GEN v, T, N, F, b0, b1, b2, bF, a1, aF, A, r, GF, tau, B, faB, u, S;
    2446         805 :   long omb0, prec = nbits2prec(bitprec);
    2447         805 :   pari_sp av = avma;
    2448             : 
    2449         805 :   if (typ(chi) != t_COL) chi = znconreylog(G,chi);
    2450         805 :   T = znchartoprimitive(G, chi);
    2451         805 :   GF  = gel(T,1);
    2452         805 :   chi = gel(T,2); /* now primitive */
    2453         805 :   N = znstar_get_N(G);
    2454         805 :   F = znstar_get_N(GF);
    2455         805 :   if (equalii(N,F)) b1 = bF = gen_1;
    2456             :   else
    2457             :   {
    2458         231 :     v = Z_ppio(diviiexact(N,F), F);
    2459         231 :     bF = gel(v,2); /* (N/F, F^oo) */
    2460         231 :     b1 = gel(v,3); /* cofactor */
    2461             :   }
    2462         805 :   if (!a) a = a1 = aF = gen_1;
    2463             :   else
    2464             :   {
    2465         756 :     if (typ(a) != t_INT) pari_err_TYPE("znchargauss",a);
    2466         756 :     a = modii(a, N);
    2467         756 :     v = Z_ppio(a, F);
    2468         756 :     aF = gel(v,2);
    2469         756 :     a1 = gel(v,3);
    2470             :   }
    2471         805 :   if (!equalii(aF, bF)) { avma = av; return gen_0; }
    2472         588 :   b0 = Z_radical(b1, &omb0);
    2473         588 :   b2 = diviiexact(b1, b0);
    2474         588 :   A = dvmdii(a1, b2, &r);
    2475         588 :   if (r != gen_0) { avma = av; return gen_0; }
    2476         588 :   B = gcdii(A,b0); faB = Z_factor(B); /* squarefree */
    2477         588 :   S = eulerphi(mkvec2(B,faB));
    2478         588 :   if (odd(omb0 + lg(gel(faB,1))-1)) S = negi(S); /* moebius(b0/B) * phi(B) */
    2479         588 :   S = mulii(S, mulii(aF,b2));
    2480         588 :   tau = znchargauss_i(GF, chi, bitprec);
    2481         588 :   u = Fp_div(b0, A, F);
    2482         588 :   if (!equali1(u))
    2483             :   {
    2484         252 :     GEN ord = zncharorder(GF, chi), z = rootsof1_cx(ord, prec);
    2485         252 :     tau = gmul(tau, znchareval(GF, chi, u, mkvec2(z,ord)));
    2486             :   }
    2487         588 :   return gerepileupto(av, gmul(tau, S));
    2488             : }

Generated by: LCOV version 1.13