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 - lfunutils.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.12.1 lcov report (development 24038-ebe36f6c4) Lines: 1387 1533 90.5 %
Date: 2019-07-23 05:53:17 Functions: 129 140 92.1 %
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: Applications                      **/
      17             : /**                                                                **/
      18             : /********************************************************************/
      19             : 
      20             : #include "pari.h"
      21             : #include "paripriv.h"
      22             : 
      23             : static GEN
      24        8687 : tag(GEN x, long t) { return mkvec2(mkvecsmall(t), x); }
      25             : 
      26             : /* v a t_VEC of length > 1 */
      27             : static int
      28       36858 : is_tagged(GEN v)
      29             : {
      30       36858 :   GEN T = gel(v,1);
      31       36858 :   return (typ(T)==t_VEC && lg(T)==3 && typ(gel(T,1))==t_VECSMALL);
      32             : }
      33             : static void
      34       36858 : checkldata(GEN ldata)
      35             : {
      36             :   GEN vga, w, N;
      37             : #if 0 /* assumed already checked and true */
      38             :   long l = lg(ldata);
      39             :   if (typ(ldata)!=t_VEC || l < 7 || l > 8 || !is_tagged(ldata))
      40             :     pari_err_TYPE("checkldata", ldata);
      41             : #endif
      42       36858 :   vga = ldata_get_gammavec(ldata);
      43       36858 :   if (typ(vga) != t_VEC) pari_err_TYPE("checkldata [gammavec]",vga);
      44       36858 :   w = gel(ldata, 4); /* FIXME */
      45       36858 :   switch(typ(w))
      46             :   {
      47       36634 :     case t_INT: break;
      48          28 :     case t_FRAC: break;
      49         196 :     case t_VEC: if (lg(w) == 3 && typ(gel(w,1)) == t_INT) break;
      50           0 :     default: pari_err_TYPE("checkldata [weight]",w);
      51             :   }
      52       36858 :   N = ldata_get_conductor(ldata);
      53       36858 :   if (typ(N) != t_INT) pari_err_TYPE("checkldata [conductor]",N);
      54       36858 : }
      55             : 
      56             : /* data may be either an object (polynomial, elliptic curve, etc...)
      57             :  * or a description vector [an,sd,Vga,k,conductor,rootno,{poles}]. */
      58             : GEN
      59         973 : lfuncreate(GEN data)
      60             : {
      61         973 :   long lx = lg(data);
      62         973 :   if (typ(data)==t_VEC && (lx == 7 || lx == 8))
      63             :   {
      64             :     GEN ldata;
      65         378 :     if (is_tagged(data)) ldata = gcopy(data);
      66             :     else
      67             :     { /* tag first component as t_LFUN_GENERIC */
      68         280 :       ldata = gcopy(data);
      69         280 :       gel(ldata, 1) = tag(gel(ldata,1), t_LFUN_GENERIC);
      70         280 :       if (typ(gel(ldata, 2))!=t_INT)
      71          21 :         gel(ldata, 2) = tag(gel(ldata,2), t_LFUN_GENERIC);
      72             :     }
      73         378 :     checkldata(ldata); return ldata;
      74             :   }
      75         595 :   return lfunmisc_to_ldata(data);
      76             : }
      77             : 
      78             : /********************************************************************/
      79             : /**                     Simple constructors                        **/
      80             : /********************************************************************/
      81             : 
      82             : static GEN
      83           0 : vecan_conj(GEN an, long n, long prec)
      84             : {
      85           0 :   GEN p1 = ldata_vecan(gel(an,1), n, prec);
      86           0 :   return typ(p1) == t_VEC? conj_i(p1): p1;
      87             : }
      88             : 
      89             : static GEN
      90         224 : vecan_mul(GEN an, long n, long prec)
      91             : {
      92         224 :   GEN p1 = ldata_vecan(gel(an,1), n, prec);
      93         224 :   GEN p2 = ldata_vecan(gel(an,2), n, prec);
      94         224 :   if (typ(p1) == t_VECSMALL) p1 = vecsmall_to_vec(p1);
      95         224 :   if (typ(p2) == t_VECSMALL) p2 = vecsmall_to_vec(p2);
      96         224 :   return dirmul(p1, p2);
      97             : }
      98             : 
      99             : static GEN
     100          42 : lfunconvol(GEN a1, GEN a2)
     101          42 : { return tag(mkvec2(a1, a2), t_LFUN_MUL); }
     102             : 
     103             : static GEN
     104         623 : vecan_div(GEN an, long n, long prec)
     105             : {
     106         623 :   GEN p1 = ldata_vecan(gel(an,1), n, prec);
     107         623 :   GEN p2 = ldata_vecan(gel(an,2), n, prec);
     108         623 :   if (typ(p1) == t_VECSMALL) p1 = vecsmall_to_vec(p1);
     109         623 :   if (typ(p2) == t_VECSMALL) p2 = vecsmall_to_vec(p2);
     110         623 :   return dirdiv(p1, p2);
     111             : }
     112             : 
     113             : static GEN
     114          49 : lfunconvolinv(GEN a1, GEN a2)
     115          49 : { return tag(mkvec2(a1,a2), t_LFUN_DIV); }
     116             : 
     117             : static GEN
     118           0 : lfunconj(GEN a1)
     119           0 : { return tag(mkvec(a1), t_LFUN_CONJ); }
     120             : 
     121             : static GEN
     122          91 : lfuncombdual(GEN fun(GEN, GEN), GEN ldata1, GEN ldata2)
     123             : {
     124          91 :   GEN a1 = ldata_get_an(ldata1), a2 = ldata_get_an(ldata2);
     125          91 :   GEN b1 = ldata_get_dual(ldata1), b2 = ldata_get_dual(ldata2);
     126          91 :   if (typ(b1)==t_INT && typ(b2)==t_INT)
     127          91 :     return utoi(signe(b1) && signe(b2));
     128             :   else
     129             :   {
     130           0 :     if (typ(b1)==t_INT) b1 = signe(b1) ? lfunconj(a1): a1;
     131           0 :     if (typ(b2)==t_INT) b2 = signe(b2) ? lfunconj(a2): a2;
     132           0 :     return fun(b1, b2);
     133             :   }
     134             : }
     135             : 
     136             : static GEN
     137         357 : vecan_twist(GEN an, long n, long prec)
     138             : {
     139         357 :   GEN p1 = ldata_vecan(gel(an,1), n, prec);
     140         357 :   GEN p2 = ldata_vecan(gel(an,2), n, prec);
     141             :   long i;
     142             :   GEN V;
     143         357 :   if (typ(p1) == t_VECSMALL) p1 = vecsmall_to_vec(p1);
     144         357 :   if (typ(p2) == t_VECSMALL) p2 = vecsmall_to_vec(p2);
     145         357 :   V = cgetg(n+1, t_VEC);
     146      223720 :   for(i = 1; i <= n ; i++)
     147      223363 :     gel(V, i) = gmul(gel(p1, i), gel(p2, i));
     148         357 :   return V;
     149             : }
     150             : 
     151             : static GEN
     152          42 : lfunmulpoles(GEN ldata1, GEN ldata2, long bitprec)
     153             : {
     154             :   long l, j;
     155          42 :   GEN k = ldata_get_k(ldata1);
     156          42 :   GEN r1 = ldata_get_residue(ldata1);
     157          42 :   GEN r2 = ldata_get_residue(ldata2), r;
     158             : 
     159          42 :   if (r1 && typ(r1) != t_VEC) r1 = mkvec(mkvec2(k, r1));
     160          42 :   if (r2 && typ(r2) != t_VEC) r2 = mkvec(mkvec2(k, r2));
     161          42 :   if (!r1)
     162             :   {
     163          14 :     if (!r2) return NULL;
     164           7 :     r1 = lfunrtopoles(r2);
     165             :   }
     166             :   else
     167             :   {
     168          28 :     r1 = lfunrtopoles(r1);
     169          28 :     if (r2) r1 = setunion_i(r1, lfunrtopoles(r2));
     170             :   }
     171          35 :   l = lg(r1); r = cgetg(l, t_VEC);
     172          70 :   for (j = 1; j < l; j++)
     173             :   {
     174          35 :     GEN be = gel(r1,j);
     175          35 :     GEN z1 = lfun(ldata1,be,bitprec), z2 = lfun(ldata2,be,bitprec);
     176          35 :     if (typ(z1) == t_SER && typ(z2) == t_SER)
     177             :     { /* pole of both, recompute to needed seriesprecision */
     178          28 :       long e = valp(z1) + valp(z2);
     179          28 :       GEN b = RgX_to_ser(deg1pol_shallow(gen_1, be, 0), 3-e);
     180          28 :       z1 = lfun(ldata1,b,bitprec);
     181          28 :       z2 = lfun(ldata2,b,bitprec);
     182             :     }
     183          35 :     gel(r,j) = mkvec2(be, gmul(z1, z2));
     184             :   }
     185          35 :   return r;
     186             : }
     187             : 
     188             : GEN
     189          42 : lfunmul(GEN ldata1, GEN ldata2, long bitprec)
     190             : {
     191          42 :   pari_sp ltop = avma;
     192             :   GEN k, r, N, Vga, eno, a1a2, b1b2, LD;
     193          42 :   ldata1 = lfunmisc_to_ldata_shallow(ldata1);
     194          42 :   ldata2 = lfunmisc_to_ldata_shallow(ldata2);
     195          42 :   k = ldata_get_k(ldata1);
     196          42 :   if (!gequal(ldata_get_k(ldata2),k))
     197           0 :     pari_err_OP("lfunmul [weight]",ldata1, ldata2);
     198          42 :   r = lfunmulpoles(ldata1, ldata2, bitprec);
     199          42 :   N = gmul(ldata_get_conductor(ldata1), ldata_get_conductor(ldata2));
     200          42 :   Vga = shallowconcat(ldata_get_gammavec(ldata1), ldata_get_gammavec(ldata2));
     201          42 :   Vga = sort(Vga);
     202          42 :   eno = gmul(ldata_get_rootno(ldata1), ldata_get_rootno(ldata2));
     203          42 :   a1a2 = lfunconvol(ldata_get_an(ldata1), ldata_get_an(ldata2));
     204          42 :   b1b2 = lfuncombdual(lfunconvol, ldata1, ldata2);
     205          42 :   LD = mkvecn(7, a1a2, b1b2, Vga, k, N, eno, r);
     206          42 :   if (!r) setlg(LD,7);
     207          42 :   return gerepilecopy(ltop, LD);
     208             : }
     209             : 
     210             : static GEN
     211          49 : lfundivpoles(GEN ldata1, GEN ldata2, long bitprec)
     212             : {
     213             :   long i, j, l;
     214          49 :   GEN k  = ldata_get_k(ldata1);
     215          49 :   GEN r1 = ldata_get_residue(ldata1);
     216          49 :   GEN r2 = ldata_get_residue(ldata2), r;
     217             : 
     218          49 :   if (r1 && typ(r1) != t_VEC) r1 = mkvec(mkvec2(k, r1));
     219          49 :   if (r2 && typ(r2) != t_VEC) r2 = mkvec(mkvec2(k, r2));
     220          49 :   if (!r1) return NULL;
     221          49 :   r1 = lfunrtopoles(r1);
     222          49 :   l = lg(r1); r = cgetg(l, t_VEC);
     223          98 :   for (i = j = 1; j < l; j++)
     224             :   {
     225          49 :     GEN be = gel(r1,j);
     226          49 :     GEN z = gdiv(lfun(ldata1,be,bitprec), lfun(ldata2,be,bitprec));
     227          49 :     if (valp(z) < 0) gel(r,i++) = mkvec2(be, z);
     228             :   }
     229          49 :   if (i == 1) return NULL;
     230          14 :   setlg(r, i); return r;
     231             : }
     232             : 
     233             : GEN
     234          49 : lfundiv(GEN ldata1, GEN ldata2, long bitprec)
     235             : {
     236          49 :   pari_sp ltop = avma;
     237             :   GEN k, r, N, v, v1, v2, eno, a1a2, b1b2, LD, eno2;
     238             :   long j, j1, j2, l1, l2;
     239          49 :   ldata1 = lfunmisc_to_ldata_shallow(ldata1);
     240          49 :   ldata2 = lfunmisc_to_ldata_shallow(ldata2);
     241          49 :   k = ldata_get_k(ldata1);
     242          49 :   if (!gequal(ldata_get_k(ldata2),k))
     243           0 :     pari_err_OP("lfundiv [weight]",ldata1, ldata2);
     244          49 :   r = lfundivpoles(ldata1, ldata2, bitprec);
     245          49 :   N = gdiv(ldata_get_conductor(ldata1), ldata_get_conductor(ldata2));
     246          49 :   if (typ(N) != t_INT) pari_err_OP("lfundiv [conductor]",ldata1, ldata2);
     247          49 :   a1a2 = lfunconvolinv(ldata_get_an(ldata1), ldata_get_an(ldata2));
     248          49 :   b1b2 = lfuncombdual(lfunconvolinv, ldata1, ldata2);
     249          49 :   eno2 = ldata_get_rootno(ldata2);
     250          49 :   eno = isintzero(eno2)? gen_0: gdiv(ldata_get_rootno(ldata1), eno2);
     251          49 :   v1 = shallowcopy(ldata_get_gammavec(ldata1));
     252          49 :   v2 = ldata_get_gammavec(ldata2);
     253          49 :   l1 = lg(v1); l2 = lg(v2);
     254         105 :   for (j2 = 1; j2 < l2; j2++)
     255             :   {
     256          84 :     for (j1 = 1; j1 < l1; j1++)
     257          84 :       if (gel(v1,j1) && gequal(gel(v1,j1), gel(v2,j2)))
     258             :       {
     259          56 :         gel(v1,j1) = NULL; break;
     260             :       }
     261          56 :     if (j1 == l1) pari_err_OP("lfundiv [Vga]",ldata1, ldata2);
     262             :   }
     263          49 :   v = cgetg(l1-l2+1, t_VEC);
     264         238 :   for (j1 = j = 1; j1 < l1; j1++)
     265         189 :     if (gel(v1, j1)) gel(v,j++) = gel(v1,j1);
     266             : 
     267          49 :   LD = mkvecn(7, a1a2, b1b2, v, k, N, eno, r);
     268          49 :   if (!r) setlg(LD,7);
     269          49 :   return gerepilecopy(ltop, LD);
     270             : }
     271             : 
     272             : static GEN
     273         140 : gamma_imagchi(GEN gam, GEN w)
     274             : {
     275         140 :   long i, j, k=1, l;
     276         140 :   GEN g = cgetg_copy(gam, &l);
     277         140 :   gam = shallowcopy(gam);
     278         420 :   for (i = l-1; i>=1; i--)
     279             :   {
     280         280 :     GEN al = gel(gam, i);
     281         280 :     if (al)
     282             :     {
     283         140 :       GEN N = gadd(w,gmul2n(real_i(al),1));
     284         140 :       if (gcmpgs(N,2) > 0)
     285             :       {
     286         140 :         GEN bl = gsubgs(al, 1);
     287         140 :         for (j=1; j < i; j++)
     288         140 :           if (gel(gam,j) && gequal(gel(gam,j), bl))
     289         140 :           { gel(gam,j) = NULL; break; }
     290         140 :         if (j==i) return NULL;
     291         140 :         gel(g, k++) = al;
     292         140 :         gel(g, k++) = bl;
     293           0 :       } else if (gequal0(N))
     294           0 :         gel(g, k++) = gaddgs(al, 1);
     295           0 :       else if (gequal1(N))
     296           0 :         gel(g, k++) = gsubgs(al, 1);
     297           0 :       else return NULL;
     298             :     }
     299             :   }
     300         140 :   return sort(g);
     301             : }
     302             : 
     303             : GEN
     304         511 : lfuntwist(GEN ldata1, GEN chi)
     305             : {
     306         511 :   pari_sp ltop = avma;
     307             :   GEN k, L, N, N1, N2, a, a1, a2, b, b1, b2, gam, gam1, gam2;
     308             :   GEN ldata2;
     309             :   long d1, t;
     310         511 :   ldata1 = lfunmisc_to_ldata_shallow(ldata1);
     311         511 :   ldata2 = lfunmisc_to_ldata_shallow(chi);
     312         511 :   t = ldata_get_type(ldata2);
     313         511 :   a1 = ldata_get_an(ldata1);
     314         511 :   a2 = ldata_get_an(ldata2);
     315         511 :   if (t == t_LFUN_ZETA)
     316         245 :     return gerepilecopy(ltop, ldata1);
     317         266 :   if (t != t_LFUN_CHIZ && t != t_LFUN_KRONECKER &&
     318           7 :     ( t != t_LFUN_CHIGEN || nf_get_degree(bnr_get_nf(gmael(a2,2,1))) != 1))
     319           0 :     pari_err_TYPE("lfuntwist", chi);
     320         266 :   N1 = ldata_get_conductor(ldata1);
     321         266 :   N2 = ldata_get_conductor(ldata2);
     322         266 :   if (!gequal1(gcdii(N1, N2)))
     323           0 :     pari_err_IMPL("lfuntwist (conductors not coprime)");
     324         266 :   k = ldata_get_k(ldata1);
     325         266 :   d1 = ldata_get_degree(ldata1);
     326         266 :   N = gmul(N1, gpowgs(N2, d1));
     327         266 :   gam1 = ldata_get_gammavec(ldata1);
     328         266 :   gam2 = ldata_get_gammavec(ldata2);
     329         266 :   if (gequal0(gel(gam2, 1)))
     330         126 :     gam = gam1;
     331             :   else
     332         140 :     gam = gamma_imagchi(ldata_get_gammavec(ldata1), gaddgs(k,-1));
     333         266 :   if (!gam) pari_err_IMPL("lfuntwist (gammafactors)");
     334         266 :   b1 = ldata_get_dual(ldata1);
     335         266 :   b2 = ldata_get_dual(ldata2);
     336         266 :   a = tag(mkvec2(a1, a2), t_LFUN_TWIST);
     337         266 :   if (typ(b1)==t_INT)
     338         266 :     b = signe(b1) && signe(b2) ? gen_0: gen_1;
     339             :   else
     340           0 :     b = tag(mkvec2(b1,lfunconj(a2)), t_LFUN_TWIST);
     341         266 :   L = mkvecn(6, a, b, gam, k, N, gen_0);
     342         266 :   return gerepilecopy(ltop, L);
     343             : }
     344             : 
     345             : /*****************************************************************/
     346             : /*  L-series from closure                                        */
     347             : /*****************************************************************/
     348             : static GEN
     349       81844 : localfactor(void *E, GEN p, long n)
     350             : {
     351       81844 :   GEN s = closure_callgen2((GEN)E, p, utoi(n));
     352       81844 :   return direuler_factor(s, n);
     353             : }
     354             : static GEN
     355        1855 : vecan_closure(GEN a, long L, long prec)
     356             : {
     357        1855 :   long ta = typ(a);
     358        1855 :   GEN gL, Sbad = NULL;
     359             : 
     360        1855 :   if (!L) return cgetg(1,t_VEC);
     361        1855 :   if (ta == t_VEC)
     362             :   {
     363         931 :     long l = lg(a);
     364         931 :     if (l == 1) pari_err_TYPE("vecan_closure", a);
     365         931 :     ta = typ(gel(a,1));
     366             :     /* regular vector, return it */
     367         931 :     if (ta != t_CLOSURE) return vecslice(a, 1, minss(L,l-1));
     368          77 :     if (l != 3) pari_err_TYPE("vecan_closure", a);
     369          77 :     Sbad = gel(a,2);
     370          77 :     if (typ(Sbad) != t_VEC) pari_err_TYPE("vecan_closure", a);
     371          77 :     a = gel(a,1);
     372             :   }
     373         924 :   else if (ta != t_CLOSURE) pari_err_TYPE("vecan_closure", a);
     374         987 :   push_localprec(prec);
     375         987 :   gL = stoi(L);
     376         987 :   switch(closure_arity(a))
     377             :   {
     378             :     case 2:
     379         434 :       a = direuler_bad((void*)a, localfactor, gen_2, gL,gL, Sbad);
     380         434 :       break;
     381             :     case 1:
     382         553 :       a = closure_callgen1(a, gL);
     383         553 :       if (typ(a) != t_VEC) pari_err_TYPE("vecan_closure", a);
     384         553 :       break;
     385           0 :     default: pari_err_TYPE("vecan_closure [wrong arity]", a);
     386             :       a = NULL; /*LCOV_EXCL_LINE*/
     387             :   }
     388         987 :   pop_localprec(); return a;
     389             : }
     390             : 
     391             : /*****************************************************************/
     392             : /*  L-series of Dirichlet characters.                            */
     393             : /*****************************************************************/
     394             : 
     395             : static GEN
     396        1316 : lfunzeta(void)
     397             : {
     398        1316 :   GEN zet = mkvecn(7, NULL, gen_0, NULL, gen_1, gen_1, gen_1, gen_1);
     399        1316 :   gel(zet,1) = tag(gen_1, t_LFUN_ZETA);
     400        1316 :   gel(zet,3) = mkvec(gen_0);
     401        1316 :   return zet;
     402             : }
     403             : static GEN
     404         175 : lfunzetainit(GEN dom, long der, long bitprec)
     405         175 : { return lfuninit(lfunzeta(), dom, der, bitprec); }
     406             : 
     407             : static GEN
     408         980 : vecan_Kronecker(GEN D, long n)
     409             : {
     410         980 :   GEN v = cgetg(n+1, t_VECSMALL);
     411         980 :   ulong Du = itou_or_0(D);
     412         980 :   long i, id, d = Du ? minuu(Du, n): n;
     413         980 :   for (i = 1; i <= d; i++) v[i] = krois(D,i);
     414       90167 :   for (id = i; i <= n; i++,id++) /* periodic mod d */
     415             :   {
     416       89187 :     if (id > d) id = 1;
     417       89187 :     gel(v, i) = gel(v, id);
     418             :   }
     419         980 :   return v;
     420             : }
     421             : 
     422             : static GEN
     423        1106 : lfunchiquad(GEN D)
     424             : {
     425             :   GEN r;
     426        1106 :   if (equali1(D)) return lfunzeta();
     427         749 :   if (!isfundamental(D)) pari_err_TYPE("lfunchiquad [not primitive]", D);
     428         749 :   r = mkvecn(6, NULL, gen_0, NULL, gen_1, NULL, gen_1);
     429         749 :   gel(r,1) = tag(icopy(D), t_LFUN_KRONECKER);
     430         749 :   gel(r,3) = mkvec(signe(D) < 0? gen_1: gen_0);
     431         749 :   gel(r,5) = mpabs(D);
     432         749 :   return r;
     433             : }
     434             : 
     435             : /* Begin Hecke characters. Here a character is assumed to be given by a
     436             :    vector on the generators of the ray class group clgp of CL_m(K).
     437             :    If clgp = [h,[d1,...,dk],[g1,...,gk]] with dk|...|d2|d1, a character chi
     438             :    is given by [a1,a2,...,ak] such that chi(gi)=\zeta_di^ai. */
     439             : 
     440             : /* Value of CHI on x, coprime to bnr.mod */
     441             : static GEN
     442       27188 : chigeneval(GEN logx, GEN nchi, GEN z, long prec)
     443             : {
     444       27188 :   pari_sp av = avma;
     445       27188 :   GEN d = gel(nchi,1);
     446       27188 :   GEN e = FpV_dotproduct(gel(nchi,2), logx, d);
     447       27188 :   if (!is_vec_t(typ(z)))
     448           0 :     return gerepileupto(av, gpow(z, e, prec));
     449             :   else
     450             :   {
     451       27188 :     ulong i = itou(e);
     452       27188 :     set_avma(av); return gel(z, i+1);
     453             :   }
     454             : }
     455             : 
     456             : /* return x + yz; y != 0; z = 0,1 "often"; x = 0 "often" */
     457             : static GEN
     458      217973 : gaddmul(GEN x, GEN y, GEN z)
     459             : {
     460             :   pari_sp av;
     461      217973 :   if (typ(z) == t_INT)
     462             :   {
     463      206752 :     if (!signe(z)) return x;
     464       19145 :     if (equali1(z)) return gadd(x,y);
     465             :   }
     466       25452 :   if (isintzero(x)) return gmul(y,z);
     467        7511 :   av = avma;
     468        7511 :   return gerepileupto(av, gadd(x, gmul(y,z)));
     469             : }
     470             : 
     471             : static GEN
     472         539 : vecan_chiZ(GEN an, long n, long prec)
     473             : {
     474             :   forprime_t iter;
     475         539 :   GEN G = gel(an,1);
     476         539 :   GEN nchi = gel(an,2), gord = gel(nchi,1), z;
     477         539 :   GEN gp = cgetipos(3), v = vec_ei(n, 1);
     478         539 :   GEN N = znstar_get_N(G);
     479         539 :   long ord = itos_or_0(gord);
     480         539 :   ulong Nu = itou_or_0(N);
     481         539 :   long i, id, d = Nu ? minuu(Nu, n): n;
     482             :   ulong p;
     483             : 
     484         539 :   if (ord && n > (ord>>4))
     485         539 :   {
     486         539 :     GEN w = ncharvecexpo(G, nchi);
     487         539 :     z = grootsof1(ord, prec);
     488       12376 :     for (i = 1; i <= d; i++)
     489       11837 :       if (w[i] >= 0) gel(v, i) = gel(z, w[i]+1);
     490             :   }
     491             :   else
     492             :   {
     493           0 :     z = rootsof1_cx(gord, prec);
     494           0 :     u_forprime_init(&iter, 2, d);
     495           0 :     while ((p = u_forprime_next(&iter)))
     496             :     {
     497             :       GEN ch;
     498             :       ulong k;
     499           0 :       if (!umodiu(N,p)) continue;
     500           0 :       gp[2] = p;
     501           0 :       ch = chigeneval(znconreylog(G, gp), nchi, z, prec);
     502           0 :       gel(v, p)  = ch;
     503           0 :       for (k = 2*p; k <= (ulong)d; k += p)
     504           0 :         gel(v, k) = gaddmul(gel(v, k), ch, gel(v, k/p));
     505             :     }
     506             :   }
     507      155526 :   for (id = i = d+1; i <= n; i++,id++) /* periodic mod d */
     508             :   {
     509      154987 :     if (id > d) id = 1;
     510      154987 :     gel(v, i) = gel(v, id);
     511             :   }
     512         539 :   return v;
     513             : }
     514             : 
     515             : static GEN
     516         672 : vecan_chigen(GEN an, long n, long prec)
     517             : {
     518             :   forprime_t iter;
     519         672 :   GEN bnr = gel(an,1), nf = bnr_get_nf(bnr);
     520         672 :   GEN nchi = gel(an,2), gord = gel(nchi,1), z;
     521         672 :   GEN gp = cgetipos(3), v = vec_ei(n, 1);
     522         672 :   GEN N = gel(bnr_get_mod(bnr), 1), NZ = gcoeff(N,1,1);
     523         672 :   long ord = itos_or_0(gord);
     524             :   ulong p;
     525             : 
     526         672 :   if (ord && n > (ord>>4))
     527         672 :     z = grootsof1(ord, prec);
     528             :   else
     529           0 :     z = rootsof1_cx(gord, prec);
     530             : 
     531         672 :   if (nf_get_degree(nf) == 1)
     532             :   {
     533         525 :     ulong Nu = itou_or_0(NZ);
     534         525 :     long i, id, d = Nu ? minuu(Nu, n): n;
     535         525 :     u_forprime_init(&iter, 2, d);
     536        3129 :     while ((p = u_forprime_next(&iter)))
     537             :     {
     538             :       GEN ch;
     539             :       ulong k;
     540        2079 :       if (!umodiu(NZ,p)) continue;
     541        1554 :       gp[2] = p;
     542        1554 :       ch = chigeneval(isprincipalray(bnr,gp), nchi, z, prec);
     543        1554 :       gel(v, p)  = ch;
     544        3626 :       for (k = 2*p; k <= (ulong)d; k += p)
     545        2072 :         gel(v, k) = gaddmul(gel(v, k), ch, gel(v, k/p));
     546             :     }
     547        7098 :     for (id = i = d+1; i <= n; i++,id++) /* periodic mod d */
     548             :     {
     549        6573 :       if (id > d) id = 1;
     550        6573 :       gel(v, i) = gel(v, id);
     551             :     }
     552             :   }
     553             :   else
     554             :   {
     555         147 :     GEN BOUND = stoi(n);
     556         147 :     u_forprime_init(&iter, 2, n);
     557       26355 :     while ((p = u_forprime_next(&iter)))
     558             :     {
     559             :       GEN L;
     560             :       long j;
     561       26061 :       int check = !umodiu(NZ,p);
     562       26061 :       gp[2] = p;
     563       26061 :       L = idealprimedec_limit_norm(nf, gp, BOUND);
     564       51793 :       for (j = 1; j < lg(L); j++)
     565             :       {
     566       25732 :         GEN pr = gel(L, j), ch;
     567             :         ulong k, q;
     568       25732 :         if (check && idealval(nf, N, pr)) continue;
     569       25634 :         ch = chigeneval(isprincipalray(bnr,pr), nchi, z, prec);
     570       25634 :         q = upr_norm(pr);
     571       25634 :         gel(v, q) = gadd(gel(v, q), ch);
     572      241535 :         for (k = 2*q; k <= (ulong)n; k += q)
     573      215901 :           gel(v, k) = gaddmul(gel(v, k), ch, gel(v, k/q));
     574             :       }
     575             :     }
     576             :   }
     577         672 :   return v;
     578             : }
     579             : 
     580             : static GEN lfunzetak_i(GEN T);
     581             : static GEN
     582        3626 : vec01(long r1, long r2)
     583             : {
     584        3626 :   long d = r1+r2, i;
     585        3626 :   GEN v = cgetg(d+1,t_VEC);
     586        3626 :   for (i = 1; i <= r1; i++) gel(v,i) = gen_0;
     587        3626 :   for (     ; i <= d;  i++) gel(v,i) = gen_1;
     588        3626 :   return v;
     589             : }
     590             : 
     591             : /* G is a bid of nftyp typ_BIDZ */
     592             : static GEN
     593         994 : lfunchiZ(GEN G, GEN chi)
     594             : {
     595         994 :   pari_sp av = avma;
     596         994 :   GEN sig = NULL;
     597         994 :   GEN N = bid_get_ideal(G), nchi, r;
     598             :   int real;
     599             : 
     600         994 :   if (typ(N) != t_INT) pari_err_TYPE("lfunchiZ", G);
     601         994 :   if (equali1(N)) return lfunzeta();
     602         749 :   if (typ(chi) != t_COL) chi = znconreylog(G,chi);
     603         749 :   N = znconreyconductor(G, chi, &chi);
     604         749 :   if (typ(N) != t_INT)
     605             :   {
     606           0 :     if (equali1(gel(N,1))) { set_avma(av); return lfunzeta(); }
     607           0 :     G = znstar0(N, 1);
     608           0 :     N = gel(N,1);
     609             :   }
     610             :   /* chi now primitive on G */
     611         749 :   switch(itou_or_0(zncharorder(G, chi)))
     612             :   {
     613           0 :     case 1: set_avma(av); return lfunzeta();
     614         420 :     case 2: if (zncharisodd(G,chi)) N = negi(N);
     615         420 :             return gerepileupto(av, lfunchiquad(N));
     616             :   }
     617         329 :   sig = mkvec( zncharisodd(G, chi)? gen_1: gen_0 );
     618         329 :   nchi = znconreylog_normalize(G, chi);
     619         329 :   real = abscmpiu(gel(nchi,1), 2) <= 0;
     620         329 :   r = mkvecn(6, tag(mkvec2(G,nchi), t_LFUN_CHIZ),
     621             :                 real? gen_0: gen_1, sig, gen_1, N, gen_0);
     622         329 :   return gerepilecopy(av, r);
     623             : }
     624             : 
     625             : static GEN
     626        1841 : lfunchigen(GEN bnr, GEN CHI)
     627             : {
     628        1841 :   pari_sp av = avma;
     629             :   GEN v;
     630             :   GEN N, sig, Ldchi, nf, nchi, NN;
     631             :   long r1, r2, n1;
     632             :   int real;
     633             : 
     634        1841 :   if (nftyp(bnr) == typ_BIDZ) return lfunchiZ(bnr, CHI);
     635             : 
     636         868 :   v = bnrconductor_i(bnr, CHI, 2);
     637         868 :   bnr = gel(v,2);
     638         868 :   CHI = gel(v,3); /* now CHI is primitive wrt bnr */
     639             : 
     640         868 :   nf = bnr_get_nf(bnr);
     641         868 :   N = bnr_get_mod(bnr);
     642         868 :   n1 = lg(vec01_to_indices(gel(N,2))) - 1; /* vecsum(N[2]) */
     643         868 :   N = gel(N,1);
     644         868 :   NN = mulii(idealnorm(nf, N), absi_shallow(nf_get_disc(nf)));
     645         868 :   if (equali1(NN)) return gerepileupto(av, lfunzeta());
     646         546 :   if (ZV_equal0(CHI)) return gerepilecopy(av, lfunzetak_i(bnr));
     647         539 :   nf_get_sign(nf, &r1, &r2);
     648         539 :   sig = vec01(r1+r2-n1, r2+n1);
     649         539 :   nchi = char_normalize(CHI, cyc_normalize(bnr_get_cyc(bnr)));
     650         539 :   real = abscmpiu(gel(nchi,1), 2) <= 0;
     651         539 :   Ldchi = mkvecn(6, tag(mkvec2(bnr, nchi), t_LFUN_CHIGEN),
     652             :                     real? gen_0: gen_1, sig, gen_1, NN, gen_0);
     653         539 :   return gerepilecopy(av, Ldchi);
     654             : }
     655             : 
     656             : /* Find all characters of clgp whose kernel contain group given by HNF H.
     657             :  * Set *pcnj[i] if chi[i] is not real */
     658             : static GEN
     659         329 : chigenkerfind(GEN bnr, GEN H, GEN *pcnj)
     660             : {
     661         329 :   GEN res, cnj, L = bnrchar(bnr, H, NULL), cyc = bnr_get_cyc(bnr);
     662         329 :   long i, k, l = lg(L);
     663             : 
     664         329 :   res = cgetg(l, t_VEC);
     665         329 :   *pcnj = cnj = cgetg(l, t_VECSMALL);
     666        1267 :   for (i = k = 1; i < l; i++)
     667             :   {
     668         938 :     GEN chi = gel(L,i), c = charconj(cyc, chi);
     669         938 :     long fl = ZV_cmp(c, chi);
     670         938 :     if (fl < 0) continue; /* keep one char in pair of conjugates */
     671         763 :     gel(res, k) = chi;
     672         763 :     cnj[k] = fl; k++;
     673             :   }
     674         329 :   setlg(cnj, k);
     675         329 :   setlg(res, k); return res;
     676             : }
     677             : 
     678             : static GEN
     679        1204 : lfunzetak_i(GEN T)
     680             : {
     681        1204 :   GEN Vga, N, nf, bnf = checkbnf_i(T), r = gen_0/*unknown*/;
     682             :   long r1, r2;
     683             : 
     684        1204 :   if (bnf)
     685          98 :     nf = bnf_get_nf(bnf);
     686             :   else
     687             :   {
     688        1106 :     nf = checknf_i(T);
     689        1106 :     if (!nf) nf = T = nfinit(T, DEFAULTPREC);
     690             :   }
     691        1204 :   nf_get_sign(nf,&r1,&r2);
     692        1204 :   N = absi_shallow(nf_get_disc(nf));
     693        1204 :   if (bnf)
     694             :   {
     695          98 :     GEN h = bnf_get_no(bnf);
     696          98 :     GEN R = bnf_get_reg(bnf);
     697          98 :     long prec = nf_get_prec(nf);
     698         196 :     r = gdiv(gmul(mulir(shifti(h, r1+r2), powru(mppi(prec), r2)), R),
     699          98 :              mulur(bnf_get_tuN(bnf), gsqrt(N, prec)));
     700             :   }
     701        1204 :   Vga = vec01(r1+r2,r2);
     702        1204 :   return mkvecn(7, tag(T,t_LFUN_NF), gen_0, Vga, gen_1, N, gen_1, r);
     703             : }
     704             : static GEN
     705         511 : lfunzetak(GEN T)
     706         511 : { pari_sp ltop = avma; return gerepilecopy(ltop, lfunzetak_i(T)); }
     707             : 
     708             : /* bnf = NULL: base field = Q */
     709             : GEN
     710         329 : lfunabelianrelinit(GEN nfabs, GEN bnf, GEN polrel, GEN dom, long der, long bitprec)
     711             : {
     712         329 :   pari_sp ltop = avma;
     713             :   GEN cond, chi, cnj, res, bnr, M, domain;
     714             :   long l, i;
     715         329 :   long v = -1;
     716             : 
     717         329 :   if (bnf) bnf = checkbnf(bnf);
     718             :   else
     719             :   {
     720         322 :     v = fetch_var();
     721         322 :     bnf = Buchall(pol_x(v), 0, nbits2prec(bitprec));
     722             :   }
     723         329 :   if (typ(polrel) != t_POL) pari_err_TYPE("lfunabelianrelinit", polrel);
     724         329 :   cond = rnfconductor(bnf, polrel);
     725         329 :   chi = chigenkerfind(gel(cond,2), gel(cond,3), &cnj);
     726         329 :   bnr = Buchray(bnf, gel(cond,1), nf_INIT);
     727         329 :   l = lg(chi); res = cgetg(l, t_VEC);
     728        1092 :   for (i = 1; i < l; ++i)
     729             :   {
     730         763 :     GEN L = lfunchigen(bnr, gel(chi,i));
     731         763 :     gel(res, i) = lfuninit(L, dom, der, bitprec);
     732             :   }
     733         329 :   if (v >= 0) delete_var();
     734         329 :   M = mkvec3(res, const_vecsmall(l-1, 1), cnj);
     735         329 :   domain = mkvec2(dom, mkvecsmall2(der, bitprec));
     736         329 :   return gerepilecopy(ltop, lfuninit_make(t_LDESC_PRODUCT, lfunzetak_i(nfabs), M, domain));
     737             : }
     738             : 
     739             : /*****************************************************************/
     740             : /*                 Dedekind zeta functions                       */
     741             : /*****************************************************************/
     742             : static GEN
     743        1267 : dirzetak0(GEN nf, ulong N)
     744             : {
     745        1267 :   GEN vect, c, c2, T = nf_get_pol(nf), index = nf_get_index(nf);
     746        1267 :   pari_sp av = avma, av2;
     747        1267 :   const ulong SQRTN = usqrt(N);
     748             :   ulong i, p, lx;
     749        1267 :   long court[] = {evaltyp(t_INT)|_evallg(3), evalsigne(1)|evallgefint(3),0};
     750             :   forprime_t S;
     751             : 
     752        1267 :   c  = cgetalloc(t_VECSMALL, N+1);
     753        1267 :   c2 = cgetalloc(t_VECSMALL, N+1);
     754        1267 :   c2[1] = c[1] = 1; for (i=2; i<=N; i++) c[i] = 0;
     755        1267 :   u_forprime_init(&S, 2, N); av2 = avma;
     756      163737 :   while ( (p = u_forprime_next(&S)) )
     757             :   {
     758      161203 :     set_avma(av2);
     759      161203 :     if (umodiu(index, p)) /* p does not divide index */
     760      161000 :       vect = gel(Flx_degfact(ZX_to_Flx(T,p), p),1);
     761             :     else
     762             :     {
     763         203 :       court[2] = p;
     764         203 :       vect = idealprimedec_degrees(nf,court);
     765             :     }
     766      161203 :     lx = lg(vect);
     767      161203 :     if (p <= SQRTN)
     768       19740 :       for (i=1; i<lx; i++)
     769             :       {
     770       13111 :         ulong qn, q = upowuu(p, vect[i]); /* Norm P[i] */
     771       13111 :         if (!q || q > N) break;
     772       11249 :         memcpy(c2 + 2, c + 2, (N-1)*sizeof(long));
     773             :         /* c2[i] <- c[i] + sum_{k = 1}^{v_q(i)} c[i/q^k] for all i <= N */
     774       23205 :         for (qn = q; qn <= N; qn *= q)
     775             :         {
     776       23205 :           ulong k0 = N/qn, k, k2; /* k2 = k*qn */
     777       23205 :           for (k = k0, k2 = k*qn; k > 0; k--, k2 -=qn) c2[k2] += c[k];
     778       23205 :           if (q > k0) break; /* <=> q*qn > N */
     779             :         }
     780       11249 :         swap(c, c2);
     781             :       }
     782             :     else /* p > sqrt(N): simpler */
     783      301252 :       for (i=1; i<lx; i++)
     784             :       {
     785             :         ulong k, k2; /* k2 = k*p */
     786      263963 :         if (vect[i] > 1) break;
     787             :         /* c2[i] <- c[i] + sum_{k = 1}^{v_q(i)} c[i/q^k] for all i <= N */
     788      148540 :         for (k = N/p, k2 = k*p; k > 0; k--, k2 -= p) c[k2] += c[k];
     789             :       }
     790             :   }
     791        1267 :   set_avma(av);
     792        1267 :   pari_free(c2); return c;
     793             : }
     794             : 
     795             : GEN
     796        1267 : dirzetak(GEN nf, GEN b)
     797             : {
     798             :   GEN z, c;
     799             :   long n;
     800             : 
     801        1267 :   if (typ(b) != t_INT) pari_err_TYPE("dirzetak",b);
     802        1267 :   if (signe(b) <= 0) return cgetg(1,t_VEC);
     803        1267 :   nf = checknf(nf);
     804        1267 :   n = itou_or_0(b); if (!n) pari_err_OVERFLOW("dirzetak");
     805        1267 :   c = dirzetak0(nf, n);
     806        1267 :   z = vecsmall_to_vec(c); pari_free(c); return z;
     807             : }
     808             : 
     809             : static GEN
     810         672 : linit_get_mat(GEN linit)
     811             : {
     812         672 :   if (linit_get_type(linit)==t_LDESC_PRODUCT)
     813         161 :     return lfunprod_get_fact(linit_get_tech(linit));
     814             :   else
     815         511 :     return mkvec3(mkvec(linit), mkvecsmall(1), mkvecsmall(0));
     816             : }
     817             : 
     818             : static GEN
     819         336 : lfunproduct(GEN ldata, GEN linit1, GEN linit2, GEN domain)
     820             : {
     821         336 :   GEN M1 = linit_get_mat(linit1);
     822         336 :   GEN M2 = linit_get_mat(linit2);
     823        1344 :   GEN M3 = mkvec3(shallowconcat(gel(M1, 1), gel(M2, 1)),
     824         672 :                   vecsmall_concat(gel(M1, 2), gel(M2, 2)),
     825         672 :                   vecsmall_concat(gel(M1, 3), gel(M2, 3)));
     826         336 :   return lfuninit_make(t_LDESC_PRODUCT, ldata, M3, domain);
     827             : }
     828             : 
     829             : static GEN
     830           0 : lfunzetakinit_raw(GEN T, GEN dom, long der, long bitprec)
     831             : {
     832           0 :   pari_sp ltop = avma;
     833           0 :   GEN ldata = lfunzetak_i(T);
     834           0 :   return gerepileupto(ltop, lfuninit(ldata, dom, der, bitprec));
     835             : }
     836             : 
     837             : static GEN
     838         336 : lfunzetakinit_quotient(GEN nf, GEN polk, GEN dom, long der, long bitprec)
     839             : {
     840         336 :   pari_sp av = avma;
     841             :   GEN ak, an, nfk, Vga, ldata, N, Lk, LKk, domain;
     842             :   long r1k, r2k, r1, r2;
     843             : 
     844         336 :   nf_get_sign(nf,&r1,&r2);
     845         336 :   nfk = nfinit(polk, nbits2prec(bitprec));
     846         336 :   Lk = lfunzetakinit(nfk, dom, der, 0, bitprec); /* zeta_k */
     847         336 :   nf_get_sign(nfk,&r1k,&r2k);
     848         336 :   Vga = vec01((r1+r2) - (r1k+r2k), r2-r2k);
     849         336 :   N = absi_shallow(diviiexact(nf_get_disc(nf), nf_get_disc(nfk)));
     850         336 :   ak = nf_get_degree(nf)==1 ? tag(gen_1, t_LFUN_ZETA): tag(nfk, t_LFUN_NF);
     851         336 :   an = tag(mkvec2(tag(nf,t_LFUN_NF), ak), t_LFUN_DIV);
     852         336 :   ldata = mkvecn(6, an, gen_0, Vga, gen_1, N, gen_1);
     853         336 :   LKk = lfuninit(ldata, dom, der, bitprec); /* zeta_K/zeta_k */
     854         336 :   domain = mkvec2(dom, mkvecsmall2(der, bitprec));
     855         336 :   return gerepilecopy(av, lfunproduct(lfunzetak_i(nf), Lk, LKk, domain));
     856             : }
     857             : 
     858             : static GEN
     859             : lfunzetakinit_artin(GEN nf, GEN gal, GEN dom, long der, long bitprec);
     860             : 
     861             : static GEN
     862         343 : lfunzetakinit_Galois(GEN nf, GEN gal, GEN dom, long der, long bitprec)
     863             : {
     864         343 :   GEN grp = galois_group(gal);
     865         343 :   if (group_isabelian(grp))
     866         322 :     return lfunabelianrelinit(nf, NULL, gal_get_pol(gal), dom, der, bitprec);
     867          21 :   else return lfunzetakinit_artin(nf, gal, dom, der, bitprec);
     868             : }
     869             : 
     870             : GEN
     871         854 : lfunzetakinit(GEN NF, GEN dom, long der, long flag, long bitprec)
     872             : {
     873         854 :   GEN nf = checknf(NF);
     874             :   GEN G, nfs, sbg;
     875         854 :   long lf, d = nf_get_degree(nf);
     876         854 :   if (d == 1) return lfunzetainit(dom, der, bitprec);
     877         679 :   G = galoisinit(nf, NULL);
     878         679 :   if (!isintzero(G))
     879         343 :     return lfunzetakinit_Galois(nf, G, dom, der, bitprec);
     880         336 :   nfs = nfsubfields(nf, 0); lf = lg(nfs)-1;
     881         336 :   sbg = gmael(nfs,lf-1,1);
     882         336 :   if (flag && d > 4*degpol(sbg))
     883           0 :     return lfunzetakinit_raw(nf, dom, der, bitprec);
     884         336 :   return lfunzetakinit_quotient(nf, sbg, dom, der, bitprec);
     885             : }
     886             : 
     887             : /***************************************************************/
     888             : /*             Elliptic Curves and Modular Forms               */
     889             : /***************************************************************/
     890             : 
     891             : static GEN
     892         168 : lfunellnf(GEN e)
     893             : {
     894         168 :   pari_sp av = avma;
     895         168 :   GEN ldata = cgetg(7, t_VEC), nf = ellnf_get_nf(e);
     896         168 :   GEN N = gel(ellglobalred(e), 1);
     897         168 :   long n = nf_get_degree(nf);
     898         168 :   gel(ldata, 1) = tag(e, t_LFUN_ELL);
     899         168 :   gel(ldata, 2) = gen_0;
     900         168 :   gel(ldata, 3) = vec01(n, n);
     901         168 :   gel(ldata, 4) = gen_2;
     902         168 :   gel(ldata, 5) = mulii(idealnorm(nf,N), sqri(nf_get_disc(nf)));
     903         168 :   gel(ldata, 6) = stoi(ellrootno_global(e));
     904         168 :   return gerepilecopy(av, ldata);
     905             : }
     906             : 
     907             : static GEN
     908        1092 : lfunellQ(GEN e)
     909             : {
     910        1092 :   pari_sp av = avma;
     911        1092 :   GEN ldata = cgetg(7, t_VEC);
     912        1092 :   gel(ldata, 1) = tag(ellanal_globalred(e, NULL), t_LFUN_ELL);
     913        1092 :   gel(ldata, 2) = gen_0;
     914        1092 :   gel(ldata, 3) = mkvec2(gen_0, gen_1);
     915        1092 :   gel(ldata, 4) = gen_2;
     916        1092 :   gel(ldata, 5) = ellQ_get_N(e);
     917        1092 :   gel(ldata, 6) = stoi(ellrootno_global(e));
     918        1092 :   return gerepilecopy(av, ldata); /* ellanal_globalred not gerepile-safe */
     919             : }
     920             : 
     921             : static GEN
     922        1260 : lfunell(GEN e)
     923             : {
     924        1260 :   long t = ell_get_type(e);
     925        1260 :   switch(t)
     926             :   {
     927        1092 :     case t_ELL_Q: return lfunellQ(e);
     928         168 :     case t_ELL_NF:return lfunellnf(e);
     929             :   }
     930           0 :   pari_err_TYPE("lfun",e);
     931             :   return NULL; /*LCOV_EXCL_LINE*/
     932             : }
     933             : 
     934             : static GEN
     935         140 : ellsympow_gamma(long m)
     936             : {
     937         140 :   GEN V = cgetg(m+2, t_VEC);
     938         140 :   long i = 1, j;
     939         140 :   if (!odd(m)) gel(V, i++) = stoi(-2*(m>>2));
     940         364 :   for (j = (m+1)>>1; j > 0; i+=2, j--)
     941             :   {
     942         224 :     gel(V,i)   = stoi(1-j);
     943         224 :     gel(V,i+1) = stoi(1-j+1);
     944             :   }
     945         140 :   return V;
     946             : }
     947             : 
     948             : static GEN
     949       72543 : ellsympow_trace(GEN p, GEN t, long m)
     950             : {
     951       72543 :   long k, n = m >> 1;
     952       72543 :   GEN tp = gpowers0(sqri(t), n, odd(m)? t: NULL);
     953       72541 :   GEN pp = gen_1, b = gen_1, r = gel(tp,n+1);
     954      204293 :   for(k=1; k<=n; k++)
     955             :   {
     956             :     GEN s;
     957      131948 :     pp = mulii(pp, p);
     958      131402 :     b  = diviuexact(muliu(b, (m-(2*k-1))*(m-(2*k-2))), k*(m-(k-1)));
     959      131193 :     s = mulii(mulii(b, gel(tp,1+n-k)), pp);
     960      131514 :     r = odd(k) ? subii(r, s): addii(r, s);
     961             :   }
     962       72345 :   return r;
     963             : }
     964             : 
     965             : static GEN
     966        2765 : ellsympow_abelian(GEN p, GEN ap, long m, long o)
     967             : {
     968        2765 :   pari_sp av = avma;
     969        2765 :   long i, M, n = (m+1)>>1;
     970             :   GEN pk, tv, pn, pm, F, v;
     971        2765 :   if (!odd(o))
     972             :   {
     973           0 :     if (odd(m)) return pol_1(0);
     974           0 :     M = m >> 1; o >>= 1;
     975             :   }
     976             :   else
     977        2765 :     M = m * ((o+1) >> 1);
     978        2765 :   pk = gpowers(p,n); pn = gel(pk,n+1);
     979        2765 :   tv = cgetg(m+2,t_VEC);
     980        2765 :   gel(tv, 1) = gen_2;
     981        2765 :   gel(tv, 2) = ap;
     982        9534 :   for (i = 3; i <= m+1; i++)
     983        6769 :     gel(tv,i) = subii(mulii(ap,gel(tv,i-1)), mulii(p,gel(tv,i-2)));
     984        2765 :   pm = odd(m)? mulii(gel(pk,n), pn): sqri(pn); /* cheap p^m */
     985        2765 :   F = deg2pol_shallow(pm, gen_0, gen_1, 0);
     986        2765 :   v = odd(m) ? pol_1(0): deg1pol_shallow(negi(pn), gen_1, 0);
     987        8393 :   for (i = M % o; i < n; i += o) /* o | m-2*i */
     988             :   {
     989        5628 :     gel(F,3) = negi(mulii(gel(tv,m-2*i+1), gel(pk,i+1)));
     990        5628 :     v = ZX_mul(v, F);
     991             :   }
     992        2765 :   return gerepilecopy(av, v);
     993             : }
     994             : 
     995             : static GEN
     996       75326 : ellsympow(GEN E, ulong m, GEN p, long n)
     997             : {
     998       75326 :   pari_sp av = avma;
     999       75326 :   GEN ap = ellap(E, p);
    1000       75305 :   if (n <= 2)
    1001             :   {
    1002       72540 :     GEN t = ellsympow_trace(p, ap, m);
    1003       72343 :     return deg1pol_shallow(t, gen_1, 0);
    1004             :   }
    1005             :   else
    1006        2765 :     return gerepileupto(av, RgXn_inv_i(ellsympow_abelian(p, ap, m, 1), n));
    1007             : }
    1008             : 
    1009             : GEN
    1010        4970 : direllsympow_worker(GEN P, ulong X, GEN E, ulong m)
    1011             : {
    1012        4970 :   pari_sp av = avma;
    1013        4970 :   long i, l = lg(P);
    1014        4970 :   GEN W = cgetg(l, t_VEC);
    1015       80305 :   for(i = 1; i < l; i++)
    1016             :   {
    1017       75329 :     ulong p = uel(P,i);
    1018       75329 :     long d = ulogint(X, p) + 1; /* minimal d such that p^d > X */
    1019       75326 :     gel(W,i) = ellsympow(E, m, utoi(uel(P,i)), d);
    1020             :   }
    1021        4976 :   return gerepilecopy(av, mkvec2(P,W));
    1022             : }
    1023             : 
    1024             : static GEN
    1025         315 : vecan_ellsympow(GEN an, long n)
    1026             : {
    1027         315 :   GEN nn = utoi(n), crvm = gel(an,1), bad = gel(an,2);
    1028         315 :   GEN worker = strtoclosure("_direllsympow_worker",2, gel(crvm,1), gel(crvm,2));
    1029         315 :   return pardireuler(worker, gen_2, nn, nn, bad);
    1030             : }
    1031             : 
    1032             : static long
    1033         196 : ellsympow_betam(long o, long m)
    1034             : {
    1035         196 :   const long c3[]={3, -1, 1};
    1036         196 :   const long c12[]={6, -2, 2, 0, 4, -4};
    1037         196 :   const long c24[]={12, -2, -4, 6, 4, -10};
    1038         196 :   if (!odd(o) && odd(m)) return 0;
    1039         161 :   switch(o)
    1040             :   {
    1041           0 :     case 1:  return m+1;
    1042          14 :     case 2:  return m+1;
    1043          84 :     case 3:  case 6: return (m+c3[m%3])/3;
    1044           0 :     case 4:  return m%4 == 0 ? (m+2)/2: m/2;
    1045          21 :     case 8:  return m%4 == 0 ? (m+4)/4: (m-2)/4;
    1046          35 :     case 12: return (m+c12[(m%12)/2])/6;
    1047           7 :     case 24: return (m+c24[(m%12)/2])/12;
    1048             :   }
    1049           0 :   return 0;
    1050             : }
    1051             : 
    1052             : static long
    1053          98 : ellsympow_epsm(long o, long m)
    1054             : {
    1055          98 :   return m+1-ellsympow_betam(o, m);
    1056             : }
    1057             : 
    1058             : static GEN
    1059          98 : ellsympow_multred(GEN E, GEN p, long m, long vN, long *cnd, long *w)
    1060             : {
    1061          98 :   if (vN==1) /* mult red*/
    1062             :   {
    1063          98 :     *cnd = m;
    1064          98 :     *w = odd(m) ? ellrootno(E, p): 1;
    1065          98 :     if (odd(m) && signe(ellap(E,p))==-1)
    1066           0 :       return deg1pol_shallow(gen_1, gen_1, 0);
    1067             :     else
    1068          98 :       return deg1pol_shallow(gen_m1, gen_1, 0);
    1069             :   } else
    1070             :   {
    1071           0 :     *cnd = odd(m)? m+1: m;
    1072           0 :     *w  = odd(m) && odd((m+1)>>1) ? ellrootno(E, p): 1;
    1073           0 :     if (odd(m) && equaliu(p,2))
    1074           0 :       *cnd += ((m+1)>>1)*(vN-2);
    1075           0 :     return odd(m)? pol_1(0): deg1pol_shallow(gen_m1, gen_1, 0);
    1076             :   }
    1077             : }
    1078             : 
    1079             : static GEN
    1080          98 : ellsympow_nonabelian(GEN p, long m, long bet)
    1081             : {
    1082          98 :  GEN pm = powiu(p, m), F;
    1083          98 :  if (odd(m)) return gpowgs(deg2pol_shallow(pm, gen_0, gen_1, 0), bet>>1);
    1084          63 :  F = gpowgs(deg2pol_shallow(negi(pm), gen_0, gen_1, 0), bet>>1);
    1085          63 :  if (!odd(bet)) return F;
    1086          28 :  if (m%4==2)
    1087          21 :    return gmul(F, deg1pol_shallow(powiu(p, m>>1), gen_1, 0));
    1088             :  else
    1089           7 :    return gmul(F, deg1pol_shallow(negi(powiu(p, m>>1)), gen_1, 0));
    1090             : }
    1091             : 
    1092             : static long
    1093           0 : safe_Z_pvalrem(GEN n, GEN p, GEN *pr)
    1094           0 : { return signe(n)==0? -1: Z_pvalrem(n, p, pr); }
    1095             : 
    1096             : static GEN
    1097           0 : c4c6_ap(GEN c4, GEN c6, GEN p)
    1098             : {
    1099           0 :   GEN N = Fp_ellcard(Fp_muls(c4,-27,p), Fp_muls(c6, -54, p), p);
    1100           0 :   return subii(addis(p, 1), N);
    1101             : }
    1102             : 
    1103             : static GEN
    1104           0 : ellsympow_abelian_twist(GEN E, GEN p, long m, long o)
    1105             : {
    1106             :   GEN ap;
    1107           0 :   GEN c4 = ell_get_c4(E), c6 = ell_get_c6(E);
    1108             :   GEN c4t, c6t;
    1109           0 :   long v4 = safe_Z_pvalrem(c4, p, &c4t);
    1110           0 :   long v6 = safe_Z_pvalrem(c6, p, &c6t);
    1111           0 :   if (v6>=0 && (v4==-1 || 3*v4>=2*v6)) c6 = c6t;
    1112           0 :   if (v4>=0 && (v6==-1 || 3*v4<=2*v6)) c4 = c4t;
    1113           0 :   ap = c4c6_ap(c4, c6, p);
    1114           0 :   return ellsympow_abelian(p, ap, m, o);
    1115             : }
    1116             : 
    1117             : static GEN
    1118           0 : ellsympow_goodred(GEN E, GEN p, long m, long *cnd, long *w)
    1119             : {
    1120           0 :   long o = 12/cgcd(12, Z_pval(ell_get_disc(E), p));
    1121           0 :   long bet = ellsympow_betam(o, m);
    1122           0 :   long eps = m + 1 - bet;
    1123           0 :   *w = odd(m) && odd(eps>>1) ? ellrootno(E,p): 1;
    1124           0 :   *cnd = eps;
    1125           0 :   if (umodiu(p, o) == 1)
    1126           0 :     return ellsympow_abelian_twist(E, p, m, o);
    1127             :   else
    1128           0 :     return ellsympow_nonabelian(p, m, bet);
    1129             : }
    1130             : 
    1131             : static long
    1132          70 : ellsympow_inertia3(GEN E, long vN)
    1133             : {
    1134          70 :   long vD = Z_lval(ell_get_disc(E), 3);
    1135          70 :   if (vN==2) return vD%2==0 ? 2: 4;
    1136          70 :   if (vN==4) return vD%4==0 ? 3: 6;
    1137          70 :   if (vN==3 || vN==5) return 12;
    1138           0 :   return 0;
    1139             : }
    1140             : 
    1141             : static long
    1142          70 : ellsympow_deltam3(long o, long m, long vN)
    1143             : {
    1144          70 :   if (o==3 || o==6) return ellsympow_epsm(3, m);
    1145          70 :   if (o==12 && vN ==3) return (ellsympow_epsm(3, m))/2;
    1146           0 :   if (o==12 && vN ==5) return (ellsympow_epsm(3, m))*3/2;
    1147           0 :   return 0;
    1148             : }
    1149             : 
    1150             : static long
    1151           0 : ellsympow_isabelian3(GEN E)
    1152             : {
    1153           0 :   ulong c4 = umodiu(ell_get_c4(E),81), c6 = umodiu(ell_get_c6(E), 243);
    1154           0 :   return (c4 == 27 || (c4%27==9 && (c6==108 || c6==135)));
    1155             : }
    1156             : 
    1157             : static long
    1158          35 : ellsympow_rootno3(GEN E, GEN p, long o, long m)
    1159             : {
    1160          35 :   const long  w6p[]={1,-1,-1,-1,1,1};
    1161          35 :   const long  w6n[]={-1,1,-1,1,-1,1};
    1162          35 :   const long w12p[]={1,1,-1,1,1,1};
    1163          35 :   const long w12n[]={-1,-1,-1,-1,-1,1};
    1164          35 :   long w = ellrootno(E, p), mm = (m%12)>>1;
    1165          35 :   switch(o)
    1166             :   {
    1167           0 :     case 2: return m%4== 1 ? -1: 1;
    1168           0 :     case 6:  return w == 1 ? w6p[mm]: w6n[mm];
    1169          35 :     case 12: return w == 1 ? w12p[mm]: w12n[mm];
    1170           0 :     default: return 1;
    1171             :   }
    1172             : }
    1173             : 
    1174             : static GEN
    1175          70 : ellsympow_goodred3(GEN E, GEN F, GEN p, long m, long vN, long *cnd, long *w)
    1176             : {
    1177          70 :   long o = ellsympow_inertia3(E, vN);
    1178          70 :   long bet = ellsympow_betam(o, m);
    1179          70 :   *cnd = m + 1 - bet + ellsympow_deltam3(o, m, vN);
    1180          70 :   *w = odd(m)? ellsympow_rootno3(E, p, o, m): 1;
    1181          70 :   if (o==1 || o==2)
    1182           0 :     return ellsympow_abelian(p, ellap(F, p), m, o);
    1183          70 :   if ((o==3 || o==6) && ellsympow_isabelian3(F))
    1184           0 :     return ellsympow_abelian(p, p, m, o);
    1185             :   else
    1186          70 :     return ellsympow_nonabelian(p, m, bet);
    1187             : }
    1188             : 
    1189             : static long
    1190          28 : ellsympow_inertia2(GEN F, long vN)
    1191             : {
    1192          28 :   long vM = itos(gel(elllocalred(F, gen_2),1));
    1193          28 :   GEN c6 = ell_get_c6(F);
    1194          28 :   long v6 = signe(c6) ? vali(c6): 24;
    1195          28 :   if (vM==0) return vN==0 ? 1: 2;
    1196          28 :   if (vM==2) return vN==2 ? 3: 6;
    1197          14 :   if (vM==5) return 8;
    1198           7 :   if (vM==8) return v6>=9? 8: 4;
    1199           7 :   if (vM==3 || vN==7) return 24;
    1200           0 :   return 0;
    1201             : }
    1202             : 
    1203             : static long
    1204          28 : ellsympow_deltam2(long o, long m, long vN)
    1205             : {
    1206          28 :   if ((o==2 || o==6) && vN==4) return ellsympow_epsm(2, m);
    1207          28 :   if ((o==2 || o==6) && vN==6) return 2*ellsympow_epsm(2, m);
    1208          28 :   if (o==4) return 2*ellsympow_epsm(4, m)+ellsympow_epsm(2, m);
    1209          28 :   if (o==8 && vN==5) return ellsympow_epsm(8, m)+ellsympow_epsm(2, m)/2;
    1210          21 :   if (o==8 && vN==6) return ellsympow_epsm(8, m)+ellsympow_epsm(2, m);
    1211          21 :   if (o==8 && vN==8) return ellsympow_epsm(8, m)+ellsympow_epsm(4, m)+ellsympow_epsm(2, m);
    1212          21 :   if (o==24 && vN==3) return (2*ellsympow_epsm(8, m)+ellsympow_epsm(2, m))/6;
    1213          14 :   if (o==24 && vN==4) return (ellsympow_epsm(8, m)+ellsympow_epsm(2, m)*2)/3;
    1214          14 :   if (o==24 && vN==6) return (ellsympow_epsm(8, m)+ellsympow_epsm(2, m)*5)/3;
    1215          14 :   if (o==24 && vN==7) return (ellsympow_epsm(8, m)*10+ellsympow_epsm(2, m)*5)/6;
    1216          14 :   return 0;
    1217             : }
    1218             : 
    1219             : static long
    1220           0 : ellsympow_isabelian2(GEN F)
    1221           0 : { return umodi2n(ell_get_c4(F),7) == 96; }
    1222             : 
    1223             : static long
    1224           0 : ellsympow_rootno2(GEN E, long vN, long m, long bet)
    1225             : {
    1226           0 :   long eps2 = (m + 1 - bet)>>1;
    1227           0 :   long eta = odd(vN) && m%8==3 ? -1 : 1;
    1228           0 :   long w2 = odd(eps2) ? ellrootno(E, gen_2): 1;
    1229           0 :   return eta == w2 ? 1 : -1;
    1230             : }
    1231             : 
    1232             : static GEN
    1233          28 : ellsympow_goodred2(GEN E, GEN F, GEN p, long m, long vN, long *cnd, long *w)
    1234             : {
    1235          28 :   long o = ellsympow_inertia2(F, vN);
    1236          28 :   long bet = ellsympow_betam(o, m);
    1237          28 :   *cnd = m + 1 - bet + ellsympow_deltam2(o, m, vN);
    1238          28 :   *w = odd(m) ? ellsympow_rootno2(E, vN, m, bet): 1;
    1239          28 :   if (o==1 || o==2)
    1240           0 :     return ellsympow_abelian(p, ellap(F, p), m, o);
    1241          28 :   if (o==4 && ellsympow_isabelian2(F))
    1242           0 :     return ellsympow_abelian(p, p, m, o);
    1243             :   else
    1244          28 :     return ellsympow_nonabelian(p, m, bet);
    1245             : }
    1246             : 
    1247             : static GEN
    1248         189 : ellminimaldotwist(GEN E, GEN *pD, long prec)
    1249             : {
    1250         189 :   GEN D = ellminimaltwistcond(E), Et = elltwist(E, D), Etmin;
    1251         189 :   if (pD) *pD = D;
    1252         189 :   Et = ellinit(Et, NULL, prec);
    1253         189 :   Etmin = ellminimalmodel(Et, NULL);
    1254         189 :   obj_free(Et); return Etmin;
    1255             : }
    1256             : 
    1257             : /* Based on
    1258             : Symmetric powers of elliptic curve L-functions,
    1259             : Phil Martin and Mark Watkins, ANTS VII
    1260             : <http://magma.maths.usyd.edu.au/users/watkins/papers/antsVII.pdf>
    1261             : with thanks to Mark Watkins. BA20180402
    1262             : */
    1263             : static GEN
    1264         140 : lfunellsympow(GEN e, ulong m)
    1265             : {
    1266         140 :   pari_sp av = avma;
    1267             :   GEN B, N, Nfa, pr, ex, ld, bad, ejd, et, pole;
    1268         140 :   long i, l, mero, w = (m&7)==1 || (m&7)==3 ? -1: 1;
    1269         140 :   checkell_Q(e);
    1270         140 :   e = ellminimalmodel(e, NULL);
    1271         140 :   ejd = Q_denom(ell_get_j(e));
    1272         140 :   mero = m==0 || (m%4==0 && ellQ_get_CM(e)<0);
    1273         140 :   ellQ_get_Nfa(e, &N, &Nfa);
    1274         140 :   pr = gel(Nfa,1);
    1275         140 :   ex = gel(Nfa,2); l = lg(pr);
    1276         140 :   if (ugcd(umodiu(N,6), 6) == 1)
    1277          21 :     et = NULL;
    1278             :   else
    1279         119 :     et = ellminimaldotwist(e, NULL, DEFAULTPREC);
    1280         140 :   B = gen_1;
    1281         140 :   bad = cgetg(l, t_VEC);
    1282         336 :   for (i=1; i<l; i++)
    1283             :   {
    1284         196 :     long vN = itos(gel(ex,i));
    1285         196 :     GEN p = gel(pr,i), eul;
    1286             :     long cnd, wp;
    1287         196 :     if (dvdii(ejd, p))
    1288          98 :       eul = ellsympow_multred(e, p, m, vN, &cnd, &wp);
    1289          98 :     else if (equaliu(p, 2))
    1290          28 :       eul = ellsympow_goodred2(e, et, p, m, vN, &cnd, &wp);
    1291          70 :     else if (equaliu(p, 3))
    1292          70 :       eul = ellsympow_goodred3(e, et, p, m, vN, &cnd, &wp);
    1293             :     else
    1294           0 :       eul = ellsympow_goodred(e, p, m, &cnd, &wp);
    1295         196 :     gel(bad, i) = mkvec2(p, ginv(eul));
    1296         196 :     B = mulii(B, powiu(p,cnd));
    1297         196 :     w *= wp;
    1298             :   }
    1299         140 :   pole = mero ? mkvec(mkvec2(stoi(1+(m>>1)),gen_0)): NULL;
    1300         280 :   ld = mkvecn(mero? 7: 6, tag(mkvec2(mkvec2(e,utoi(m)),bad), t_LFUN_SYMPOW_ELL),
    1301         140 :         gen_0, ellsympow_gamma(m), stoi(m+1), B, stoi(w), pole);
    1302         140 :   if (et) obj_free(et);
    1303         140 :   return gerepilecopy(av, ld);
    1304             : }
    1305             : 
    1306             : GEN
    1307          70 : lfunsympow(GEN ldata, ulong m)
    1308             : {
    1309          70 :   ldata = lfunmisc_to_ldata_shallow(ldata);
    1310          70 :   if (ldata_get_type(ldata) != t_LFUN_ELL)
    1311           0 :     pari_err_IMPL("lfunsympow");
    1312          70 :   return lfunellsympow(gel(ldata_get_an(ldata), 2), m);
    1313             : }
    1314             : 
    1315             : static GEN
    1316          28 : lfunmfspec_i(GEN lmisc, long bit)
    1317             : {
    1318             :   GEN linit, ldataf, v, ve, vo, om, op, B, dom;
    1319             :   long k, k2, j;
    1320             : 
    1321          28 :   ldataf = lfunmisc_to_ldata_shallow(lmisc);
    1322          28 :   if (!gequal(ldata_get_gammavec(ldataf), mkvec2(gen_0,gen_1)))
    1323           0 :     pari_err_TYPE("lfunmfspec", lmisc);
    1324          28 :   k = gtos(ldata_get_k(ldataf));
    1325          28 :   if (k == 1) return mkvec2(cgetg(1, t_VEC), gen_1);
    1326          21 :   dom = mkvec3(dbltor(k/2.), dbltor((k-2)/2.), gen_0);
    1327          21 :   if (is_linit(lmisc) && linit_get_type(lmisc) == t_LDESC_INIT
    1328           0 :       && sdomain_isincl((double)k, dom, lfun_get_dom(linit_get_tech(lmisc))))
    1329           0 :     linit = lmisc;
    1330             :   else
    1331          21 :     linit = lfuninit(ldataf, dom, 0, bit);
    1332          21 :   B = int2n(bit/4);
    1333          21 :   v = cgetg(k, t_VEC);
    1334          21 :   for (j = 1; j < k; j++) gel(v,j) = lfunlambda(linit, utoi(j), bit);
    1335          21 :   om = gel(v,1);
    1336          21 :   if (odd(k)) return mkvec2(bestappr(gdiv(v, om), B), om);
    1337             : 
    1338           7 :   k2 = k/2;
    1339           7 :   ve = cgetg(k2, t_VEC);
    1340           7 :   vo = cgetg(k2+1, t_VEC);
    1341           7 :   gel(vo,1) = om;
    1342          42 :   for (j = 1; j < k2; j++)
    1343             :   {
    1344          35 :     gel(ve,j) = gel(v,2*j);
    1345          35 :     gel(vo,j+1) = gel(v,2*j+1);
    1346             :   }
    1347           7 :   if (k2 == 1) { om = gen_1;    op = gel(v,1); }
    1348           7 :   else         { om = gel(v,2); op = gel(v,3); }
    1349           7 :   if (maxss(gexpo(imag_i(om)), gexpo(imag_i(op))) > -bit/2)
    1350           0 :     pari_err_TYPE("lfunmfspec", lmisc);
    1351           7 :   ve = gdiv(ve, om);
    1352           7 :   vo = gdiv(vo, op);
    1353           7 :   return mkvec4(bestappr(ve,B), bestappr(vo,B), om, op);
    1354             : }
    1355             : GEN
    1356          28 : lfunmfspec(GEN lmisc, long bit)
    1357             : {
    1358          28 :   pari_sp av = avma;
    1359          28 :   return gerepilecopy(av, lfunmfspec_i(lmisc, bit));
    1360             : }
    1361             : 
    1362             : static long
    1363          28 : ellsymsq_bad2(GEN c4, GEN c6, long e)
    1364             : {
    1365          28 :   switch (e)
    1366             :   {
    1367          14 :     case 2: return 1;
    1368           7 :     case 3: return 0;
    1369           7 :     case 5: return 0;
    1370           0 :     case 7: return 0;
    1371             :     case 8:
    1372           0 :       if (dvdiu(c6,512)) return 0;
    1373           0 :       return umodiu(c4,128)==32 ? 1 : -1;
    1374           0 :     default: return 0;
    1375             :   }
    1376             : }
    1377             : static long
    1378          14 : ellsymsq_bad3(GEN c4, GEN c6, long e)
    1379             : {
    1380             :   long c6_243, c4_81;
    1381          14 :   switch (e)
    1382             :   {
    1383           0 :     case 2: return 1;
    1384          14 :     case 3: return 0;
    1385           0 :     case 5: return 0;
    1386             :     case 4:
    1387           0 :       c4_81 = umodiu(c4,81);
    1388           0 :       if (c4_81 == 27) return -1;
    1389           0 :       if (c4_81%27 != 9) return 1;
    1390           0 :       c6_243 = umodiu(c6,243);
    1391           0 :       return (c6_243==108 || c6_243==135)? -1: 1;
    1392           0 :     default: return 0;
    1393             :   }
    1394             : }
    1395             : static int
    1396           0 : c4c6_testp(GEN c4, GEN c6, GEN p)
    1397           0 : { GEN p2 = sqri(p); return (dvdii(c6,p2) && !dvdii(c4,p2)); }
    1398             : /* assume e = v_p(N) >= 2 */
    1399             : static long
    1400          42 : ellsymsq_badp(GEN c4, GEN c6, GEN p, long e)
    1401             : {
    1402          42 :   if (absequaliu(p, 2)) return ellsymsq_bad2(c4, c6, e);
    1403          14 :   if (absequaliu(p, 3)) return ellsymsq_bad3(c4, c6, e);
    1404           0 :   switch(umodiu(p, 12UL))
    1405             :   {
    1406           0 :     case 1: return -1;
    1407           0 :     case 5: return c4c6_testp(c4,c6,p)? -1: 1;
    1408           0 :     case 7: return c4c6_testp(c4,c6,p)?  1:-1;
    1409           0 :     default:return 1; /* p%12 = 11 */
    1410             :   }
    1411             : }
    1412             : static GEN
    1413          70 : lfunellsymsqmintwist(GEN e)
    1414             : {
    1415          70 :   pari_sp av = avma;
    1416             :   GEN N, Nfa, P, E, V, c4, c6, ld;
    1417             :   long i, l, k;
    1418          70 :   checkell_Q(e);
    1419          70 :   e = ellminimalmodel(e, NULL);
    1420          70 :   ellQ_get_Nfa(e, &N, &Nfa);
    1421          70 :   c4 = ell_get_c4(e);
    1422          70 :   c6 = ell_get_c6(e);
    1423          70 :   P = gel(Nfa,1); l = lg(P);
    1424          70 :   E = gel(Nfa,2);
    1425          70 :   V = cgetg(l, t_VEC);
    1426         196 :   for (i=k=1; i<l; i++)
    1427             :   {
    1428         126 :     GEN p = gel(P,i);
    1429         126 :     long a, e = itos(gel(E,i));
    1430         126 :     if (e == 1) continue;
    1431          42 :     a = ellsymsq_badp(c4, c6, p, e);
    1432          42 :     gel(V,k++) = mkvec2(p, stoi(a));
    1433             :   }
    1434          70 :   setlg(V, k);
    1435          70 :   ld = lfunellsympow(e, 2);
    1436          70 :   return gerepilecopy(av, mkvec2(ld, V));
    1437             : }
    1438             : 
    1439             : static GEN
    1440          70 : mfpeters(GEN ldata2, GEN fudge, GEN N, long k, long bitprec)
    1441             : {
    1442          70 :   GEN t, L = real_i(lfun(ldata2, stoi(k), bitprec));
    1443          70 :   long prec = nbits2prec(bitprec);
    1444          70 :   t = powrs(mppi(prec), k+1); shiftr_inplace(t, 2*k-1); /* Pi/2 * (4Pi)^k */
    1445          70 :   return gmul(gdiv(gmul(mulii(N,mpfact(k-1)), fudge), t), L);
    1446             : }
    1447             : 
    1448             : /* Assume E to be twist-minimal */
    1449             : static GEN
    1450          70 : lfunellmfpetersmintwist(GEN E, long bitprec)
    1451             : {
    1452          70 :   pari_sp av = avma;
    1453          70 :   GEN symsq, veceuler, N = ellQ_get_N(E), fudge = gen_1;
    1454          70 :   long j, k = 2;
    1455          70 :   symsq = lfunellsymsqmintwist(E);
    1456          70 :   veceuler = gel(symsq,2);
    1457         112 :   for (j = 1; j < lg(veceuler); j++)
    1458             :   {
    1459          42 :     GEN v = gel(veceuler,j), p = gel(v,1), q = powis(p,1-k);
    1460          42 :     long s = signe(gel(v,2));
    1461          42 :     if (s) fudge = gmul(fudge, s==1 ? gaddsg(1, q): gsubsg(1, q));
    1462             :   }
    1463          70 :   return gerepileupto(av, mfpeters(gel(symsq,1),fudge,N,k,bitprec));
    1464             : }
    1465             : 
    1466             : /* From Christophe Delaunay, http://delaunay.perso.math.cnrs.fr/these.pdf */
    1467             : static GEN
    1468          70 : elldiscfix(GEN E, GEN Et, GEN D)
    1469             : {
    1470          70 :   GEN N = ellQ_get_N(E), Nt = ellQ_get_N(Et);
    1471          70 :   GEN P = gel(absZ_factor(D), 1);
    1472          70 :   GEN f = gen_1;
    1473          70 :   long i, l = lg(P);
    1474         133 :   for (i=1; i < l; i++)
    1475             :   {
    1476          63 :     GEN r, p = gel(P,i);
    1477          63 :     long v = Z_pval(N, p), vt = Z_pval(Nt, p);
    1478          63 :     if (v <= vt) continue;
    1479             :     /* v > vt */
    1480          49 :     if (absequaliu(p, 2))
    1481             :     {
    1482          28 :       if (vt == 0 && v >= 4)
    1483           0 :         r = shifti(subsi(9, sqri(ellap(Et, p))), v-3);  /* 9=(2+1)^2 */
    1484          28 :       else if (vt == 1)
    1485           7 :         r = gmul2n(utoipos(3), v-3);  /* not in Z if v=2 */
    1486          21 :       else if (vt >= 2)
    1487          21 :         r = int2n(v-vt);
    1488             :       else
    1489           0 :         r = gen_1; /* vt = 0, 1 <= v <= 3 */
    1490             :     }
    1491          21 :     else if (vt >= 1)
    1492          14 :       r = gdiv(subiu(sqri(p), 1), p);
    1493             :     else
    1494           7 :       r = gdiv(mulii(subiu(p, 1), subii(sqri(addiu(p, 1)), sqri(ellap(Et, p)))), p);
    1495          49 :     f = gmul(f, r);
    1496             :   }
    1497          70 :   return f;
    1498             : }
    1499             : 
    1500             : GEN
    1501          70 : lfunellmfpeters(GEN E, long bitprec)
    1502             : {
    1503          70 :   pari_sp ltop = avma;
    1504          70 :   GEN D, Et = ellminimaldotwist(E, &D, nbits2prec(bitprec));
    1505          70 :   GEN nor = lfunellmfpetersmintwist(Et, bitprec);
    1506          70 :   GEN nor2 = gmul(nor, elldiscfix(E, Et, D));
    1507          70 :   obj_free(Et); return gerepileupto(ltop, nor2);
    1508             : }
    1509             : 
    1510             : /*************************************************************/
    1511             : /*               Genus 2 curves                              */
    1512             : /*************************************************************/
    1513             : 
    1514             : static void
    1515      181453 : Flv_diffnext(GEN d, ulong p)
    1516             : {
    1517      181453 :   long j, n = lg(d)-1;
    1518     1267110 :   for(j = n; j>=2; j--)
    1519     1085575 :     uel(d,j) = Fl_add(uel(d,j), uel(d,j-1), p);
    1520      181535 : }
    1521             : 
    1522             : static GEN
    1523        1764 : Flx_difftable(GEN P, ulong p)
    1524             : {
    1525        1764 :   long i, n = degpol(P);
    1526        1764 :   GEN V = cgetg(n+2, t_VECSMALL);
    1527        1764 :   uel(V, n+1) = Flx_constant(P);
    1528       12331 :   for(i = n; i >= 1; i--)
    1529             :   {
    1530       10568 :     P = Flx_diff1(P, p);
    1531       10565 :     uel(V, i) = Flx_constant(P);
    1532             :   }
    1533        1763 :   return V;
    1534             : }
    1535             : 
    1536             : static long
    1537        1762 : Flx_genus2trace_naive(GEN H, ulong p)
    1538             : {
    1539        1762 :   pari_sp av = avma;
    1540             :   ulong i, j;
    1541        1762 :   long a, n = degpol(H);
    1542        1762 :   GEN k = const_vecsmall(p, -1), d;
    1543        1764 :   k[1] = 0;
    1544       96190 :   for (i=1, j=1; i < p; i += 2, j = Fl_add(j, i, p))
    1545       94426 :     k[j+1] = 1;
    1546        1764 :   a = n == 5 ? 0: k[1+Flx_lead(H)];
    1547        1764 :   d = Flx_difftable(H, p);
    1548      183360 :   for (i=0; i < p; i++)
    1549             :   {
    1550      181598 :     a += k[1+uel(d,n+1)];
    1551      181598 :     Flv_diffnext(d, p);
    1552             :   }
    1553        1762 :   set_avma(av);
    1554        1763 :   return a;
    1555             : }
    1556             : 
    1557             : static GEN
    1558        1889 : dirgenus2(GEN Q, GEN p, long n)
    1559             : {
    1560        1889 :   pari_sp av = avma;
    1561             :   GEN f;
    1562        1889 :   if (n > 2)
    1563         126 :     f = RgX_recip(hyperellcharpoly(gmul(Q,gmodulo(gen_1, p))));
    1564             :   else
    1565             :   {
    1566        1763 :     ulong pp = itou(p);
    1567        1763 :     GEN Qp = ZX_to_Flx(Q, pp);
    1568        1762 :     long t = Flx_genus2trace_naive(Qp, pp);
    1569        1763 :     f = deg1pol_shallow(stoi(t), gen_1, 0);
    1570             :   }
    1571        1888 :   return gerepileupto(av, RgXn_inv_i(f, n));
    1572             : }
    1573             : 
    1574             : GEN
    1575         756 : dirgenus2_worker(GEN P, ulong X, GEN Q)
    1576             : {
    1577         756 :   pari_sp av = avma;
    1578         756 :   long i, l = lg(P);
    1579         756 :   GEN V = cgetg(l, t_VEC);
    1580        2644 :   for(i = 1; i < l; i++)
    1581             :   {
    1582        1888 :     ulong p = uel(P,i);
    1583        1888 :     long d = ulogint(X, p) + 1; /* minimal d such that p^d > X */
    1584        1889 :     gel(V,i) = dirgenus2(Q, utoi(uel(P,i)), d);
    1585             :   }
    1586         756 :   return gerepilecopy(av, mkvec2(P,V));
    1587             : }
    1588             : 
    1589             : static GEN
    1590         168 : vecan_genus2(GEN an, long L)
    1591             : {
    1592         168 :   GEN Q = gel(an,1), bad = gel(an, 2);
    1593         168 :   GEN worker = strtoclosure("_dirgenus2_worker", 1, Q);
    1594         168 :   return pardireuler(worker, gen_2, stoi(L), NULL, bad);
    1595             : }
    1596             : 
    1597             : static GEN
    1598          42 : genus2_redmodel(GEN P, GEN p)
    1599             : {
    1600          42 :   GEN M = FpX_factor(P, p);
    1601          42 :   GEN F = gel(M,1), E = gel(M,2);
    1602          42 :   long i, k, r = lg(F);
    1603          42 :   GEN U = scalarpol(leading_coeff(P), varn(P));
    1604          42 :   GEN G = cgetg(r, t_COL);
    1605         140 :   for (i=1, k=0; i<r; i++)
    1606             :   {
    1607          98 :     if (E[i]>1)
    1608          49 :       gel(G,++k) = gel(F,i);
    1609          98 :     if (odd(E[i]))
    1610          63 :       U = FpX_mul(U, gel(F,i), p);
    1611             :   }
    1612          42 :   setlg(G,++k);
    1613          42 :   return mkvec2(G,U);
    1614             : }
    1615             : 
    1616             : static GEN
    1617         280 : oneminusxd(long d)
    1618             : {
    1619         280 :   return gsub(gen_1, pol_xn(d, 0));
    1620             : }
    1621             : 
    1622             : static GEN
    1623          21 : ellfromeqncharpoly(GEN P, GEN Q, GEN p)
    1624             : {
    1625             :   long v;
    1626             :   GEN E, F, t, y;
    1627          21 :   v = fetch_var();
    1628          21 :   y = pol_x(v);
    1629          21 :   F = gsub(gadd(ZX_sqr(y), gmul(y, Q)), P);
    1630          21 :   E = ellinit(ellfromeqn(F), p, DEFAULTPREC);
    1631          21 :   delete_var();
    1632          21 :   t = ellap(E, p);
    1633          21 :   obj_free(E);
    1634          21 :   return mkpoln(3, gen_1, negi(t), p);
    1635             : }
    1636             : 
    1637             : static GEN
    1638          42 : genus2_eulerfact(GEN P, GEN p)
    1639             : {
    1640          42 :   GEN Pp = FpX_red(P, p);
    1641          42 :   GEN GU = genus2_redmodel(Pp, p);
    1642          42 :   long d = 6-degpol(Pp), v = d/2, w = odd(d);
    1643             :   GEN abe, tor;
    1644          42 :   GEN ki, kp = pol_1(0), kq = pol_1(0);
    1645          42 :   GEN F = gel(GU,1), Q = gel(GU,2);
    1646          42 :   long dQ = degpol(Q), lF = lg(F)-1;
    1647             : 
    1648          42 :   abe = dQ >= 5 ? RgX_recip(hyperellcharpoly(gmul(Q,gmodulo(gen_1,p))))
    1649         105 :       : dQ >= 3 ? RgX_recip(ellfromeqncharpoly(Q,gen_0,p))
    1650          63 :                 : pol_1(0);
    1651          42 :   ki = dQ != 0 ? oneminusxd(1)
    1652          56 :               : Fp_issquare(gel(Q,2),p) ? ZX_sqr(oneminusxd(1))
    1653          14 :                                         : oneminusxd(2);
    1654          42 :   if (lF)
    1655             :   {
    1656             :     long i;
    1657          84 :     for(i=1; i <= lF; i++)
    1658             :     {
    1659          49 :       GEN Fi = gel(F, i);
    1660          49 :       long d = degpol(Fi);
    1661          49 :       GEN e = FpX_rem(Q, Fi, p);
    1662          84 :       GEN kqf = lgpol(e)==0 ? oneminusxd(d):
    1663          70 :                 FpXQ_issquare(e, Fi, p) ? ZX_sqr(oneminusxd(d))
    1664          70 :                                         : oneminusxd(2*d);
    1665          49 :       kp = gmul(kp, oneminusxd(d));
    1666          49 :       kq = gmul(kq, kqf);
    1667             :     }
    1668             :   }
    1669          42 :   if (v)
    1670             :   {
    1671           7 :     GEN kqoo = w==1 ? oneminusxd(1):
    1672           0 :                Fp_issquare(leading_coeff(Q), p)? ZX_sqr(oneminusxd(1))
    1673           0 :                                               : oneminusxd(2);
    1674           7 :     kp = gmul(kp, oneminusxd(1));
    1675           7 :     kq = gmul(kq, kqoo);
    1676             :   }
    1677          42 :   tor = RgX_div(ZX_mul(oneminusxd(1), kq), ZX_mul(ki, kp));
    1678          42 :   return ginv( ZX_mul(abe, tor) );
    1679             : }
    1680             : 
    1681             : static GEN
    1682          28 : F2x_genus2_find_trans(GEN P, GEN Q, GEN F)
    1683             : {
    1684          28 :   pari_sp av = avma;
    1685          28 :   long i, d = F2x_degree(F), v = P[1];
    1686             :   GEN M, C, V;
    1687          28 :   M = cgetg(d+1, t_MAT);
    1688          84 :   for (i=1; i<=d; i++)
    1689             :   {
    1690          56 :     GEN Mi = F2x_rem(F2x_add(F2x_shift(Q,i-1), monomial_F2x(2*i-2,v)), F);
    1691          56 :     gel(M,i) = F2x_to_F2v(Mi, d);
    1692             :   }
    1693          28 :   C = F2x_to_F2v(F2x_rem(P, F), d);
    1694          28 :   V = F2m_F2c_invimage(M, C);
    1695          28 :   return gerepileuptoleaf(av, F2v_to_F2x(V, v));
    1696             : }
    1697             : 
    1698             : static GEN
    1699          35 : F2x_genus2_trans(GEN P, GEN Q, GEN H)
    1700             : {
    1701          35 :   return F2x_add(P,F2x_add(F2x_mul(H,Q), F2x_sqr(H)));
    1702             : }
    1703             : 
    1704             : static GEN
    1705          42 : F2x_genus_redoo(GEN P, GEN Q, long k)
    1706             : {
    1707          42 :   if (F2x_degree(P)==2*k)
    1708             :   {
    1709          14 :     long c = F2x_coeff(P,2*k-1), dQ = F2x_degree(Q);
    1710          14 :     if ((dQ==k-1 && c==1) || (dQ<k-1 && c==0))
    1711           7 :      return F2x_genus2_trans(P, Q, monomial_F2x(k, P[1]));
    1712             :   }
    1713          35 :   return P;
    1714             : }
    1715             : 
    1716             : static GEN
    1717          35 : F2x_pseudodisc(GEN P, GEN Q)
    1718             : {
    1719          35 :   GEN dP = F2x_deriv(P), dQ = F2x_deriv(Q);
    1720          35 :   return F2x_gcd(Q, F2x_add(F2x_mul(P, F2x_sqr(dQ)), F2x_sqr(dP)));
    1721             : }
    1722             : 
    1723             : static GEN
    1724          14 : F2x_genus_red(GEN P, GEN Q)
    1725             : {
    1726             :   long dP, dQ;
    1727             :   GEN F, FF;
    1728          14 :   P = F2x_genus_redoo(P, Q, 3);
    1729          14 :   P = F2x_genus_redoo(P, Q, 2);
    1730          14 :   P = F2x_genus_redoo(P, Q, 1);
    1731          14 :   dP = F2x_degree(P);
    1732          14 :   dQ = F2x_degree(Q);
    1733          14 :   FF = F = F2x_pseudodisc(P,Q);
    1734          49 :   while(F2x_degree(F)>0)
    1735             :   {
    1736          21 :     GEN M = gel(F2x_factor(F),1);
    1737          21 :     long i, l = lg(M);
    1738          49 :     for(i=1; i<l; i++)
    1739             :     {
    1740          28 :       GEN R = F2x_sqr(gel(M,i));
    1741          28 :       GEN H = F2x_genus2_find_trans(P, Q, R);
    1742          28 :       P = F2x_div(F2x_genus2_trans(P, Q, H), R);
    1743          28 :       Q = F2x_div(Q, gel(M,i));
    1744             :     }
    1745          21 :     F = F2x_pseudodisc(P, Q);
    1746             :   }
    1747          14 :   return mkvec4(P,Q,FF,mkvecsmall2(dP,dQ));
    1748             : }
    1749             : 
    1750             : /* Number of solutions of x^2+b*x+c */
    1751             : static long
    1752          21 : F2xqX_quad_nbroots(GEN b, GEN c, GEN T)
    1753             : {
    1754          21 :   if (lgpol(b) > 0)
    1755             :   {
    1756          14 :     GEN d = F2xq_div(c, F2xq_sqr(b, T), T);
    1757          14 :     return F2xq_trace(d, T)? 0: 2;
    1758             :   }
    1759             :   else
    1760           7 :     return 1;
    1761             : }
    1762             : 
    1763             : static GEN
    1764          28 : genus2_redmodel2(GEN P)
    1765             : {
    1766          28 :   GEN Q = pol_0(varn(P));
    1767          28 :   GEN P2 = ZX_to_F2x(P);
    1768          77 :   while (F2x_issquare(P2))
    1769             :   {
    1770          21 :     GEN H = F2x_to_ZX(F2x_sqrt(P2));
    1771          21 :     GEN P1 = ZX_sub(P, ZX_sqr(H));
    1772          21 :     GEN Q1 = ZX_add(Q, ZX_mulu(H, 2));
    1773          21 :     if ((signe(P1)==0 ||  ZX_lval(P1,2)>=2)
    1774          21 :      && (signe(Q1)==0 ||  ZX_lval(Q1,2)>=1))
    1775             :     {
    1776          21 :       P = ZX_shifti(P1, -2);
    1777          21 :       Q = ZX_shifti(Q1, -1);
    1778          21 :       P2= ZX_to_F2x(P);
    1779             :     } else break;
    1780             :   }
    1781          28 :   return mkvec2(P,Q);
    1782             : }
    1783             : 
    1784             : static GEN
    1785          14 : genus2_eulerfact2(GEN PQ)
    1786             : {
    1787          14 :   GEN V = F2x_genus_red(ZX_to_F2x(gel(PQ, 1)), ZX_to_F2x(gel(PQ, 2)));
    1788          14 :   GEN P = gel(V, 1), Q = gel(V, 2);
    1789          14 :   GEN F = gel(V, 3), v = gel(V, 4);
    1790             :   GEN abe, tor;
    1791          14 :   GEN ki, kp = pol_1(0), kq = pol_1(0);
    1792          14 :   long dP = F2x_degree(P), dQ = F2x_degree(Q), d = maxss(dP, 2*dQ);
    1793          14 :   if (!lgpol(F)) return pol_1(0);
    1794          21 :   ki = dQ!=0 || dP>0 ? oneminusxd(1):
    1795           7 :       dP==-1 ? ZX_sqr(oneminusxd(1)): oneminusxd(2);
    1796          28 :   abe = d>=5? RgX_recip(hyperellcharpoly(gmul(PQ,gmodulss(1,2)))):
    1797          14 :         d>=3? RgX_recip(ellfromeqncharpoly(F2x_to_ZX(P), F2x_to_ZX(Q), gen_2)):
    1798             :         pol_1(0);
    1799          14 :   if (lgpol(F))
    1800             :   {
    1801          14 :     GEN M = gel(F2x_factor(F), 1);
    1802          14 :     long i, lF = lg(M)-1;
    1803          35 :     for(i=1; i <= lF; i++)
    1804             :     {
    1805          21 :       GEN Fi = gel(M, i);
    1806          21 :       long d = F2x_degree(Fi);
    1807          21 :       long nb  = F2xqX_quad_nbroots(F2x_rem(Q, Fi), F2x_rem(P, Fi), Fi);
    1808          35 :       GEN kqf = nb==1 ? oneminusxd(d):
    1809           0 :                 nb==2 ? ZX_sqr(oneminusxd(d))
    1810          14 :                       : oneminusxd(2*d);
    1811          21 :       kp = gmul(kp, oneminusxd(d));
    1812          21 :       kq = gmul(kq, kqf);
    1813             :     }
    1814             :   }
    1815          14 :   if (maxss(v[1],2*v[2])<5)
    1816             :   {
    1817          14 :     GEN kqoo = v[1]>2*v[2] ? oneminusxd(1):
    1818           7 :                v[1]<2*v[2] ? ZX_sqr(oneminusxd(1))
    1819           7 :                            : oneminusxd(2);
    1820           7 :     kp = gmul(kp, oneminusxd(1));
    1821           7 :     kq = gmul(kq, kqoo);
    1822             :   }
    1823          14 :   tor = RgX_div(ZX_mul(oneminusxd(1),kq), ZX_mul(ki, kp));
    1824          14 :   return ginv( ZX_mul(abe, tor) );
    1825             : }
    1826             : 
    1827             : GEN
    1828          28 : lfungenus2(GEN G)
    1829             : {
    1830          28 :   pari_sp ltop = avma;
    1831             :   GEN Ldata;
    1832          28 :   GEN gr = genus2red(G, NULL);
    1833          28 :   GEN N  = gel(gr, 1), M = gel(gr, 2), Q = gel(gr, 3), L = gel(gr, 4);
    1834          28 :   GEN PQ = genus2_redmodel2(Q);
    1835             :   GEN e;
    1836          28 :   long i, lL = lg(L), ram2;
    1837          28 :   ram2 = absequaliu(gmael(M,1,1),2);
    1838          28 :   if (ram2 && equalis(gmael(M,2,1),-1))
    1839          14 :     pari_warn(warner,"unknown valuation of conductor at 2");
    1840          28 :   e = cgetg(lL+(ram2?0:1), t_VEC);
    1841          42 :   gel(e,1) = mkvec2(gen_2, ram2 ? genus2_eulerfact2(PQ)
    1842          14 :            : ginv( RgX_recip(hyperellcharpoly(gmul(PQ,gmodulss(1,2))))) );
    1843          70 :   for(i = ram2? 2: 1; i < lL; i++)
    1844             :   {
    1845          42 :     GEN Li = gel(L, i);
    1846          42 :     GEN p = gel(Li, 1);
    1847          42 :     gel(e, ram2 ? i: i+1) = mkvec2(p, genus2_eulerfact(Q,p));
    1848             :   }
    1849          28 :   Ldata = mkvecn(6, tag(mkvec2(Q,e), t_LFUN_GENUS2),
    1850             :       gen_0, mkvec4(gen_0, gen_0, gen_1, gen_1), gen_2, N, gen_0);
    1851          28 :   return gerepilecopy(ltop, Ldata);
    1852             : }
    1853             : 
    1854             : /*************************************************************/
    1855             : /*                        ETA QUOTIENTS                      */
    1856             : /* An eta quotient is a matrix with 2 columns [m, r_m] with  */
    1857             : /* m >= 1 representing f(\tau)=\prod_m\eta(m\tau)^{r_m}.     */
    1858             : /*************************************************************/
    1859             : 
    1860             : /* eta(x^v) + O(x^L) */
    1861             : GEN
    1862         714 : eta_ZXn(long v, long L)
    1863             : {
    1864         714 :   long n, k = 0, v2 = 2*v, bn = v, cn = 0;
    1865             :   GEN P;
    1866         714 :   if (!L) return zeropol(0);
    1867         714 :   P = cgetg(L+2,t_POL); P[1] = evalsigne(1);
    1868         714 :   for(n = 0; n < L; n++) gel(P,n+2) = gen_0;
    1869        2653 :   for(n = 0;; n++, bn += v2, cn += v)
    1870        1939 :   { /* k = v * (3*n-1) * n / 2; bn = v * (2*n+1); cn = v * n */
    1871             :     long k2;
    1872        2653 :     gel(P, k+2) = odd(n)? gen_m1: gen_1;
    1873        2653 :     k2 = k+cn; if (k2 >= L) break;
    1874        2429 :     k = k2;
    1875             :     /* k = v * (3*n+1) * n / 2 */;
    1876        2429 :     gel(P, k+2) = odd(n)? gen_m1: gen_1;
    1877        2429 :     k2 = k+bn; if (k2 >= L) break;
    1878        1939 :     k = k2;
    1879             :   }
    1880         714 :   setlg(P, k+3); return P;
    1881             : }
    1882             : GEN
    1883         259 : eta_product_ZXn(GEN eta, long L)
    1884             : {
    1885         259 :   pari_sp av = avma;
    1886         259 :   GEN P = NULL, D = gel(eta,1), R = gel(eta,2);
    1887         259 :   long i, l = lg(D);
    1888         952 :   for (i = 1; i < l; ++i)
    1889             :   {
    1890         693 :     GEN Q = eta_ZXn(D[i], L);
    1891         693 :     long r = R[i];
    1892         693 :     if (r < 0) { Q = RgXn_inv_i(Q, L); r = -r; }
    1893         693 :     if (r != 1) Q = RgXn_powu_i(Q, r, L);
    1894         693 :     P = P? ZXn_mul(P, Q, L): Q;
    1895         693 :     if (gc_needed(av,1) && i > 1)
    1896             :     {
    1897           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"eta_product_ZXn");
    1898           0 :       P = gerepilecopy(av, P);
    1899             :     }
    1900             :   }
    1901         259 :   return P;
    1902             : }
    1903             : static GEN
    1904         126 : vecan_eta(GEN an, long L)
    1905             : {
    1906         126 :   GEN v = RgX_to_RgC(eta_product_ZXn(an, L), L);
    1907         126 :   settyp(v, t_VEC); return v;
    1908             : }
    1909             : 
    1910             : /* return 1 if cuspidal, 0 if holomorphic, -1 otherwise */
    1911             : static int
    1912         175 : etacuspidal(GEN N, GEN k, GEN B, GEN R, GEN NB)
    1913             : {
    1914         175 :   long i, j, lD, l, cusp = 1;
    1915             :   GEN D;
    1916         175 :   if (gsigne(k) < 0) return -1;
    1917         168 :   D = divisors(N); lD = lg(D); l = lg(B);
    1918        1162 :   for (i = 1; i < lD; i++)
    1919             :   {
    1920         994 :     GEN t = gen_0, d = gel(D,i);
    1921             :     long s;
    1922        3157 :     for (j = 1; j < l; j++)
    1923        2163 :       t = addii(t, mulii(gel(NB,j), mulii(gel(R,j), sqri(gcdii(d, gel(B,j))))));
    1924         994 :     s = signe(t);
    1925         994 :     if (s < 0) return -1;
    1926         994 :     if (!s) cusp = 0;
    1927             :   }
    1928         168 :   return cusp;
    1929             : }
    1930             : /* u | 24, level = u*N0, N0 = lcm(B), NB[i] = N0/B[i] */
    1931             : static int
    1932          42 : etaselfdual(GEN B, GEN R, GEN NB, ulong u)
    1933             : {
    1934          42 :   long i, l = lg(B);
    1935         140 :   for (i = 1; i < l; ++i)
    1936             :   {
    1937          98 :     GEN t = muliu(gel(NB,i), u); /* *pN / B[i] */
    1938          98 :     long j = ZV_search(B, t);
    1939          98 :     if (!j || !equalii(gel(R,i),gel(R,j))) return 0;
    1940             :   }
    1941          42 :   return 1;
    1942             : }
    1943             : /* return Nebentypus of eta quotient, k2 = 2*k integral */
    1944             : static GEN
    1945         140 : etachar(GEN B, GEN R, GEN k2)
    1946             : {
    1947         140 :   long i, l = lg(B);
    1948         140 :   GEN P = gen_1;
    1949         140 :   for (i = 1; i < l; ++i) if (mpodd(gel(R,i))) P = mulii(P, gel(B,i));
    1950         140 :   switch(Mod4(k2))
    1951             :   {
    1952          84 :     case 0: break;
    1953          28 :     case 2:  P = negi(P); break;
    1954          28 :     default: P = shifti(P, 1); break;
    1955             :   }
    1956         140 :   return coredisc(P);
    1957             : }
    1958             : /* Return 0 if not on gamma_0(N). Sets conductor, modular weight, character,
    1959             :  * canonical matrix, v_q(eta), sd = 1 iff self-dual, cusp = 1 iff cuspidal
    1960             :  * [0 if holomorphic at all cusps, else -1] */
    1961             : long
    1962         189 : etaquotype(GEN *peta, GEN *pN, GEN *pk, GEN *CHI, long *pv, long *sd,
    1963             :            long *cusp)
    1964             : {
    1965         189 :   GEN B, R, S, T, N, NB, eta = *peta;
    1966             :   long l, i, u, S24;
    1967             : 
    1968         189 :   if (lg(eta) != 3) pari_err_TYPE("lfunetaquo", eta);
    1969         189 :   switch(typ(eta))
    1970             :   {
    1971          77 :     case t_VEC: eta = mkmat2(mkcol(gel(eta,1)), mkcol(gel(eta,2))); break;
    1972         112 :     case t_MAT: break;
    1973           0 :     default: pari_err_TYPE("lfunetaquo", eta);
    1974             :   }
    1975         189 :   if (!RgV_is_ZVpos(gel(eta,1)) || !RgV_is_ZV(gel(eta,2)))
    1976           0 :     pari_err_TYPE("lfunetaquo", eta);
    1977         189 :   *peta = eta = famat_reduce(eta);
    1978         189 :   B = gel(eta,1); l = lg(B); /* sorted in increasing order */
    1979         189 :   R = gel(eta,2);
    1980         189 :   N = glcm0(B, NULL);
    1981         189 :   NB = cgetg(l, t_VEC);
    1982         189 :   for (i = 1; i < l; i++) gel(NB,i) = diviiexact(N, gel(B,i));
    1983         189 :   S = gen_0; T = gen_0; u = 0;
    1984         546 :   for (i = 1; i < l; ++i)
    1985             :   {
    1986         357 :     GEN b = gel(B,i), r = gel(R,i);
    1987         357 :     S = addii(S, mulii(b, r));
    1988         357 :     T = addii(T, r);
    1989         357 :     u += umodiu(r,24) * umodiu(gel(NB,i), 24);
    1990             :   }
    1991         189 :   S = divis_rem(S, 24, &S24);
    1992         189 :   if (S24) return 0; /* non-integral valuation at oo */
    1993         182 :   u %= 24; if (u < 0) u += 24;
    1994         182 :   u = 24/ugcd(24, u);
    1995         182 :   *pN = muliu(N, u); /* level */
    1996         182 :   *pk = gmul2n(T,-1); /* weight */
    1997         182 :   *pv = itos(S); /* valuation */
    1998         182 :   if (cusp) *cusp = etacuspidal(*pN, *pk, B, R, NB);
    1999         182 :   if (sd) *sd = etaselfdual(B, R, NB, u);
    2000         182 :   if (CHI) *CHI = etachar(B, R, T);
    2001         182 :   return 1;
    2002             : }
    2003             : 
    2004             : GEN
    2005          42 : lfunetaquo(GEN eta0)
    2006             : {
    2007          42 :   pari_sp ltop = avma;
    2008          42 :   GEN Ldata, N, BR, k, eta = eta0;
    2009             :   long v, sd, cusp;
    2010          42 :   if (!etaquotype(&eta, &N, &k, NULL, &v, &sd, &cusp))
    2011           0 :     pari_err_TYPE("lfunetaquo", eta0);
    2012          42 :   if (!cusp) pari_err_IMPL("noncuspidal eta quotient");
    2013          42 :   if (v != 1) pari_err_IMPL("valuation != 1");
    2014          42 :   if (!sd) pari_err_IMPL("non self-dual eta quotient");
    2015          42 :   if (typ(k) != t_INT) pari_err_TYPE("lfunetaquo [non-integral weight]", eta0);
    2016          42 :   BR = mkvec2(ZV_to_zv(gel(eta,1)), ZV_to_zv(gel(eta,2)));
    2017          42 :   Ldata = mkvecn(6, tag(BR,t_LFUN_ETA), gen_0, mkvec2(gen_0,gen_1), k,N, gen_1);
    2018          42 :   return gerepilecopy(ltop, Ldata);
    2019             : }
    2020             : 
    2021             : static GEN
    2022         455 : vecan_qf(GEN Q, long L)
    2023             : {
    2024         455 :   GEN v, w = qfrep0(Q, utoi(L), 1);
    2025             :   long i;
    2026         455 :   v = cgetg(L+1, t_VEC);
    2027         455 :   for (i = 1; i <= L; i++) gel(v,i) = utoi(2 * w[i]);
    2028         455 :   return v;
    2029             : }
    2030             : 
    2031             : long
    2032         154 : qfiseven(GEN M)
    2033             : {
    2034         154 :   long i, l = lg(M);
    2035         357 :   for (i=1; i<l; i++)
    2036         280 :     if (mpodd(gcoeff(M,i,i))) return 0;
    2037          77 :   return 1;
    2038             : }
    2039             : 
    2040             : GEN
    2041          28 : lfunqf(GEN M, long prec)
    2042             : {
    2043          28 :   pari_sp ltop = avma;
    2044             :   long n;
    2045             :   GEN k, D, d, Mi, Ldata, poles, res0, res2, eno, dual;
    2046             : 
    2047          28 :   if (typ(M) != t_MAT) pari_err_TYPE("lfunqf", M);
    2048          28 :   if (!RgM_is_ZM(M))   pari_err_TYPE("lfunqf [not integral]", M);
    2049          28 :   n = lg(M)-1;
    2050          28 :   k = sstoQ(n,2);
    2051          28 :   M = Q_primpart(M);
    2052          28 :   Mi = ZM_inv(M, &d); /* d M^(-1) */
    2053          28 :   if (!qfiseven(M)) { M = gmul2n(M, 1); d = shifti(d,1); }
    2054          28 :   if (!qfiseven(Mi)){ Mi= gmul2n(Mi,1); d = shifti(d,1); }
    2055             :   /* det(Mi) = d^n/det(M), D^2 = det(Mi)/det(M) */
    2056          28 :   D = gdiv(gpow(d,k,prec), ZM_det(M));
    2057          28 :   if (!issquareall(D, &eno)) eno = gsqrt(D, prec);
    2058          28 :   dual = gequal1(D) ? gen_0: tag(Mi, t_LFUN_QF);
    2059          28 :   res0 = RgX_to_ser(deg1pol_shallow(gen_m2, gen_0, 0), 3);
    2060          28 :   setvalp(res0, -1);
    2061          28 :   res2 = RgX_to_ser(deg1pol_shallow(gmulgs(eno,2), gen_0, 0), 3);
    2062          28 :   setvalp(res2, -1);
    2063          28 :   poles = mkcol2(mkvec2(k,res2), mkvec2(gen_0,res0));
    2064          28 :   Ldata = mkvecn(7, tag(M, t_LFUN_QF), dual,
    2065             :        mkvec2(gen_0, gen_1), k, d, eno, poles);
    2066          28 :   return gerepilecopy(ltop, Ldata);
    2067             : }
    2068             : 
    2069             : /********************************************************************/
    2070             : /**  Artin L function, based on a GP script by Charlotte Euvrard   **/
    2071             : /********************************************************************/
    2072             : 
    2073             : static GEN
    2074         119 : artin_charfromgens(GEN G, GEN M)
    2075             : {
    2076         119 :   GEN R, V, ord = gal_get_orders(G), grp = gal_get_group(G);
    2077         119 :   long i, j, k, n = lg(ord)-1, m = lg(grp)-1;
    2078             : 
    2079         119 :   if (lg(M)-1 != n) pari_err_DIM("lfunartin");
    2080         119 :   R = cgetg(m+1, t_VEC);
    2081         119 :   gel(R, 1) = matid(lg(gel(M, 1))-1);
    2082         357 :   for (i = 1, k = 1; i <= n; ++i)
    2083             :   {
    2084         238 :     long c = k*(ord[i] - 1);
    2085         238 :     gel(R, ++k) = gel(M, i);
    2086         238 :     for (j = 2; j <= c; ++j) gel(R, ++k) = gmul(gel(R,j), gel(M,i));
    2087             :   }
    2088         119 :   V = cgetg(m+1, t_VEC);
    2089         119 :   for (i = 1; i <= m; i++) gel(V, gel(grp,i)[1]) = gtrace(gel(R,i));
    2090         119 :   return V;
    2091             : }
    2092             : 
    2093             : static GEN
    2094          70 : galois_get_conj(GEN G)
    2095             : {
    2096          70 :   GEN grp = gal_get_group(G);
    2097          70 :   long i, k, r = lg(grp)-1;
    2098          70 :   GEN b = zero_F2v(r);
    2099         224 :   for (k = 2; k <= r; ++k)
    2100             :   {
    2101         224 :     GEN g = gel(grp,k);
    2102         224 :     if (!F2v_coeff(b,g[1]) && g[g[1]]==1)
    2103             :     {
    2104          70 :       pari_sp av = avma;
    2105          70 :       GEN F = galoisfixedfield(G, g, 1, -1);
    2106          70 :       if (ZX_sturmpart(F, NULL) > 0) { set_avma(av); return g; }
    2107           0 :       for (i = 1; i<=r; i++)
    2108             :       {
    2109           0 :         GEN h = gel(grp, i);
    2110           0 :         long t = h[1];
    2111           0 :         while (h[t]!=1) t = h[t];
    2112           0 :         F2v_set(b, h[g[t]]);
    2113             :       }
    2114           0 :       set_avma(av);
    2115             :     }
    2116             :   }
    2117           0 :   pari_err_BUG("galois_get_conj");
    2118           0 :   return NULL;
    2119             : }
    2120             : 
    2121        7700 : static GEN  cyclotoi(GEN v) { return simplify_shallow(lift_shallow(v)); }
    2122        2891 : static long cyclotos(GEN v) { return gtos(cyclotoi(v)); }
    2123        2821 : static long char_dim(GEN ch) { return cyclotos(gel(ch,1)); }
    2124             : 
    2125             : static GEN
    2126        1379 : artin_gamma(GEN N, GEN G, GEN ch)
    2127             : {
    2128        1379 :   long a, t, d = char_dim(ch);
    2129        1379 :   if (nf_get_r2(N) == 0) return vec01(d, 0);
    2130          70 :   a = galois_get_conj(G)[1];
    2131          70 :   t = cyclotos(gel(ch,a));
    2132          70 :   return vec01((d+t) / 2, (d-t) / 2);
    2133             : }
    2134             : 
    2135             : static long
    2136        3213 : artin_dim(GEN ind, GEN ch)
    2137             : {
    2138        3213 :   long n = lg(ch)-1;
    2139        3213 :   GEN elts = group_elts(ind, n);
    2140        3213 :   long i, d = lg(elts)-1;
    2141        3213 :   GEN s = gen_0;
    2142       18123 :   for(i=1; i<=d; i++)
    2143       14910 :     s = gadd(s, gel(ch, gel(elts,i)[1]));
    2144        3213 :   return gtos(gdivgs(cyclotoi(s), d));
    2145             : }
    2146             : 
    2147             : static GEN
    2148         623 : artin_ind(GEN elts, GEN ch, GEN p)
    2149             : {
    2150         623 :   long i, d = lg(elts)-1;
    2151         623 :   GEN s = gen_0;
    2152        2149 :   for(i=1; i<=d; i++)
    2153        1526 :     s = gadd(s, gel(ch, gmul(gel(elts,i),p)[1]));
    2154         623 :   return gdivgs(s, d);
    2155             : }
    2156             : 
    2157             : static GEN
    2158        2282 : artin_ram(GEN nf, GEN gal, GEN aut, GEN pr, GEN ramg, GEN ch, long d)
    2159             : {
    2160        2282 :   pari_sp av = avma;
    2161             :   long i, v, n;
    2162             :   GEN p, q, V, elts;
    2163        2282 :   if (d==0) return pol_1(0);
    2164         616 :   n = degpol(gal_get_pol(gal));
    2165         616 :   q = p = idealramfrobenius_aut(nf, gal, pr, ramg, aut);
    2166         616 :   elts = group_elts(gel(ramg,2), n);
    2167         616 :   v = fetch_var_higher();
    2168         616 :   V = cgetg(d+3, t_POL);
    2169         616 :   V[1] = evalsigne(1)|evalvarn(v);
    2170         616 :   gel(V,2) = gen_0;
    2171        1239 :   for(i=1; i<=d; i++)
    2172             :   {
    2173         623 :     gel(V,i+2) = gdivgs(artin_ind(elts, ch, q), -i);
    2174         623 :     q = gmul(q, p);
    2175             :   }
    2176         616 :   delete_var();
    2177         616 :   V = RgXn_exp(V,d+1);
    2178         616 :   setvarn(V,0); return gerepileupto(av, ginv(V));
    2179             : }
    2180             : 
    2181             : /* [Artin conductor, vec of [p, Lp]] */
    2182             : static GEN
    2183        1379 : artin_badprimes(GEN N, GEN G, GEN aut, GEN ch)
    2184             : {
    2185        1379 :   pari_sp av = avma;
    2186        1379 :   long i, d = char_dim(ch);
    2187        1379 :   GEN P = gel(absZ_factor(nf_get_disc(N)), 1);
    2188        1379 :   long lP = lg(P);
    2189        1379 :   GEN B = cgetg(lP, t_VEC), C = cgetg(lP, t_VEC);
    2190             : 
    2191        3661 :   for (i = 1; i < lP; ++i)
    2192             :   {
    2193        2282 :     GEN p = gel(P, i), pr = idealprimedec_galois(N, p);
    2194        2282 :     GEN J = idealramgroups_aut(N, G, pr, aut);
    2195        2282 :     GEN G0 = gel(J,2); /* inertia group */
    2196        2282 :     long lJ = lg(J);
    2197        2282 :     long sdec = artin_dim(G0, ch);
    2198        2282 :     long ndec = group_order(G0);
    2199        2282 :     long j, v = ndec * (d - sdec);
    2200        3213 :     for (j = 3; j < lJ; ++j)
    2201             :     {
    2202         931 :       GEN Jj = gel(J, j);
    2203         931 :       long s = artin_dim(Jj, ch);
    2204         931 :       v += group_order(Jj) * (d - s);
    2205             :     }
    2206        2282 :     gel(C, i) = powiu(p, v/ndec);
    2207        2282 :     gel(B, i) = mkvec2(p, artin_ram(N, G, aut, pr, J, ch, sdec));
    2208             :   }
    2209        1379 :   return gerepilecopy(av, mkvec2(ZV_prod(C), B));
    2210             : }
    2211             : 
    2212             : /* p does not divide nf.index */
    2213             : static GEN
    2214       51278 : idealfrobenius_easy(GEN nf, GEN gal, GEN aut, GEN T, GEN p)
    2215             : {
    2216       51278 :   long i, l = lg(aut), f = degpol(T);
    2217       51277 :   GEN D, Dzk, DzkT, DXp, grp = gal_get_group(gal);
    2218       51277 :   pari_sp av = avma;
    2219       51277 :   if (f==1) return gel(grp,1);
    2220       49513 :   Dzk = nf_get_zkprimpart(nf);
    2221       49512 :   D = modii(nf_get_zkden(nf), p);
    2222       49505 :   DzkT = RgV_to_RgM(FqV_red(Dzk, T, p), f);
    2223       49514 :   DXp = RgX_to_RgC(FpX_Frobenius(T, p), f);
    2224       49510 :   if (!equali1(D)) DXp = FpC_Fp_mul(DXp, D, p);
    2225      324667 :   for(i=1; i < l; i++)
    2226             :   {
    2227      324667 :     GEN g = gel(grp,i);
    2228      324667 :     if (perm_order(g)==f)
    2229             :     {
    2230      166222 :       GEN A = FpM_FpC_mul(DzkT, gel(aut,g[1]), p);
    2231      166185 :       if (ZV_equal(A, DXp)) {set_avma(av); return g; }
    2232             :     }
    2233             :   }
    2234             :   return NULL; /* LCOV_EXCL_LINE */
    2235             : }
    2236             : /* true nf; p divides nf.index, pr/p unramified */
    2237             : static GEN
    2238        1470 : idealfrobenius_hard(GEN nf, GEN gal, GEN aut, GEN pr)
    2239             : {
    2240        1470 :   long i, l = lg(aut), f = pr_get_f(pr);
    2241        1470 :   GEN modpr, p, T, X, Xp, pi, grp = gal_get_group(gal);
    2242        1470 :   pari_sp av = avma;
    2243        1470 :   if (f==1) return gel(grp,1);
    2244        1218 :   pi = pr_get_gen(pr);
    2245        1218 :   modpr = zkmodprinit(nf, pr);
    2246        1218 :   p = modpr_get_p(modpr);
    2247        1218 :   T = modpr_get_T(modpr);
    2248        1218 :   X = modpr_genFq(modpr);
    2249        1218 :   Xp = FpX_Frobenius(T, p);
    2250        8162 :   for (i = 1; i < l; i++)
    2251             :   {
    2252        8162 :     GEN g = gel(grp,i);
    2253        8162 :     if (perm_order(g)==f)
    2254             :     {
    2255        3808 :       GEN S = gel(aut,g[1]);
    2256        3808 :       GEN A = nf_to_Fq(nf, zk_galoisapplymod(nf,X,S,p), modpr);
    2257             :       /* sigma(X) = X^p (mod pr) and sigma(pi) in pr */
    2258        5138 :       if (ZX_equal(A, Xp) && (f == nf_get_degree(nf) ||
    2259        2548 :           ZC_prdvd(zk_galoisapplymod(nf,pi,S,p),pr))) { set_avma(av); return g; }
    2260             :     }
    2261             :   }
    2262             :   return NULL; /* LCOV_EXCL_LINE */
    2263             : }
    2264             : 
    2265             : static GEN
    2266       52745 : dirartin(GEN nf, GEN G, GEN V, GEN aut, GEN p, long n)
    2267             : {
    2268       52745 :   pari_sp av = avma;
    2269             :   GEN pr, frob;
    2270             :   /* pick one maximal ideal in the conjugacy class above p */
    2271       52745 :   GEN T = nf_get_pol(nf);
    2272       52743 :   if (!dvdii(nf_get_index(nf), p))
    2273             :   { /* simple case */
    2274       51265 :     GEN F = FpX_factor(T, p), P = gmael(F,1,1);
    2275       51278 :     frob = idealfrobenius_easy(nf, G, aut, P, p);
    2276             :   }
    2277             :   else
    2278             :   {
    2279        1470 :     pr = idealprimedec_galois(nf,p);
    2280        1470 :     frob = idealfrobenius_hard(nf, G, aut, pr);
    2281             :   }
    2282       52750 :   set_avma(av); return RgXn_inv(gel(V, frob[1]), n);
    2283             : }
    2284             : 
    2285             : GEN
    2286       15340 : dirartin_worker(GEN P, ulong X, GEN nf, GEN G, GEN V, GEN aut)
    2287             : {
    2288       15340 :   pari_sp av = avma;
    2289       15340 :   long i, l = lg(P);
    2290       15340 :   GEN W = cgetg(l, t_VEC);
    2291       68088 :   for(i = 1; i < l; i++)
    2292             :   {
    2293       52745 :     ulong p = uel(P,i);
    2294       52745 :     long d = ulogint(X, p) + 1; /* minimal d such that p^d > X */
    2295       52745 :     gel(W,i) = dirartin(nf, G, V, aut, utoi(uel(P,i)), d);
    2296             :   }
    2297       15343 :   return gerepilecopy(av, mkvec2(P,W));
    2298             : }
    2299             : 
    2300             : static GEN
    2301        2905 : vecan_artin(GEN an, long L, long prec)
    2302             : {
    2303        2905 :   GEN A, Sbad = gel(an,5);
    2304        2905 :   long n = itos(gel(an,6)), isreal = lg(an)<8 ? 0: !itos(gel(an,7));
    2305        2905 :   GEN worker = strtoclosure("_dirartin_worker", 4, gel(an,1), gel(an,2), gel(an,3), gel(an,4));
    2306        2905 :   A = lift_shallow(pardireuler(worker, gen_2, stoi(L), NULL, Sbad));
    2307        2905 :   A = RgXV_RgV_eval(A, grootsof1(n, prec));
    2308        2905 :   if (isreal) A = real_i(A);
    2309        2905 :   return A;
    2310             : }
    2311             : 
    2312             : static GEN
    2313        2856 : char_expand(GEN conj, GEN ch)
    2314             : {
    2315        2856 :   long i, l = lg(conj);
    2316        2856 :   GEN V = cgetg(l, t_COL);
    2317        2856 :   for (i=1; i<l; i++) gel(V,i) = gel(ch, conj[i]);
    2318        2856 :   return V;
    2319             : }
    2320             : 
    2321             : static GEN
    2322        1596 : handle_zeta(long n, GEN ch, long *m)
    2323             : {
    2324             :   GEN c;
    2325        1596 :   long t, i, l = lg(ch);
    2326        1596 :   GEN dim = cyclotoi(vecsum(ch));
    2327        1596 :   if (typ(dim) != t_INT)
    2328           0 :     pari_err_DOMAIN("lfunartin","chi","is not a", strtoGENstr("character"), ch);
    2329        1596 :   t = itos(dim);
    2330        1596 :   if (t < 0 || t % n)
    2331           0 :     pari_err_DOMAIN("lfunartin","chi","is not a", strtoGENstr("character"), ch);
    2332        1596 :   if (t == 0) { *m = 0; return ch; }
    2333         224 :   *m = t / n;
    2334         224 :   c = cgetg(l, t_COL);
    2335        2065 :   for (i=1; i<l; i++)
    2336        1841 :     gel(c,i) = gsubgs(gel(ch,i), *m);
    2337         224 :   return c;
    2338             : }
    2339             : 
    2340             : static int
    2341        6496 : cyclo_is_real(GEN v, GEN ix)
    2342             : {
    2343        6496 :   pari_sp av = avma;
    2344        6496 :   GEN w = poleval(lift_shallow(v), ix);
    2345        6496 :   return gc_bool(av, gequal(w, v));
    2346             : }
    2347             : 
    2348             : static int
    2349        1379 : char_is_real(GEN ch, GEN mod)
    2350             : {
    2351        1379 :   long i, l = lg(ch);
    2352        1379 :   GEN ix = QXQ_inv(pol_x(varn(mod)), mod);
    2353        7014 :   for (i=1; i<l; i++)
    2354        6496 :     if (!cyclo_is_real(gel(ch,i), ix)) return 0;
    2355         518 :   return 1;
    2356             : }
    2357             : 
    2358             : GEN
    2359        1610 : lfunartin(GEN nf, GEN gal, GEN ch, long o, long bitprec)
    2360             : {
    2361        1610 :   pari_sp av = avma;
    2362        1610 :   GEN bc, V, aut, mod, Ldata = NULL, chx, cc, conj, repr;
    2363             :   long tmult, var;
    2364        1610 :   nf = checknf(nf);
    2365        1610 :   checkgal(gal);
    2366        1610 :   var = gvar(ch);
    2367        1610 :   if (var == 0) pari_err_PRIORITY("lfunartin",ch,"=",0);
    2368        1610 :   if (var < 0) var = 1;
    2369        1610 :   if (!is_vec_t(typ(ch))) pari_err_TYPE("lfunartin", ch);
    2370        1610 :   cc = group_to_cc(gal);
    2371        1610 :   conj = gel(cc,2);
    2372        1610 :   repr = gel(cc,3);
    2373        1610 :   mod = mkpolmod(gen_1, polcyclo(o, var));
    2374        1610 :   if (lg(ch)>1 && typ(gel(ch,1))==t_MAT)
    2375         119 :     chx = artin_charfromgens(gal, gmul(ch,mod));
    2376             :   else
    2377             :   {
    2378        1491 :     if (lg(repr) != lg(ch)) pari_err_DIM("lfunartin");
    2379        1477 :     chx = char_expand(conj, gmul(ch,mod));
    2380             :   }
    2381        1596 :   chx = handle_zeta(nf_get_degree(nf), chx, &tmult);
    2382        1596 :   ch = shallowextract(chx, repr);
    2383        1596 :   if (!gequal0(chx))
    2384             :   {
    2385        1379 :     GEN real = char_is_real(chx, gel(mod,1))? gen_0: gen_1;
    2386        1379 :     aut = nfgaloispermtobasis(nf, gal);
    2387        1379 :     V = gmul(char_expand(conj, galoischarpoly(gal, ch, o)), mod);
    2388        1379 :     bc = artin_badprimes(nf, gal, aut, chx);
    2389        4137 :     Ldata = mkvecn(6,
    2390        1379 :       tag(mkcoln(7, nf, gal, V, aut, gel(bc, 2), stoi(o), real), t_LFUN_ARTIN),
    2391        1379 :       real, artin_gamma(nf, gal, chx), gen_1, gel(bc,1), gen_0);
    2392             :   }
    2393        1596 :   if (tmult==0 && Ldata==NULL) pari_err_TYPE("lfunartin",ch);
    2394        1596 :   if (tmult)
    2395             :   {
    2396             :     long i;
    2397         224 :     if (Ldata==NULL) { Ldata = lfunzeta(); tmult--; }
    2398         231 :     for(i=1; i<=tmult; i++)
    2399           7 :       Ldata = lfunmul(Ldata, gen_1, bitprec);
    2400             :   }
    2401        1596 :   return gerepilecopy(av, Ldata);
    2402             : }
    2403             : 
    2404             : static GEN
    2405          21 : lfunzetakinit_artin(GEN nf, GEN gal, GEN dom, long der, long bitprec)
    2406             : {
    2407          21 :   pari_sp ltop = avma;
    2408          21 :   GEN To = galoischartable(gal), T = gel(To, 1);
    2409          21 :   long o = itos(gel(To, 2));
    2410             :   GEN F, E, M, domain;
    2411          21 :   long i, l = lg(T);
    2412          21 :   F = cgetg(l, t_VEC);
    2413          21 :   E = cgetg(l, t_VECSMALL);
    2414          84 :   for (i = 1; i < l; ++i)
    2415             :   {
    2416          63 :     GEN L = lfunartin(nf, gal, gel(T,i), o, bitprec);
    2417          63 :     gel(F, i) = lfuninit(L, dom, der, bitprec);
    2418          63 :     E[i] = char_dim(gel(T,i));
    2419             :   }
    2420          21 :   domain = mkvec2(dom, mkvecsmall2(der, bitprec));
    2421          21 :   M = mkvec3(F, E, zero_zv(l-1));
    2422          21 :   return gerepilecopy(ltop, lfuninit_make(t_LDESC_PRODUCT, lfunzetak_i(nf),
    2423             :                                           M, domain));
    2424             : }
    2425             : 
    2426             : /********************************************************************/
    2427             : /**                    High-level Constructors                     **/
    2428             : /********************************************************************/
    2429             : enum { t_LFUNMISC_POL, t_LFUNMISC_CHIQUAD, t_LFUNMISC_CHICONREY,
    2430             :        t_LFUNMISC_CHIGEN, t_LFUNMISC_ELLINIT, t_LFUNMISC_ETAQUO };
    2431             : static long
    2432        3556 : lfundatatype(GEN data)
    2433             : {
    2434             :   long l;
    2435        3556 :   switch(typ(data))
    2436             :   {
    2437         686 :     case t_INT: return t_LFUNMISC_CHIQUAD;
    2438          21 :     case t_INTMOD: return t_LFUNMISC_CHICONREY;
    2439         399 :     case t_POL: return t_LFUNMISC_POL;
    2440             :     case t_VEC:
    2441        2450 :       if (checknf_i(data)) return t_LFUNMISC_POL;
    2442        2338 :       l = lg(data);
    2443        2338 :       if (l == 17) return t_LFUNMISC_ELLINIT;
    2444        1078 :       if (l == 3 && typ(gel(data,1)) == t_VEC) return t_LFUNMISC_CHIGEN;
    2445           0 :       break;
    2446             :   }
    2447           0 :   return -1;
    2448             : }
    2449             : static GEN
    2450       40036 : lfunmisc_to_ldata_i(GEN ldata, long shallow)
    2451             : {
    2452             :   long lx;
    2453       40036 :   if (is_linit(ldata)) ldata = linit_get_ldata(ldata);
    2454       40036 :   lx = lg(ldata);
    2455       40036 :   if (typ(ldata)==t_VEC && (lx == 7 || lx == 8) && is_tagged(ldata))
    2456             :   {
    2457       36480 :     if (!shallow) ldata = gcopy(ldata);
    2458       36480 :     checkldata(ldata); return ldata;
    2459             :   }
    2460        3556 :   switch (lfundatatype(ldata))
    2461             :   {
    2462         511 :     case t_LFUNMISC_POL: return lfunzetak(ldata);
    2463         686 :     case t_LFUNMISC_CHIQUAD: return lfunchiquad(ldata);
    2464             :     case t_LFUNMISC_CHICONREY:
    2465             :     {
    2466          21 :       GEN G = znstar0(gel(ldata,1), 1);
    2467          21 :       return lfunchiZ(G, gel(ldata,2));
    2468             :     }
    2469        1078 :     case t_LFUNMISC_CHIGEN: return lfunchigen(gel(ldata,1), gel(ldata,2));
    2470        1260 :     case t_LFUNMISC_ELLINIT: return lfunell(ldata);
    2471             :   }
    2472           0 :   pari_err_TYPE("lfunmisc_to_ldata",ldata);
    2473             :   return NULL; /* LCOV_EXCL_LINE */
    2474             : }
    2475             : 
    2476             : GEN
    2477         595 : lfunmisc_to_ldata(GEN ldata)
    2478         595 : { return lfunmisc_to_ldata_i(ldata, 0); }
    2479             : 
    2480             : GEN
    2481       39441 : lfunmisc_to_ldata_shallow(GEN ldata)
    2482       39441 : { return lfunmisc_to_ldata_i(ldata, 1); }
    2483             : 
    2484             : /********************************************************************/
    2485             : /**                    High-level an expansion                     **/
    2486             : /********************************************************************/
    2487             : /* van is the output of ldata_get_an: return a_1,...a_L at precision prec */
    2488             : GEN
    2489       15141 : ldata_vecan(GEN van, long L, long prec)
    2490             : {
    2491       15141 :   GEN an = gel(van, 2);
    2492       15141 :   long t = mael(van,1,1);
    2493             :   pari_timer ti;
    2494       15141 :   if (DEBUGLEVEL >= 1)
    2495           0 :     err_printf("Lfun: computing %ld coeffs, prec %ld, type %ld\n", L, prec, t);
    2496       15141 :   if (DEBUGLEVEL >= 2) timer_start(&ti);
    2497       15141 :   switch (t)
    2498             :   {
    2499             :     long n;
    2500             :     case t_LFUN_GENERIC:
    2501        1855 :       an = vecan_closure(an, L, prec);
    2502        1841 :       n = lg(an)-1;
    2503        1841 :       if (n < L)
    2504          21 :         pari_warn(warner, "#an = %ld < %ld, results may be imprecise", n, L);
    2505        1841 :       break;
    2506        2765 :     case t_LFUN_ZETA: an = const_vecsmall(L, 1); break;
    2507        1260 :     case t_LFUN_NF:  an = dirzetak(an, stoi(L)); break;
    2508             :     case t_LFUN_ELL:
    2509        1435 :       an = (ell_get_type(an) == t_ELL_Q) ? ellanQ_zv(an, L): ellan(an, L);
    2510        1435 :       break;
    2511         980 :     case t_LFUN_KRONECKER: an = vecan_Kronecker(an, L); break;
    2512         539 :     case t_LFUN_CHIZ: an = vecan_chiZ(an, L, prec); break;
    2513         672 :     case t_LFUN_CHIGEN: an = vecan_chigen(an, L, prec); break;
    2514        2905 :     case t_LFUN_ARTIN: an = vecan_artin(an, L, prec); break;
    2515         126 :     case t_LFUN_ETA: an = vecan_eta(an, L); break;
    2516         455 :     case t_LFUN_QF: an = vecan_qf(an, L); break;
    2517         623 :     case t_LFUN_DIV: an = vecan_div(an, L, prec); break;
    2518         224 :     case t_LFUN_MUL: an = vecan_mul(an, L, prec); break;
    2519           0 :     case t_LFUN_CONJ: an = vecan_conj(an, L, prec); break;
    2520         315 :     case t_LFUN_SYMPOW_ELL: an = vecan_ellsympow(an, L); break;
    2521         168 :     case t_LFUN_GENUS2: an = vecan_genus2(an, L); break;
    2522             :     case t_LFUN_MFCLOS:
    2523             :     {
    2524         462 :       GEN F = gel(an,1), E = gel(an,2), c = gel(an,3);
    2525         462 :       an = mfcoefs(F,L,1) + 1; /* skip a_0 */
    2526         462 :       an[0] = evaltyp(t_VEC)|evallg(L+1);
    2527         462 :       an = mfvecembed(E, an);
    2528         462 :       if (!isint1(c)) an = RgV_Rg_mul(an,c);
    2529         462 :       break;
    2530             :     }
    2531         357 :     case t_LFUN_TWIST: an = vecan_twist(an, L, prec); break;
    2532           0 :     default: pari_err_TYPE("ldata_vecan", van);
    2533             :   }
    2534       15127 :   if (DEBUGLEVEL >= 2) timer_printf(&ti, "ldata_vecan");
    2535       15127 :   return an;
    2536             : }

Generated by: LCOV version 1.13