Code coverage tests

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

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

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

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

Generated by: LCOV version 1.13