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

Generated by: LCOV version 1.13