Code coverage tests

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

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

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

LCOV - code coverage report
Current view: top level - basemath - modsym.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.10.0 lcov report (development 22303-eb3e11d) Lines: 2549 2615 97.5 %
Date: 2018-04-21 06:16:28 Functions: 274 275 99.6 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2011  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., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. */
      13             : 
      14             : #include "pari.h"
      15             : #include "paripriv.h"
      16             : 
      17             : /* Adapted from shp_package/moments by Robert Pollack
      18             :  * http://www.math.mcgill.ca/darmon/programs/shp/shp.html */
      19             : static GEN mskinit(ulong N, long k, long sign);
      20             : static GEN mshecke_i(GEN W, ulong p);
      21             : static GEN ZSl2_star(GEN v);
      22             : static GEN getMorphism(GEN W1, GEN W2, GEN v);
      23             : static GEN voo_act_Gl2Q(GEN g, long k);
      24             : 
      25             : /* Input: P^1(Z/NZ) (formed by create_p1mod)
      26             :    Output: # P^1(Z/NZ) */
      27             : static long
      28        6664 : p1_size(GEN p1N) { return lg(gel(p1N,1)) - 1; }
      29             : static ulong
      30     6762980 : p1N_get_N(GEN p1N) { return gel(p1N,3)[2]; }
      31             : static GEN
      32     2797228 : p1N_get_hash(GEN p1N) { return gel(p1N,2); }
      33             : static GEN
      34        1792 : p1N_get_fa(GEN p1N) { return gel(p1N,4); }
      35             : static GEN
      36        1687 : p1N_get_div(GEN p1N) { return gel(p1N,5); }
      37             : static GEN
      38     2677087 : p1N_get_invsafe(GEN p1N) { return gel(p1N,6); }
      39             : static GEN
      40      904162 : p1N_get_inverse(GEN p1N) { return gel(p1N,7); }
      41             : 
      42             : /* ms-specific accessors */
      43             : /* W = msinit, return the output of msinit_N */
      44             : static GEN
      45      656439 : get_msN(GEN W) { return lg(W) == 4? gel(W,1): W; }
      46             : static GEN
      47      750113 : msN_get_p1N(GEN W) { return gel(W,1); }
      48             : static GEN
      49      208446 : msN_get_genindex(GEN W) { return gel(W,5); }
      50             : static GEN
      51     5381873 : msN_get_E2fromE1(GEN W) { return gel(W,7); }
      52             : static GEN
      53        1393 : msN_get_annT2(GEN W) { return gel(W,8); }
      54             : static GEN
      55        1393 : msN_get_annT31(GEN W) { return gel(W,9); }
      56             : static GEN
      57        1358 : msN_get_singlerel(GEN W) { return gel(W,10); }
      58             : static GEN
      59      676956 : msN_get_section(GEN W) { return gel(W,12); }
      60             : 
      61             : static GEN
      62       83062 : ms_get_p1N(GEN W) { return msN_get_p1N(get_msN(W)); }
      63             : static long
      64       72961 : ms_get_N(GEN W) { return p1N_get_N(ms_get_p1N(W)); }
      65             : static GEN
      66        1869 : ms_get_hashcusps(GEN W) { W = get_msN(W); return gel(W,16); }
      67             : static GEN
      68       14770 : ms_get_section(GEN W) { return msN_get_section(get_msN(W)); }
      69             : static GEN
      70      204400 : ms_get_genindex(GEN W) { return msN_get_genindex(get_msN(W)); }
      71             : static long
      72      199605 : ms_get_nbgen(GEN W) { return lg(ms_get_genindex(W))-1; }
      73             : static long
      74      161602 : ms_get_nbE1(GEN W)
      75             : {
      76             :   GEN W11;
      77      161602 :   W = get_msN(W); W11 = gel(W,11);
      78      161602 :   return W11[4] - W11[3];
      79             : }
      80             : 
      81             : /* msk-specific accessors */
      82             : static long
      83          77 : msk_get_dim(GEN W) { return gmael(W,3,2)[2]; }
      84             : static GEN
      85       82838 : msk_get_basis(GEN W) { return gmael(W,3,1); }
      86             : static long
      87       45808 : msk_get_weight(GEN W) { return gmael(W,3,2)[1]; }
      88             : static long
      89       20762 : msk_get_sign(GEN W)
      90             : {
      91       20762 :   GEN t = gel(W,2);
      92       20762 :   return typ(t)==t_INT? 0: itos(gel(t,1));
      93             : }
      94             : static GEN
      95         623 : msk_get_star(GEN W) { return gmael(W,2,2); }
      96             : static GEN
      97        3626 : msk_get_starproj(GEN W) { return gmael(W,2,3); }
      98             : 
      99             : static int
     100         294 : is_Qevproj(GEN x)
     101         294 : { return typ(x) == t_VEC && lg(x) == 5 && typ(gel(x,1)) == t_MAT; }
     102             : long
     103         105 : msdim(GEN W)
     104             : {
     105         105 :   if (is_Qevproj(W)) return lg(gel(W,1)) - 1;
     106          91 :   checkms(W);
     107          91 :   if (!msk_get_sign(W)) return msk_get_dim(W);
     108          28 :   return lg(gel(msk_get_starproj(W), 1)) - 1;
     109             : }
     110             : long
     111          14 : msgetlevel(GEN W) { checkms(W); return ms_get_N(W); }
     112             : long
     113          14 : msgetweight(GEN W) { checkms(W); return msk_get_weight(W); }
     114             : long
     115          28 : msgetsign(GEN W) { checkms(W); return msk_get_sign(W); }
     116             : 
     117             : void
     118       35959 : checkms(GEN W)
     119             : {
     120       35959 :   if (typ(W) != t_VEC || lg(W) != 4)
     121           0 :     pari_err_TYPE("checkms [please apply msinit]", W);
     122       35959 : }
     123             : 
     124             : /** MODULAR TO SYM **/
     125             : 
     126             : /* q a t_FRAC or t_INT */
     127             : static GEN
     128        7630 : Q_log_init(ulong N, GEN q)
     129             : {
     130             :   long l, n;
     131             :   GEN Q;
     132             : 
     133        7630 :   q = gboundcf(q, 0);
     134        7630 :   l = lg(q);
     135        7630 :   Q = cgetg(l, t_VECSMALL);
     136        7630 :   Q[1] = 1;
     137        7630 :   for (n=2; n <l; n++) Q[n] = umodiu(gel(q,n), N);
     138       17087 :   for (n=3; n < l; n++)
     139        9457 :     Q[n] = Fl_add(Fl_mul(Q[n], Q[n-1], N), Q[n-2], N);
     140        7630 :   return Q;
     141             : }
     142             : 
     143             : /** INIT MODSYM STRUCTURE, WEIGHT 2 **/
     144             : 
     145             : /* num = [Gamma : Gamma_0(N)] = N * Prod_{p|N} (1+p^-1) */
     146             : static ulong
     147        1687 : count_Manin_symbols(ulong N, GEN P)
     148             : {
     149        1687 :   long i, l = lg(P);
     150        1687 :   ulong num = N;
     151        1687 :   for (i = 1; i < l; i++) { ulong p = P[i]; num *= p+1; num /= p; }
     152        1687 :   return num;
     153             : }
     154             : /* returns the list of "Manin symbols" (c,d) in (Z/NZ)^2, (c,d,N) = 1
     155             :  * generating H^1(X_0(N), Z) */
     156             : static GEN
     157        1687 : generatemsymbols(ulong N, ulong num, GEN divN)
     158             : {
     159        1687 :   GEN ret = cgetg(num+1, t_VEC);
     160        1687 :   ulong c, d, curn = 0;
     161             :   long i, l;
     162             :   /* generate Manin-symbols in two lists: */
     163             :   /* list 1: (c:1) for 0 <= c < N */
     164        1687 :   for (c = 0; c < N; c++) gel(ret, ++curn) = mkvecsmall2(c, 1);
     165        1687 :   if (N == 1) return ret;
     166             :   /* list 2: (c:d) with 1 <= c < N, c | N, 0 <= d < N, gcd(d,N) > 1, gcd(c,d)=1.
     167             :    * Furthermore, d != d0 (mod N/c) with c,d0 already in the list */
     168        1666 :   l = lg(divN) - 1;
     169             :   /* c = 1 first */
     170        1666 :   gel(ret, ++curn) = mkvecsmall2(1,0);
     171      152320 :   for (d = 2; d < N; d++)
     172      150654 :     if (ugcd(d,N) != 1UL)
     173       54586 :       gel(ret, ++curn) = mkvecsmall2(1,d);
     174             :   /* omit c = 1 (first) and c = N (last) */
     175        4781 :   for (i=2; i < l; i++)
     176             :   {
     177             :     ulong Novc, d0;
     178        3115 :     c = divN[i];
     179        3115 :     Novc = N / c;
     180       84672 :     for (d0 = 2; d0 <= Novc; d0++)
     181             :     {
     182       81557 :       ulong k, d = d0;
     183       81557 :       if (ugcd(d, Novc) == 1UL) continue;
     184      131257 :       for (k = 0; k < c; k++, d += Novc)
     185      119763 :         if (ugcd(c,d) == 1UL)
     186             :         {
     187       18592 :           gel(ret, ++curn) = mkvecsmall2(c,d);
     188       18592 :           break;
     189             :         }
     190             :     }
     191             :   }
     192        1666 :   if (curn != num) pari_err_BUG("generatemsymbols [wrong number of symbols]");
     193        1666 :   return ret;
     194             : }
     195             : 
     196             : static GEN
     197        1687 : inithashmsymbols(ulong N, GEN symbols)
     198             : {
     199        1687 :   GEN H = zerovec(N);
     200        1687 :   long k, l = lg(symbols);
     201             :   /* skip the (c:1), 0 <= c < N and (1:0) */
     202       74865 :   for (k=N+2; k < l; k++)
     203             :   {
     204       73178 :     GEN s = gel(symbols, k);
     205       73178 :     ulong c = s[1], d = s[2], Novc = N/c;
     206       73178 :     if (gel(H,c) == gen_0) gel(H,c) = const_vecsmall(Novc+1,0);
     207       73178 :     if (c != 1) { d %= Novc; if (!d) d = Novc; }
     208       73178 :     mael(H, c, d) = k;
     209             :   }
     210        1687 :   return H;
     211             : }
     212             : 
     213             : /** Helper functions for Sl2(Z) / Gamma_0(N) **/
     214             : /* M a 2x2 ZM in SL2(Z) */
     215             : static GEN
     216     1207451 : SL2_inv(GEN M)
     217             : {
     218     1207451 :   GEN a = gcoeff(M,1,1), b = gcoeff(M,1,2);
     219     1207451 :   GEN c = gcoeff(M,2,1), d = gcoeff(M,2,2);
     220     1207451 :   return mkmat22(d,negi(b), negi(c),a);
     221             : }
     222             : /* SL2_inv(M)[2] */
     223             : static GEN
     224        3409 : SL2_inv2(GEN M)
     225             : {
     226        3409 :   GEN a = gcoeff(M,1,1), b = gcoeff(M,1,2);
     227        3409 :   return mkcol2(negi(b),a);
     228             : }
     229             : /* M a 2x2 mat2 in SL2(Z) */
     230             : static GEN
     231      660681 : sl2_inv(GEN M)
     232             : {
     233      660681 :   long a=coeff(M,1,1), b=coeff(M,1,2), c=coeff(M,2,1), d=coeff(M,2,2);
     234      660681 :   return mkvec2(mkvecsmall2(d, -c), mkvecsmall2(-b, a));
     235             : }
     236             : /* Return the mat2 [a,b; c,d], not a zm to avoid GP problems */
     237             : static GEN
     238     1306879 : mat2(long a, long b, long c, long d)
     239     1306879 : { return mkvec2(mkvecsmall2(a,c), mkvecsmall2(b,d)); }
     240             : static GEN
     241      319683 : mat2_to_ZM(GEN M)
     242             : {
     243      319683 :   GEN A = gel(M,1), B = gel(M,2);
     244      319683 :   retmkmat2(mkcol2s(A[1],A[2]), mkcol2s(B[1],B[2]));
     245             : }
     246             : 
     247             : /* Input: a = 2-vector = path = {r/s,x/y}
     248             :  * Output: either [r,x;s,y] or [-r,x;-s,y], whichever has determinant > 0 */
     249             : static GEN
     250       80178 : path_to_ZM(GEN a)
     251             : {
     252       80178 :   GEN v = gel(a,1), w = gel(a,2);
     253       80178 :   long r = v[1], s = v[2], x = w[1], y = w[2];
     254       80178 :   if (cmpii(mulss(r,y), mulss(x,s)) < 0) { r = -r; s = -s; }
     255       80178 :   return mkmat22(stoi(r),stoi(x),stoi(s),stoi(y));
     256             : }
     257             : static GEN
     258      805315 : path_to_zm(GEN a)
     259             : {
     260      805315 :   GEN v = gel(a,1), w = gel(a,2);
     261      805315 :   long r = v[1], s = v[2], x = w[1], y = w[2];
     262      805315 :   if (cmpii(mulss(r,y), mulss(x,s)) < 0) { r = -r; s = -s; }
     263      805315 :   return mat2(r,x,s,y);
     264             : }
     265             : /* path from c1 to c2 */
     266             : static GEN
     267      463596 : mkpath(GEN c1, GEN c2) { return mat2(c1[1], c2[1], c1[2], c2[2]); }
     268             : static long
     269      654913 : cc(GEN M) { GEN v = gel(M,1); return v[2]; }
     270             : static long
     271      654913 : dd(GEN M) { GEN v = gel(M,2); return v[2]; }
     272             : 
     273             : /*Input: a,b = 2 paths, N = integer
     274             :  *Output: 1 if the a,b are \Gamma_0(N)-equivalent; 0 otherwise */
     275             : static int
     276       75306 : gamma_equiv(GEN a, GEN b, ulong N)
     277             : {
     278       75306 :   pari_sp av = avma;
     279       75306 :   GEN m = path_to_zm(a);
     280       75306 :   GEN n = path_to_zm(b);
     281       75306 :   GEN d = subii(mulss(cc(m),dd(n)), mulss(dd(m),cc(n)));
     282       75306 :   ulong res = umodiu(d, N);
     283       75306 :   avma = av; return res == 0;
     284             : }
     285             : /* Input: a,b = 2 paths that are \Gamma_0(N)-equivalent, N = integer
     286             :  * Output: M in \Gamma_0(N) such that Mb=a */
     287             : static GEN
     288       39725 : gamma_equiv_matrix(GEN a, GEN b)
     289             : {
     290       39725 :   GEN m = path_to_ZM(a);
     291       39725 :   GEN n = path_to_ZM(b);
     292       39725 :   return ZM_mul(m, SL2_inv(n));
     293             : }
     294             : 
     295             : /*************/
     296             : /* P^1(Z/NZ) */
     297             : /*************/
     298             : /* a != 0 in Z/NZ. Return v in (Z/NZ)^* such that av = gcd(a, N) (mod N)*/
     299             : static ulong
     300      354165 : Fl_inverse(ulong a, ulong N) { ulong g; return Fl_invgen(a,N,&g); }
     301             : 
     302             : /* Input: N = integer
     303             :  * Output: creates P^1(Z/NZ) = [symbols, H, N]
     304             :  *   symbols: list of vectors [x,y] that give a set of representatives
     305             :  *            of P^1(Z/NZ)
     306             :  *   H: an M by M grid whose value at the r,c-th place is the index of the
     307             :  *      "standard representative" equivalent to [r,c] occuring in the first
     308             :  *      list. If gcd(r,c,N) > 1 the grid has value 0. */
     309             : static GEN
     310        1687 : create_p1mod(ulong N)
     311             : {
     312        1687 :   GEN fa = factoru(N), div = divisorsu_fact(fa);
     313        1687 :   ulong i, nsym = count_Manin_symbols(N, gel(fa,1));
     314        1687 :   GEN symbols = generatemsymbols(N, nsym, div);
     315        1687 :   GEN H = inithashmsymbols(N,symbols);
     316        1687 :   GEN invsafe = cgetg(N, t_VECSMALL), inverse = cgetg(N, t_VECSMALL);
     317      154007 :   for (i = 1; i < N; i++)
     318             :   {
     319      152320 :     invsafe[i] = Fl_invsafe(i,N);
     320      152320 :     inverse[i] = Fl_inverse(i,N);
     321             :   }
     322        1687 :   return mkvecn(7, symbols, H, utoipos(N), fa, div, invsafe, inverse);
     323             : }
     324             : 
     325             : /* Let (c : d) in P1(Z/NZ).
     326             :  * If c = 0 return (0:1). If d = 0 return (1:0).
     327             :  * Else replace by (cu : du), where u in (Z/NZ)^* such that C := cu = gcd(c,N).
     328             :  * In create_p1mod(), (c : d) is represented by (C:D) where D = du (mod N/c)
     329             :  * is smallest such that gcd(C,D) = 1. Return (C : du mod N/c), which need
     330             :  * not belong to P1(Z/NZ) ! A second component du mod N/c = 0 is replaced by
     331             :  * N/c in this case to avoid problems with array indices */
     332             : static void
     333     2797228 : p1_std_form(long *pc, long *pd, GEN p1N)
     334             : {
     335     2797228 :   ulong N = p1N_get_N(p1N);
     336             :   ulong u;
     337     2797228 :   *pc = smodss(*pc, N); if (!*pc) { *pd = 1; return; }
     338     2716749 :   *pd = smodss(*pd, N); if (!*pd) { *pc = 1; return; }
     339     2677087 :   u = p1N_get_invsafe(p1N)[*pd];
     340     2677087 :   if (u) { *pc = Fl_mul(*pc,u,N); *pd = 1; return; } /* (d,N) = 1 */
     341             : 
     342      904162 :   u = p1N_get_inverse(p1N)[*pc];
     343      904162 :   if (u > 1) { *pc = Fl_mul(*pc,u,N); *pd = Fl_mul(*pd,u,N); }
     344             :   /* c | N */
     345      904162 :   if (*pc != 1) *pd %= (N / *pc);
     346      904162 :   if (!*pd) *pd = N / *pc;
     347             : }
     348             : 
     349             : /* Input: v = [x,y] = elt of P^1(Z/NZ) = class in Gamma_0(N) \ PSL2(Z)
     350             :  * Output: returns the index of the standard rep equivalent to v */
     351             : static long
     352     2797228 : p1_index(long x, long y, GEN p1N)
     353             : {
     354     2797228 :   ulong N = p1N_get_N(p1N);
     355     2797228 :   GEN H = p1N_get_hash(p1N);
     356             : 
     357     2797228 :   p1_std_form(&x, &y, p1N);
     358     2797228 :   if (y == 1) return x+1;
     359      943824 :   if (y == 0) return N+1;
     360      904162 :   if (mael(H,x,y) == 0) pari_err_BUG("p1_index");
     361      904162 :   return mael(H,x,y);
     362             : }
     363             : 
     364             : /* Cusps for \Gamma_0(N) */
     365             : 
     366             : /* \sum_{d | N} \phi(gcd(d, N/d)), using multiplicativity. fa = factor(N) */
     367             : ulong
     368        1736 : mfnumcuspsu_fact(GEN fa)
     369             : {
     370        1736 :   GEN P = gel(fa,1), E = gel(fa,2);
     371        1736 :   long i, l = lg(P);
     372        1736 :   ulong T = 1;
     373        4361 :   for (i = 1; i < l; i++)
     374             :   {
     375        2625 :     long e = E[i], e2 = e >> 1; /* floor(E[i] / 2) */
     376        2625 :     ulong p = P[i];
     377        2625 :     if (odd(e))
     378        2289 :       T *= 2 * upowuu(p, e2);
     379             :     else
     380         336 :       T *= (p+1) * upowuu(p, e2-1);
     381             :   }
     382        1736 :   return T;
     383             : }
     384             : ulong
     385           7 : mfnumcuspsu(ulong n)
     386             : {
     387           7 :   pari_sp av = avma;
     388           7 :   ulong t = mfnumcuspsu_fact( factoru(n) );
     389           7 :   avma = av; return t;
     390             : }
     391             : /* \sum_{d | N} \phi(gcd(d, N/d)), using multiplicativity. fa = factor(N) */
     392             : GEN
     393          14 : mfnumcusps_fact(GEN fa)
     394             : {
     395          14 :   GEN P = gel(fa,1), E = gel(fa,2), T = gen_1;
     396          14 :   long i, l = lg(P);
     397          35 :   for (i = 1; i < l; i++)
     398             :   {
     399          21 :     GEN p = gel(P,i), c;
     400          21 :     long e = itos(gel(E,i)), e2 = e >> 1; /* floor(E[i] / 2) */
     401          21 :     if (odd(e))
     402           0 :       c = shifti(powiu(p, e2), 1);
     403             :     else
     404          21 :       c = mulii(addiu(p,1), powiu(p, e2-1));
     405          21 :     T = T? mulii(T, c): c;
     406             :   }
     407          14 :   return T? T: gen_1;
     408             : }
     409             : GEN
     410          21 : mfnumcusps(GEN n)
     411             : {
     412          21 :   pari_sp av = avma;
     413          21 :   GEN F = check_arith_pos(n,"mfnumcusps");
     414          21 :   if (!F)
     415             :   {
     416          14 :     if (lgefint(n) == 3) return utoi( mfnumcuspsu(n[2]) );
     417           7 :     F = absZ_factor(n);
     418             :   }
     419          14 :   return gerepileuptoint(av, mfnumcusps_fact(F));
     420             : }
     421             : 
     422             : 
     423             : /* to each cusp in \Gamma_0(N) P1(Q), represented by p/q, we associate a
     424             :  * unique index. Canonical representative: (1:0) or (p:q) with q | N, q < N,
     425             :  * p defined modulo d := gcd(N/q,q), (p,d) = 1.
     426             :  * Return [[N, nbcusps], H, cusps]*/
     427             : static GEN
     428        1687 : inithashcusps(GEN p1N)
     429             : {
     430        1687 :   ulong N = p1N_get_N(p1N);
     431        1687 :   GEN div = p1N_get_div(p1N), H = zerovec(N+1);
     432        1687 :   long k, ind, l = lg(div), ncusp = mfnumcuspsu_fact(p1N_get_fa(p1N));
     433        1687 :   GEN cusps = cgetg(ncusp+1, t_VEC);
     434             : 
     435        1687 :   gel(H,1) = mkvecsmall2(0/*empty*/, 1/* first cusp: (1:0) */);
     436        1687 :   gel(cusps, 1) = mkvecsmall2(1,0);
     437        1687 :   ind = 2;
     438        6468 :   for (k=1; k < l-1; k++) /* l-1: remove q = N */
     439             :   {
     440        4781 :     ulong p, q = div[k], d = ugcd(q, N/q);
     441        4781 :     GEN h = const_vecsmall(d+1,0);
     442        4781 :     gel(H,q+1) = h ;
     443       12698 :     for (p = 0; p < d; p++)
     444        7917 :       if (ugcd(p,d) == 1)
     445             :       {
     446        6195 :         h[p+1] = ind;
     447        6195 :         gel(cusps, ind) = mkvecsmall2(p,q);
     448        6195 :         ind++;
     449             :       }
     450             :   }
     451        1687 :   return mkvec3(mkvecsmall2(N,ind-1), H, cusps);
     452             : }
     453             : /* c = [p,q], (p,q) = 1, return a canonical representative for
     454             :  * \Gamma_0(N)(p/q) */
     455             : static GEN
     456      203770 : cusp_std_form(GEN c, GEN S)
     457             : {
     458      203770 :   long p, N = gel(S,1)[1], q = smodss(c[2], N);
     459             :   ulong u, d;
     460      203770 :   if (q == 0) return mkvecsmall2(1, 0);
     461      201845 :   p = smodss(c[1], N);
     462      201845 :   u = Fl_inverse(q, N);
     463      201845 :   q = Fl_mul(q,u, N);
     464      201845 :   d = ugcd(q, N/q);
     465      201845 :   return mkvecsmall2(Fl_div(p % d,u % d, d), q);
     466             : }
     467             : /* c = [p,q], (p,q) = 1, return the index of the corresponding cusp.
     468             :  * S from inithashcusps */
     469             : static ulong
     470      203770 : cusp_index(GEN c, GEN S)
     471             : {
     472             :   long p, q;
     473      203770 :   GEN H = gel(S,2);
     474      203770 :   c = cusp_std_form(c, S);
     475      203770 :   p = c[1]; q = c[2];
     476      203770 :   if (!mael(H,q+1,p+1)) pari_err_BUG("cusp_index");
     477      203770 :   return mael(H,q+1,p+1);
     478             : }
     479             : 
     480             : /* M a square invertible ZM, return a ZM iM such that iM M = M iM = d.Id */
     481             : static GEN
     482        3052 : ZM_inv_denom(GEN M)
     483             : {
     484        3052 :   GEN diM, iM = ZM_inv(M, &diM);
     485        3052 :   return mkvec2(iM, diM);
     486             : }
     487             : /* return M^(-1) v, dinv = ZM_inv_denom(M) OR Qevproj_init(M) */
     488             : static GEN
     489      745115 : ZC_apply_dinv(GEN dinv, GEN v)
     490             : {
     491             :   GEN x, c, iM;
     492      745115 :   if (lg(dinv) == 3)
     493             :   {
     494      666421 :     iM = gel(dinv,1);
     495      666421 :     c = gel(dinv,2);
     496             :   }
     497             :   else
     498             :   { /* Qevproj_init */
     499       78694 :     iM = gel(dinv,2);
     500       78694 :     c = gel(dinv,3);
     501      157388 :     v = typ(v) == t_MAT? rowpermute(v, gel(dinv,4))
     502       78694 :                        : vecpermute(v, gel(dinv,4));
     503             :   }
     504      745115 :   x = RgM_RgC_mul(iM, v);
     505      745115 :   if (!isint1(c)) x = RgC_Rg_div(x, c);
     506      745115 :   return x;
     507             : }
     508             : 
     509             : /* M an n x d ZM of rank d (basis of a Q-subspace), n >= d.
     510             :  * Initialize a projector on M */
     511             : GEN
     512        4592 : Qevproj_init(GEN M)
     513             : {
     514             :   GEN v, perm, MM, iM, diM;
     515        4592 :   v = ZM_indexrank(M); perm = gel(v,1);
     516        4592 :   MM = rowpermute(M, perm); /* square invertible */
     517        4592 :   iM = ZM_inv(MM, &diM);
     518        4592 :   return mkvec4(M, iM, diM, perm);
     519             : }
     520             : 
     521             : /* same with typechecks */
     522             : static GEN
     523         714 : Qevproj_init0(GEN M)
     524             : {
     525         714 :   switch(typ(M))
     526             :   {
     527             :     case t_VEC:
     528         658 :       if (lg(M) == 5) return M;
     529           0 :       break;
     530             :     case t_COL:
     531          49 :       M = mkmat(M);/*fall through*/
     532             :     case t_MAT:
     533          56 :       M = Q_primpart(M);
     534          56 :       RgM_check_ZM(M,"Qevproj_init");
     535          56 :       return Qevproj_init(M);
     536             :   }
     537           0 :   pari_err_TYPE("Qevproj_init",M);
     538           0 :   return NULL;
     539             : }
     540             : 
     541             : /* T an n x n QM, pro = Qevproj_init(M), pro2 = Qevproj_init(M2); TM \subset M2.
     542             :  * Express these column vectors on M2's basis */
     543             : static GEN
     544        3724 : Qevproj_apply2(GEN T, GEN pro, GEN pro2)
     545             : {
     546        3724 :   GEN M = gel(pro,1), iM = gel(pro2,2), ciM = gel(pro2,3), perm = gel(pro2,4);
     547        3724 :   return RgM_Rg_div(RgM_mul(iM, RgM_mul(rowpermute(T,perm), M)), ciM);
     548             : }
     549             : /* T an n x n QM, stabilizing d-dimensional Q-vector space spanned by the
     550             :  * d columns of M, pro = Qevproj_init(M). Return dxd matrix of T acting on M */
     551             : GEN
     552        3122 : Qevproj_apply(GEN T, GEN pro) { return Qevproj_apply2(T, pro, pro); }
     553             : /* Qevproj_apply(T,pro)[,k] */
     554             : GEN
     555         819 : Qevproj_apply_vecei(GEN T, GEN pro, long k)
     556             : {
     557         819 :   GEN M = gel(pro,1), iM = gel(pro,2), ciM = gel(pro,3), perm = gel(pro,4);
     558         819 :   GEN v = RgM_RgC_mul(iM, RgM_RgC_mul(rowpermute(T,perm), gel(M,k)));
     559         819 :   return RgC_Rg_div(v, ciM);
     560             : }
     561             : 
     562             : static GEN
     563        1183 : QM_ker_r(GEN M) { return ZM_ker(Q_primpart(M)); }
     564             : static GEN
     565         959 : QM_image(GEN A)
     566             : {
     567         959 :   A = vec_Q_primpart(A);
     568         959 :   return vecpermute(A, ZM_indeximage(A));
     569             : }
     570             : 
     571             : static int
     572         420 : cmp_dim(void *E, GEN a, GEN b)
     573             : {
     574             :   long k;
     575             :   (void)E;
     576         420 :   a = gel(a,1);
     577         420 :   b = gel(b,1); k = lg(a)-lg(b);
     578         420 :   return k? ((k > 0)? 1: -1): 0;
     579             : }
     580             : 
     581             : /* FIXME: could use ZX_roots for deglim = 1 */
     582             : static GEN
     583         329 : ZX_factor_limit(GEN T, long deglim, long *pl)
     584             : {
     585         329 :   GEN fa = ZX_factor(T), P, E;
     586             :   long i, l;
     587         329 :   P = gel(fa,1); *pl = l = lg(P);
     588         329 :   if (deglim <= 0) return fa;
     589         224 :   E = gel(fa,2);
     590         567 :   for (i = 1; i < l; i++)
     591         406 :     if (degpol(gel(P,i)) > deglim) break;
     592         224 :   setlg(P,i);
     593         224 :   setlg(E,i); return fa;
     594             : }
     595             : 
     596             : /* Decompose the subspace H (Qevproj format) in simple subspaces.
     597             :  * Eg for H = msnew */
     598             : static GEN
     599         259 : mssplit_i(GEN W, GEN H, long deglim)
     600             : {
     601         259 :   ulong p, N = ms_get_N(W);
     602             :   long first, dim;
     603             :   forprime_t S;
     604         259 :   GEN T1 = NULL, T2 = NULL, V;
     605         259 :   dim = lg(gel(H,1))-1;
     606         259 :   V = vectrunc_init(dim+1);
     607         259 :   if (!dim) return V;
     608         252 :   (void)u_forprime_init(&S, 2, ULONG_MAX);
     609         252 :   vectrunc_append(V, H);
     610         252 :   first = 1; /* V[1..first-1] contains simple subspaces */
     611         630 :   while ((p = u_forprime_next(&S)))
     612             :   {
     613             :     GEN T;
     614             :     long j, lV;
     615         378 :     if (N % p == 0) continue;
     616         322 :     if (T1 && T2) {
     617          21 :       T = RgM_add(T1,T2);
     618          21 :       T2 = NULL;
     619             :     } else {
     620         301 :       T2 = T1;
     621         301 :       T1 = T = mshecke(W, p, NULL);
     622             :     }
     623         322 :     lV = lg(V);
     624         651 :     for (j = first; j < lV; j++)
     625             :     {
     626         329 :       pari_sp av = avma;
     627             :       long lP;
     628         329 :       GEN Vj = gel(V,j), P = gel(Vj,1);
     629         329 :       GEN TVj = Qevproj_apply(T, Vj); /* c T | V_j */
     630         329 :       GEN ch = QM_charpoly_ZX(TVj), fa = ZX_factor_limit(ch,deglim, &lP);
     631         329 :       GEN F = gel(fa, 1), E = gel(fa, 2);
     632         329 :       long k, lF = lg(F);
     633         329 :       if (lF == 2 && lP == 2)
     634             :       {
     635         336 :         if (isint1(gel(E,1)))
     636             :         { /* simple subspace */
     637         168 :           swap(gel(V,first), gel(V,j));
     638         168 :           first++;
     639             :         }
     640             :         else
     641           0 :           avma = av;
     642             :       }
     643         161 :       else if (lF == 1) /* discard V[j] */
     644           7 :       { swap(gel(V,j), gel(V,lg(V)-1)); setlg(V, lg(V)-1); }
     645             :       else
     646             :       { /* can split Vj */
     647             :         GEN pows;
     648         154 :         long D = 1;
     649         616 :         for (k = 1; k < lF; k++)
     650             :         {
     651         462 :           long d = degpol(gel(F,k));
     652         462 :           if (d > D) D = d;
     653             :         }
     654             :         /* remove V[j] */
     655         154 :         swap(gel(V,j), gel(V,lg(V)-1)); setlg(V, lg(V)-1);
     656         154 :         pows = RgM_powers(TVj, minss((long)2*sqrt((double)D), D));
     657         616 :         for (k = 1; k < lF; k++)
     658             :         {
     659         462 :           GEN f = gel(F,k);
     660         462 :           GEN K = QM_ker_r( RgX_RgMV_eval(f, pows)) ; /* Ker f(TVj) */
     661         462 :           GEN p = vec_Q_primpart( RgM_mul(P, K) );
     662         462 :           vectrunc_append(V, Qevproj_init(p));
     663         462 :           if (lg(K) == 2 || isint1(gel(E,k)))
     664             :           { /* simple subspace */
     665         385 :             swap(gel(V,first), gel(V, lg(V)-1));
     666         385 :             first++;
     667             :           }
     668             :         }
     669         154 :         if (j < first) j = first;
     670             :       }
     671             :     }
     672         322 :     if (first >= lg(V)) {
     673         252 :       gen_sort_inplace(V, NULL, cmp_dim, NULL);
     674         252 :       return V;
     675             :     }
     676             :   }
     677           0 :   pari_err_BUG("subspaces not found");
     678           0 :   return NULL;
     679             : }
     680             : GEN
     681         259 : mssplit(GEN W, GEN H, long deglim)
     682             : {
     683         259 :   pari_sp av = avma;
     684         259 :   checkms(W);
     685         259 :   if (!msk_get_sign(W))
     686           0 :     pari_err_DOMAIN("mssplit","abs(sign)","!=",gen_1,gen_0);
     687         259 :   if (!H) H = msnew(W);
     688         259 :   H = Qevproj_init0(H);
     689         259 :   return gerepilecopy(av, mssplit_i(W,H,deglim));
     690             : }
     691             : 
     692             : /* proV = Qevproj_init of a Hecke simple subspace, return [ a_n, n <= B ] */
     693             : static GEN
     694         245 : msqexpansion_i(GEN W, GEN proV, ulong B)
     695             : {
     696         245 :   ulong p, N = ms_get_N(W), sqrtB;
     697         245 :   long i, d, k = msk_get_weight(W);
     698             :   forprime_t S;
     699         245 :   GEN T1=NULL, T2=NULL, TV=NULL, ch=NULL, v, dTiv, Tiv, diM, iM, L;
     700         245 :   switch(B)
     701             :   {
     702           0 :     case 0: return cgetg(1,t_VEC);
     703           0 :     case 1: return mkvec(gen_1);
     704             :   }
     705         245 :   (void)u_forprime_init(&S, 2, ULONG_MAX);
     706         602 :   while ((p = u_forprime_next(&S)))
     707             :   {
     708             :     GEN T;
     709         357 :     if (N % p == 0) continue;
     710         266 :     if (T1 && T2)
     711             :     {
     712           0 :       T = RgM_add(T1,T2);
     713           0 :       T2 = NULL;
     714             :     }
     715             :     else
     716             :     {
     717         266 :       T2 = T1;
     718         266 :       T1 = T = mshecke(W, p, NULL);
     719             :     }
     720         266 :     TV = Qevproj_apply(T, proV); /* T | V */
     721         266 :     ch = QM_charpoly_ZX(TV);
     722         266 :     if (ZX_is_irred(ch)) break;
     723          21 :     ch = NULL;
     724             :   }
     725         245 :   if (!ch) pari_err_BUG("q-Expansion not found");
     726             :   /* T generates the Hecke algebra (acting on V) */
     727         245 :   d = degpol(ch);
     728         245 :   v = vec_ei(d, 1); /* take v = e_1 */
     729         245 :   Tiv = cgetg(d+1, t_MAT); /* Tiv[i] = T^(i-1)v */
     730         245 :   gel(Tiv, 1) = v;
     731         245 :   for (i = 2; i <= d; i++) gel(Tiv, i) = RgM_RgC_mul(TV, gel(Tiv,i-1));
     732         245 :   Tiv = Q_remove_denom(Tiv, &dTiv);
     733         245 :   iM = ZM_inv(Tiv, &diM);
     734         245 :   if (dTiv) diM = gdiv(diM, dTiv);
     735         245 :   L = const_vec(B,NULL);
     736         245 :   sqrtB = usqrt(B);
     737         245 :   gel(L,1) = d > 1? mkpolmod(gen_1,ch): gen_1;
     738        2471 :   for (p = 2; p <= B; p++)
     739             :   {
     740        2226 :     pari_sp av = avma;
     741             :     GEN T, u, Tv, ap, P;
     742             :     ulong m;
     743        2226 :     if (gel(L,p)) continue;  /* p not prime */
     744         819 :     T = mshecke(W, p, NULL);
     745         819 :     Tv = Qevproj_apply_vecei(T, proV, 1); /* Tp.v */
     746             :     /* Write Tp.v = \sum u_i T^i v */
     747         819 :     u = RgC_Rg_div(RgM_RgC_mul(iM, Tv), diM);
     748         819 :     ap = gerepilecopy(av, RgV_to_RgX(u, 0));
     749         819 :     if (d > 1)
     750         399 :       ap = mkpolmod(ap,ch);
     751             :     else
     752         420 :       ap = simplify_shallow(ap);
     753         819 :     gel(L,p) = ap;
     754         819 :     if (!(N % p))
     755             :     { /* p divides the level */
     756         147 :       ulong C = B/p;
     757         546 :       for (m=1; m<=C; m++)
     758         399 :         if (gel(L,m)) gel(L,m*p) = gmul(gel(L,m), ap);
     759         147 :       continue;
     760             :     }
     761         672 :     P = powuu(p,k-1);
     762         672 :     if (p <= sqrtB) {
     763         119 :       ulong pj, oldpj = 1;
     764         546 :       for (pj = p; pj <= B; oldpj=pj, pj *= p)
     765             :       {
     766         427 :         GEN apj = (pj==p)? ap
     767         427 :                          : gsub(gmul(ap,gel(L,oldpj)), gmul(P,gel(L,oldpj/p)));
     768         427 :         gel(L,pj) = apj;
     769        3136 :         for (m = B/pj; m > 1; m--)
     770        2709 :           if (gel(L,m) && m%p) gel(L,m*pj) = gmul(gel(L,m), apj);
     771             :       }
     772             :     } else {
     773         553 :       gel(L,p) = ap;
     774        1092 :       for (m = B/p; m > 1; m--)
     775         539 :         if (gel(L,m)) gel(L,m*p) = gmul(gel(L,m), ap);
     776             :     }
     777             :   }
     778         245 :   return L;
     779             : }
     780             : GEN
     781         245 : msqexpansion(GEN W, GEN proV, ulong B)
     782             : {
     783         245 :   pari_sp av = avma;
     784         245 :   checkms(W);
     785         245 :   proV = Qevproj_init0(proV);
     786         245 :   return gerepilecopy(av, msqexpansion_i(W,proV,B));
     787             : }
     788             : 
     789             : static GEN
     790         266 : Qevproj_apply0(GEN T, GEN pro)
     791             : {
     792         266 :   GEN iM = gel(pro,2), perm = gel(pro,4);
     793         266 :   return vec_Q_primpart(ZM_mul(iM, rowpermute(T,perm)));
     794             : }
     795             : /* T a ZC or ZM */
     796             : GEN
     797        1057 : Qevproj_down(GEN T, GEN pro)
     798             : {
     799        1057 :   GEN iM = gel(pro,2), ciM = gel(pro,3), perm = gel(pro,4);
     800        1057 :   if (typ(T) == t_COL)
     801        1057 :     return RgC_Rg_div(ZM_ZC_mul(iM, vecpermute(T,perm)), ciM);
     802             :   else
     803           0 :     return RgM_Rg_div(ZM_mul(iM, rowpermute(T,perm)), ciM);
     804             : }
     805             : 
     806             : static GEN
     807         357 : Qevproj_star(GEN W, GEN H)
     808             : {
     809         357 :   long s = msk_get_sign(W);
     810         357 :   if (s)
     811             :   { /* project on +/- component */
     812         266 :     GEN A = RgM_mul(msk_get_star(W), H);
     813         266 :     A = (s > 0)? gadd(A, H): gsub(A, H);
     814             :     /* Im(star + sign) = Ker(star - sign) */
     815         266 :     H = QM_image(A);
     816         266 :     H = Qevproj_apply0(H, msk_get_starproj(W));
     817             :   }
     818         357 :   return H;
     819             : }
     820             : 
     821             : static GEN
     822        2604 : Tp_matrices(ulong p)
     823             : {
     824        2604 :   GEN v = cgetg(p+2, t_VEC);
     825             :   ulong i;
     826        2604 :   for (i = 1; i <= p; i++) gel(v,i) = mat2(1, i-1, 0, p);
     827        2604 :   gel(v,i) = mat2(p, 0, 0, 1);
     828        2604 :   return v;
     829             : }
     830             : static GEN
     831         973 : Up_matrices(ulong p)
     832             : {
     833         973 :   GEN v = cgetg(p+1, t_VEC);
     834             :   ulong i;
     835         973 :   for (i = 1; i <= p; i++) gel(v,i) = mat2(1, i-1, 0, p);
     836         973 :   return v;
     837             : }
     838             : 
     839             : /* M = N/p. Classes of Gamma_0(M) / Gamma_O(N) when p | M */
     840             : static GEN
     841         168 : NP_matrices(ulong M, ulong p)
     842             : {
     843         168 :   GEN v = cgetg(p+1, t_VEC);
     844             :   ulong i;
     845         168 :   for (i = 1; i <= p; i++) gel(v,i) = mat2(1, 0, (i-1)*M, 1);
     846         168 :   return v;
     847             : }
     848             : /* M = N/p. Extra class of Gamma_0(M) / Gamma_O(N) when p \nmid M */
     849             : static GEN
     850          84 : NP_matrix_extra(ulong M, ulong p)
     851             : {
     852          84 :   long w,z, d = cbezout(p, -M, &w, &z);
     853          84 :   if (d != 1) return NULL;
     854          84 :   return mat2(w,z,M,p);
     855             : }
     856             : static GEN
     857          98 : WQ_matrix(long N, long Q)
     858             : {
     859          98 :   long w,z, d = cbezout(Q, N/Q, &w, &z);
     860          98 :   if (d != 1) return NULL;
     861          98 :   return mat2(Q,1,-N*z,Q*w);
     862             : }
     863             : 
     864             : GEN
     865         280 : msnew(GEN W)
     866             : {
     867         280 :   pari_sp av = avma;
     868         280 :   GEN S = mscuspidal(W, 0);
     869         280 :   ulong N = ms_get_N(W);
     870         280 :   long s = msk_get_sign(W), k = msk_get_weight(W);
     871         280 :   if (N > 1 && (!uisprime(N) || (k == 12 || k > 14)))
     872             :   {
     873         105 :     GEN p1N = ms_get_p1N(W), P = gel(p1N_get_fa(p1N), 1);
     874         105 :     long i, nP = lg(P)-1;
     875         105 :     GEN v = cgetg(2*nP + 1, t_COL);
     876         105 :     S = gel(S,1); /* Q basis */
     877         273 :     for (i = 1; i <= nP; i++)
     878             :     {
     879         168 :       pari_sp av = avma, av2;
     880         168 :       long M = N/P[i];
     881         168 :       GEN T1,Td, Wi = mskinit(M, k, s);
     882         168 :       GEN v1 = NP_matrices(M, P[i]);
     883         168 :       GEN vd = Up_matrices(P[i]);
     884             :       /* p^2 \nmid N */
     885         168 :       if (M % P[i])
     886             :       {
     887          84 :         v1 = shallowconcat(v1, mkvec(NP_matrix_extra(M,P[i])));
     888          84 :         vd = shallowconcat(vd, mkvec(WQ_matrix(N,P[i])));
     889             :       }
     890         168 :       T1 = getMorphism(W, Wi, v1);
     891         168 :       Td = getMorphism(W, Wi, vd);
     892         168 :       if (s)
     893             :       {
     894         154 :         T1 = Qevproj_apply2(T1, msk_get_starproj(W), msk_get_starproj(Wi));
     895         154 :         Td = Qevproj_apply2(Td, msk_get_starproj(W), msk_get_starproj(Wi));
     896             :       }
     897         168 :       av2 = avma;
     898         168 :       T1 = RgM_mul(T1,S);
     899         168 :       Td = RgM_mul(Td,S);  /* multiply by S = restrict to mscusp */
     900         168 :       gerepileallsp(av, av2, 2, &T1, &Td);
     901         168 :       gel(v,2*i-1) = T1;
     902         168 :       gel(v,2*i)   = Td;
     903             :     }
     904         105 :     S = ZM_mul(S, QM_ker_r(matconcat(v))); /* Snew */
     905         105 :     S = Qevproj_init(vec_Q_primpart(S));
     906             :   }
     907         280 :   return gerepilecopy(av, S);
     908             : }
     909             : 
     910             : /* Solve the Manin relations for a congruence subgroup \Gamma by constructing
     911             :  * a well-formed fundamental domain for the action of \Gamma on upper half
     912             :  * space. See
     913             :  * Pollack and Stevens, Overconvergent modular symbols and p-adic L-functions
     914             :  * Annales scientifiques de l'ENS 44, fascicule 1 (2011), 1-42
     915             :  * http://math.bu.edu/people/rpollack/Papers/Overconvergent_modular_symbols_and_padic_Lfunctions.pdf
     916             :  *
     917             :  * FIXME: Implemented for \Gamma = \Gamma_0(N) only. */
     918             : 
     919             : #if 0 /* Pollack-Stevens shift their paths so as to solve equations of the
     920             :          form f(z+1) - f(z) = g. We don't (to avoid mistakes) so we will
     921             :          have to solve eqs of the form f(z-1) - f(z) = g */
     922             : /* c = a/b; as a t_VECSMALL [a,b]; return c-1 as a t_VECSMALL */
     923             : static GEN
     924             : Shift_left_cusp(GEN c) { long a=c[1], b=c[2]; return mkvecsmall2(a - b, b); }
     925             : /* c = a/b; as a t_VECSMALL [a,b]; return c+1 as a t_VECSMALL */
     926             : static GEN
     927             : Shift_right_cusp(GEN c) { long a=c[1], b=c[2]; return mkvecsmall2(a + b, b); }
     928             : /*Input: path = [r,s] (thought of as a geodesic between these points)
     929             :  *Output: The path shifted by one to the left, i.e. [r-1,s-1] */
     930             : static GEN
     931             : Shift_left(GEN path)
     932             : {
     933             :   GEN r = gel(path,1), s = gel(path,2);
     934             :   return mkvec2(Shift_left_cusp(r), Shift_left_cusp(s)); }
     935             : /*Input: path = [r,s] (thought of as a geodesic between these points)
     936             :  *Output: The path shifted by one to the right, i.e. [r+1,s+1] */
     937             : GEN
     938             : Shift_right(GEN path)
     939             : {
     940             :   GEN r = gel(path,1), s = gel(path,2);
     941             :   return mkvec2(Shift_right_cusp(r), Shift_right_cusp(s)); }
     942             : #endif
     943             : 
     944             : /* linked lists */
     945             : typedef struct list_t { GEN data; struct list_t *next; } list_t;
     946             : static list_t *
     947       77700 : list_new(GEN x)
     948             : {
     949       77700 :   list_t *L = (list_t*)stack_malloc(sizeof(list_t));
     950       77700 :   L->data = x;
     951       77700 :   L->next = NULL; return L;
     952             : }
     953             : static void
     954       76034 : list_insert(list_t *L, GEN x)
     955             : {
     956       76034 :   list_t *l = list_new(x);
     957       76034 :   l->next = L->next;
     958       76034 :   L->next = l;
     959       76034 : }
     960             : 
     961             : /*Input: N > 1, p1N = P^1(Z/NZ)
     962             :  *Output: a connected fundamental domain for the action of \Gamma_0(N) on
     963             :  *  upper half space.  When \Gamma_0(N) is torsion free, the domain has the
     964             :  *  property that all of its vertices are cusps.  When \Gamma_0(N) has
     965             :  *  three-torsion, 2 extra triangles need to be added.
     966             :  *
     967             :  * The domain is constructed by beginning with the triangle with vertices 0,1
     968             :  * and oo.  Each adjacent triangle is successively tested to see if it contains
     969             :  * points not \Gamma_0(N) equivalent to some point in our region.  If a
     970             :  * triangle contains new points, it is added to the region.  This process is
     971             :  * continued until the region can no longer be extended (and still be a
     972             :  * fundamental domain) by added an adjacent triangle.  The list of cusps
     973             :  * between 0 and 1 are then returned
     974             :  *
     975             :  * Precisely, the function returns a list such that the elements of the list
     976             :  * with odd index are the cusps in increasing order.  The even elements of the
     977             :  * list are either an "x" or a "t".  A "t" represents that there is an element
     978             :  * of order three such that its fixed point is in the triangle directly
     979             :  * adjacent to the our region with vertices given by the cusp before and after
     980             :  * the "t".  The "x" represents that this is not the case. */
     981             : enum { type_X, type_DO /* ? */, type_T };
     982             : static GEN
     983        1666 : form_list_of_cusps(ulong N, GEN p1N)
     984             : {
     985        1666 :   pari_sp av = avma;
     986        1666 :   long i, position, nbC = 2;
     987             :   GEN v, L;
     988             :   list_t *C, *c;
     989             :   /* Let t be the index of a class in PSL2(Z) / \Gamma in our fixed enumeration
     990             :    * v[t] != 0 iff it is the class of z tau^r for z a previous alpha_i
     991             :    * or beta_i.
     992             :    * For \Gamma = \Gamma_0(N), the enumeration is given by p1_index.
     993             :    * We write cl(gamma) = the class of gamma mod \Gamma */
     994        1666 :   v = const_vecsmall(p1_size(p1N), 0);
     995        1666 :   i = p1_index( 0, 1, p1N); v[i] = 1;
     996        1666 :   i = p1_index( 1,-1, p1N); v[i] = 2;
     997        1666 :   i = p1_index(-1, 0, p1N); v[i] = 3;
     998             :   /* the value is unused [debugging]: what matters is whether it is != 0 */
     999        1666 :   position = 4;
    1000             :   /* at this point, Fund = R, v contains the classes of Id, tau, tau^2 */
    1001             : 
    1002        1666 :   C  = list_new(mkvecsmall3(0,1, type_X));
    1003        1666 :   list_insert(C, mkvecsmall3(1,1,type_DO));
    1004             :   /* C is a list of triples[a,b,t], where c = a/b is a cusp, and t is the type
    1005             :    * of the path between c and the PREVIOUS cusp in the list, coded as
    1006             :    *   type_DO = "?", type_X = "x", type_T = "t"
    1007             :    * Initially, C = [0/1,"?",1/1]; */
    1008             : 
    1009             :   /* loop through the current set of cusps C and check to see if more cusps
    1010             :    * should be added */
    1011             :   for (;;)
    1012             :   {
    1013        9149 :     int done = 1;
    1014      369530 :     for (c = C; c; c = c->next)
    1015             :     {
    1016             :       GEN cusp1, cusp2, gam;
    1017             :       long pos, b1, b2, b;
    1018             : 
    1019      369530 :       if (!c->next) break;
    1020      360381 :       cusp1 = c->data; /* = a1/b1 */
    1021      360381 :       cusp2 = (c->next)->data; /* = a2/b2 */
    1022      360381 :       if (cusp2[3] != type_DO) continue;
    1023             : 
    1024             :       /* gam (oo -> 0) = (cusp2 -> cusp1), gam in PSL2(Z) */
    1025      150402 :       gam = path_to_zm(mkpath(cusp2, cusp1)); /* = [a2,a1;b2,b1] */
    1026             :       /* we have normalized the cusp representation so that a1 b2 - a2 b1 = 1 */
    1027      150402 :       b1 = coeff(gam,2,1); b2 = coeff(gam,2,2);
    1028             :       /* gam.1  = (a1 + a2) / (b1 + b2) */
    1029      150402 :       b = b1 + b2;
    1030             :       /* Determine whether the adjacent triangle *below* (cusp1->cusp2)
    1031             :        * should be added */
    1032      150402 :       pos = p1_index(b1,b2, p1N); /* did we see cl(gam) before ? */
    1033      150402 :       if (v[pos])
    1034       75306 :         cusp2[3] = type_X; /* NO */
    1035             :       else
    1036             :       { /* YES */
    1037             :         ulong B1, B2;
    1038       75096 :         v[pos] = position;
    1039       75096 :         i = p1_index(-(b1+b2), b1, p1N); v[i] = position+1;
    1040       75096 :         i = p1_index(b2, -(b1+b2), p1N); v[i] = position+2;
    1041             :         /* add cl(gam), cl(gam*TAU), cl(gam*TAU^2) to v */
    1042       75096 :         position += 3;
    1043             :         /* gam tau gam^(-1) in \Gamma ? */
    1044       75096 :         B1 = smodss(b1, N);
    1045       75096 :         B2 = smodss(b2, N);
    1046       75096 :         if ((Fl_sqr(B2,N) + Fl_sqr(B1,N) + Fl_mul(B1,B2,N)) % N == 0)
    1047         728 :           cusp2[3] = type_T;
    1048             :         else
    1049             :         {
    1050       74368 :           long a1 = coeff(gam, 1,1), a2 = coeff(gam, 1,2);
    1051       74368 :           long a = a1 + a2; /* gcd(a,b) = 1 */
    1052       74368 :           list_insert(c, mkvecsmall3(a,b,type_DO));
    1053       74368 :           c = c->next;
    1054       74368 :           nbC++;
    1055       74368 :           done = 0;
    1056             :         }
    1057             :       }
    1058             :     }
    1059        9149 :     if (done) break;
    1060        7483 :   }
    1061        1666 :   L = cgetg(nbC+1, t_VEC); i = 1;
    1062        1666 :   for (c = C; c; c = c->next) gel(L,i++) = c->data;
    1063        1666 :   return gerepilecopy(av, L);
    1064             : }
    1065             : 
    1066             : /* W an msN. M in PSL2(Z). Return index of M in P1^(Z/NZ) = Gamma0(N) \ PSL2(Z),
    1067             :  * and M0 in Gamma_0(N) such that M = M0 * M', where M' = chosen
    1068             :  * section( PSL2(Z) -> P1^(Z/NZ) ). */
    1069             : static GEN
    1070      499023 : Gamma0N_decompose(GEN W, GEN M, long *index)
    1071             : {
    1072      499023 :   GEN p1N = msN_get_p1N(W), W3 = gel(W,3), section = msN_get_section(W);
    1073             :   GEN A;
    1074      499023 :   ulong N = p1N_get_N(p1N);
    1075      499023 :   ulong c = umodiu(gcoeff(M,2,1), N);
    1076      499023 :   ulong d = umodiu(gcoeff(M,2,2), N);
    1077      499023 :   long s, ind = p1_index(c, d, p1N); /* as an elt of P1(Z/NZ) */
    1078      499023 :   *index = W3[ind]; /* as an elt of F, E2, ... */
    1079      499023 :   M = ZM_zm_mul(M, sl2_inv(gel(section,ind)));
    1080             :   /* normalize mod +/-Id */
    1081      499023 :   A = gcoeff(M,1,1);
    1082      499023 :   s = signe(A);
    1083      499023 :   if (s < 0)
    1084      237363 :     M = ZM_neg(M);
    1085      261660 :   else if (!s)
    1086             :   {
    1087         322 :     GEN C = gcoeff(M,2,1);
    1088         322 :     if (signe(C) < 0) M = ZM_neg(M);
    1089             :   }
    1090      499023 :   return M;
    1091             : }
    1092             : /* W an msN; as above for a path. Return [[ind], M] */
    1093             : static GEN
    1094      158732 : path_Gamma0N_decompose(GEN W, GEN path)
    1095             : {
    1096      158732 :   GEN p1N = msN_get_p1N(W);
    1097      158732 :   GEN p1index_to_ind = gel(W,3);
    1098      158732 :   GEN section = msN_get_section(W);
    1099      158732 :   GEN M = path_to_zm(path);
    1100      158732 :   long p1index = p1_index(cc(M), dd(M), p1N);
    1101      158732 :   long ind = p1index_to_ind[p1index];
    1102      158732 :   GEN M0 = ZM_zm_mul(mat2_to_ZM(M), sl2_inv(gel(section,p1index)));
    1103      158732 :   return mkvec2(mkvecsmall(ind), M0);
    1104             : }
    1105             : 
    1106             : /*Form generators of H_1(X_0(N),{cusps},Z)
    1107             : *
    1108             : *Input: N = integer > 1, p1N = P^1(Z/NZ)
    1109             : *Output: [cusp_list,E,F,T2,T3,E1] where
    1110             : *  cusps_list = list of cusps describing fundamental domain of
    1111             : *    \Gamma_0(N).
    1112             : *  E = list of paths in the boundary of the fundamental domains and oriented
    1113             : *    clockwise such that they do not contain a point
    1114             : *    fixed by an element of order 2 and they are not an edge of a
    1115             : *    triangle containing a fixed point of an element of order 3
    1116             : *  F = list of paths in the interior of the domain with each
    1117             : *    orientation appearing separately
    1118             : * T2 = list of paths in the boundary of domain containing a point fixed
    1119             : *    by an element of order 2 (oriented clockwise)
    1120             : * T3 = list of paths in the boundard of domain which are the edges of
    1121             : *    some triangle containing a fixed point of a matrix of order 3 (both
    1122             : *    orientations appear)
    1123             : * E1 = a sublist of E such that every path in E is \Gamma_0(N)-equivalent to
    1124             : *    either an element of E1 or the flip (reversed orientation) of an element
    1125             : *    of E1.
    1126             : * (Elements of T2 are \Gamma_0(N)-equivalent to their own flip.)
    1127             : *
    1128             : * sec = a list from 1..#p1N of matrices describing a section of the map
    1129             : *   SL_2(Z) to P^1(Z/NZ) given by [a,b;c,d]-->[c,d].
    1130             : *   Given our fixed enumeration of P^1(Z/NZ), the j-th element of the list
    1131             : *   represents the image of the j-th element of P^1(Z/NZ) under the section. */
    1132             : 
    1133             : /* insert path in set T */
    1134             : static void
    1135      228830 : set_insert(hashtable *T, GEN path)
    1136      228830 : { hash_insert(T, path,  (void*)(T->nb + 1)); }
    1137             : 
    1138             : static GEN
    1139       14994 : hash_to_vec(hashtable *h)
    1140             : {
    1141       14994 :   GEN v = cgetg(h->nb + 1, t_VEC);
    1142             :   ulong i;
    1143     1931398 :   for (i = 0; i < h->len; i++)
    1144             :   {
    1145     1916404 :     hashentry *e = h->table[i];
    1146     4211914 :     while (e)
    1147             :     {
    1148      379106 :       GEN key = (GEN)e->key;
    1149      379106 :       long index = (long)e->val;
    1150      379106 :       gel(v, index) = key;
    1151      379106 :       e = e->next;
    1152             :     }
    1153             :   }
    1154       14994 :   return v;
    1155             : }
    1156             : 
    1157             : static long
    1158      116739 : path_to_p1_index(GEN path, GEN p1N)
    1159             : {
    1160      116739 :   GEN M = path_to_zm(path);
    1161      116739 :   return p1_index(cc(M), dd(M), p1N);
    1162             : }
    1163             : 
    1164             : /* Pollack-Stevens sets */
    1165             : typedef struct PS_sets_t {
    1166             :   hashtable *F, *T2, *T31, *T32, *E1, *E2;
    1167             :   GEN E2fromE1, stdE1;
    1168             : } PS_sets_t;
    1169             : 
    1170             : static hashtable *
    1171       13734 : set_init(long max)
    1172       13734 : { return hash_create(max, (ulong(*)(void*))&hash_GEN,
    1173             :                           (int(*)(void*,void*))&gidentical, 1); }
    1174             : /* T = E2fromE1[i] = [c, gamma] */
    1175             : static ulong
    1176     5426141 : E2fromE1_c(GEN T) { return itou(gel(T,1)); }
    1177             : static GEN
    1178      580237 : E2fromE1_Zgamma(GEN T) { return gel(T,2); }
    1179             : static GEN
    1180       38913 : E2fromE1_gamma(GEN T) { return gcoeff(gel(T,2),1,1); }
    1181             : 
    1182             : static void
    1183       77826 : insert_E(GEN path, PS_sets_t *S, GEN p1N)
    1184             : {
    1185       77826 :   GEN rev = vecreverse(path);
    1186       77826 :   long std = path_to_p1_index(rev, p1N);
    1187       77826 :   GEN v = gel(S->stdE1, std);
    1188       77826 :   if (v)
    1189             :   { /* [s, p1], where E1[s] is the path p1 = vecreverse(path) mod \Gamma */
    1190       38913 :     GEN gamma, p1 = gel(v,2);
    1191       38913 :     long r, s = itou(gel(v,1));
    1192             : 
    1193       38913 :     set_insert(S->E2, path);
    1194       38913 :     r = S->E2->nb;
    1195       38913 :     if (gel(S->E2fromE1, r) != gen_0) pari_err_BUG("insert_E");
    1196             : 
    1197       38913 :     gamma = gamma_equiv_matrix(rev, p1);
    1198             :     /* reverse(E2[r]) = gamma * E1[s] */
    1199       38913 :     gel(S->E2fromE1, r) = mkvec2(utoipos(s), to_famat_shallow(gamma,gen_m1));
    1200             :   }
    1201             :   else
    1202             :   {
    1203       38913 :     set_insert(S->E1, path);
    1204       38913 :     std = path_to_p1_index(path, p1N);
    1205       38913 :     gel(S->stdE1, std) = mkvec2(utoipos(S->E1->nb), path);
    1206             :   }
    1207       77826 : }
    1208             : 
    1209             : static GEN
    1210        6664 : cusp_infinity(void) { return mkvecsmall2(1,0); }
    1211             : 
    1212             : static void
    1213        1666 : form_E_F_T(ulong N, GEN p1N, GEN *pC, PS_sets_t *S)
    1214             : {
    1215        1666 :   GEN C, cusp_list = form_list_of_cusps(N, p1N);
    1216        1666 :   long nbgen = lg(cusp_list)-1, nbmanin = p1_size(p1N), r, s, i;
    1217             :   hashtable *F, *T2, *T31, *T32, *E1, *E2;
    1218             : 
    1219        1666 :   *pC = C = cgetg(nbgen+1, t_VEC);
    1220       79366 :   for (i = 1; i <= nbgen; i++)
    1221             :   {
    1222       77700 :     GEN c = gel(cusp_list,i);
    1223       77700 :     gel(C,i) = mkvecsmall2(c[1], c[2]);
    1224             :   }
    1225        1666 :   S->F  = F  = set_init(nbmanin);
    1226        1666 :   S->E1 = E1 = set_init(nbgen);
    1227        1666 :   S->E2 = E2 = set_init(nbgen);
    1228        1666 :   S->T2 = T2 = set_init(nbgen);
    1229        1666 :   S->T31 = T31 = set_init(nbgen);
    1230        1666 :   S->T32 = T32 = set_init(nbgen);
    1231             : 
    1232             :   /* T31 represents the three torsion paths going from left to right */
    1233             :   /* T32 represents the three torsion paths going from right to left */
    1234       77700 :   for (r = 1; r < nbgen; r++)
    1235             :   {
    1236       76034 :     GEN c2 = gel(cusp_list,r+1);
    1237       76034 :     if (c2[3] == type_T)
    1238             :     {
    1239         728 :       GEN c1 = gel(cusp_list,r), path = mkpath(c1,c2), path2 = vecreverse(path);
    1240         728 :       set_insert(T31, path);
    1241         728 :       set_insert(T32, path2);
    1242             :     }
    1243             :   }
    1244             : 
    1245             :   /* to record relations between E2 and E1 */
    1246        1666 :   S->E2fromE1 = zerovec(nbgen);
    1247        1666 :   S->stdE1 = const_vec(nbmanin, NULL);
    1248             : 
    1249             :   /* Assumption later: path [oo,0] is E1[1], path [1,oo] is E2[1] */
    1250             :   {
    1251        1666 :     GEN oo = cusp_infinity();
    1252        1666 :     GEN p1 = mkpath(oo, mkvecsmall2(0,1)); /* [oo, 0] */
    1253        1666 :     GEN p2 = mkpath(mkvecsmall2(1,1), oo); /* [1, oo] */
    1254        1666 :     insert_E(p1, S, p1N);
    1255        1666 :     insert_E(p2, S, p1N);
    1256             :   }
    1257             : 
    1258       77700 :   for (r = 1; r < nbgen; r++)
    1259             :   {
    1260       76034 :     GEN c1 = gel(cusp_list,r);
    1261    18062940 :     for (s = r+1; s <= nbgen; s++)
    1262             :     {
    1263    17986906 :       pari_sp av = avma;
    1264    17986906 :       GEN c2 = gel(cusp_list,s), path;
    1265    17986906 :       GEN d = subii(mulss(c1[1],c2[2]), mulss(c1[2],c2[1]));
    1266    17986906 :       avma = av;
    1267    17986906 :       if (!is_pm1(d)) continue;
    1268             : 
    1269      150402 :       path = mkpath(c1,c2);
    1270      150402 :       if (r+1 == s)
    1271             :       {
    1272       76034 :         GEN w = path;
    1273       76034 :         ulong hash = T31->hash(w); /* T31, T32 use the same hash function */
    1274       76034 :         if (!hash_search2(T31, w, hash) && !hash_search2(T32, w, hash))
    1275             :         {
    1276       75306 :           if (gamma_equiv(path, vecreverse(path), N))
    1277         812 :             set_insert(T2, path);
    1278             :           else
    1279       74494 :             insert_E(path, S, p1N);
    1280             :         }
    1281             :       } else {
    1282       74368 :         set_insert(F, mkvec2(path, mkvecsmall2(r,s)));
    1283       74368 :         set_insert(F, mkvec2(vecreverse(path), mkvecsmall2(s,r)));
    1284             :       }
    1285             :     }
    1286             :   }
    1287        1666 :   setlg(S->E2fromE1, E2->nb+1);
    1288        1666 : }
    1289             : 
    1290             : /* v = \sum n_i g_i, g_i in Sl(2,Z), return \sum n_i g_i^(-1) */
    1291             : static GEN
    1292      846321 : ZSl2_star(GEN v)
    1293             : {
    1294             :   long i, l;
    1295             :   GEN w, G;
    1296      846321 :   if (typ(v) == t_INT) return v;
    1297      846321 :   G = gel(v,1);
    1298      846321 :   w = cgetg_copy(G, &l);
    1299     2016707 :   for (i = 1; i < l; i++)
    1300             :   {
    1301     1170386 :     GEN g = gel(G,i);
    1302     1170386 :     if (typ(g) == t_MAT) g = SL2_inv(g);
    1303     1170386 :     gel(w,i) = g;
    1304             :   }
    1305      846321 :   return ZG_normalize(mkmat2(w, gel(v,2)));
    1306             : }
    1307             : 
    1308             : /* Input: h = set of unimodular paths, p1N = P^1(Z/NZ) = Gamma_0(N)\PSL2(Z)
    1309             :  * Output: Each path is converted to a matrix and then an element of P^1(Z/NZ)
    1310             :  * Append the matrix to W[12], append the index that represents
    1311             :  * these elements of P^1 (the classes mod Gamma_0(N) via our fixed
    1312             :  * enumeration to W[2]. */
    1313             : static void
    1314        9996 : paths_decompose(GEN W, hashtable *h, int flag)
    1315             : {
    1316        9996 :   GEN p1N = ms_get_p1N(W), section = ms_get_section(W);
    1317        9996 :   GEN v = hash_to_vec(h);
    1318        9996 :   long i, l = lg(v);
    1319      238826 :   for (i = 1; i < l; i++)
    1320             :   {
    1321      228830 :     GEN e = gel(v,i);
    1322      228830 :     GEN M = path_to_zm(flag? gel(e,1): e);
    1323      228830 :     long index = p1_index(cc(M), dd(M), p1N);
    1324      228830 :     vecsmalltrunc_append(gel(W,2), index);
    1325      228830 :     gel(section, index) = M;
    1326             :   }
    1327        9996 : }
    1328             : static void
    1329        1666 : fill_W2_W12(GEN W, PS_sets_t *S)
    1330             : {
    1331        1666 :   GEN p1N = msN_get_p1N(W);
    1332        1666 :   long n = p1_size(p1N);
    1333        1666 :   gel(W, 2) = vecsmalltrunc_init(n+1);
    1334        1666 :   gel(W,12) = cgetg(n+1, t_VEC);
    1335             :   /* F contains [path, [index cusp1, index cusp2]]. Others contain paths only */
    1336        1666 :   paths_decompose(W, S->F, 1);
    1337        1666 :   paths_decompose(W, S->E2, 0);
    1338        1666 :   paths_decompose(W, S->T32, 0);
    1339        1666 :   paths_decompose(W, S->E1, 0);
    1340        1666 :   paths_decompose(W, S->T2, 0);
    1341        1666 :   paths_decompose(W, S->T31, 0);
    1342        1666 : }
    1343             : 
    1344             : /* x t_VECSMALL, corresponds to a map x(i) = j, where 1 <= j <= max for all i
    1345             :  * Return y s.t. y[j] = i or 0 (not in image) */
    1346             : static GEN
    1347        3332 : reverse_list(GEN x, long max)
    1348             : {
    1349        3332 :   GEN y = const_vecsmall(max, 0);
    1350        3332 :   long r, lx = lg(x);
    1351        3332 :   for (r = 1; r < lx; r++) y[ x[r] ] = r;
    1352        3332 :   return y;
    1353             : }
    1354             : 
    1355             : /* go from C[a] to C[b]; return the indices of paths
    1356             :  * E.g. if a < b
    1357             :  *   (C[a]->C[a+1], C[a+1]->C[a+2], ... C[b-1]->C[b])
    1358             :  * (else reverse direction)
    1359             :  * = b - a paths */
    1360             : static GEN
    1361      144942 : F_indices(GEN W, long a, long b)
    1362             : {
    1363      144942 :   GEN v = cgetg(labs(b-a) + 1, t_VEC);
    1364      144942 :   long s, k = 1;
    1365      144942 :   if (a < b) {
    1366       72471 :     GEN index_forward = gel(W,13);
    1367       72471 :     for (s = a; s < b; s++) gel(v,k++) = gel(index_forward,s);
    1368             :   } else {
    1369       72471 :     GEN index_backward = gel(W,14);
    1370       72471 :     for (s = a; s > b; s--) gel(v,k++) = gel(index_backward,s);
    1371             :   }
    1372      144942 :   return v;
    1373             : }
    1374             : /* go from C[a] to C[b] via oo; return the indices of paths
    1375             :  * E.g. if a < b
    1376             :  *   (C[a]->C[a-1], ... C[2]->C[1],
    1377             :  *    C[1]->oo, oo-> C[end],
    1378             :  *    C[end]->C[end-1], ... C[b+1]->C[b])
    1379             :  *  a-1 + 2 + end-(b+1)+1 = end - b + a + 1 paths  */
    1380             : static GEN
    1381        3794 : F_indices_oo(GEN W, long end, long a, long b)
    1382             : {
    1383        3794 :   GEN index_oo = gel(W,15);
    1384        3794 :   GEN v = cgetg(end-labs(b-a)+1 + 1, t_VEC);
    1385        3794 :   long s, k = 1;
    1386             : 
    1387        3794 :   if (a < b) {
    1388        1897 :     GEN index_backward = gel(W,14);
    1389        1897 :     for (s = a; s > 1; s--) gel(v,k++) = gel(index_backward,s);
    1390        1897 :     gel(v,k++) = gel(index_backward,1); /* C[1] -> oo */
    1391        1897 :     gel(v,k++) = gel(index_oo,2); /* oo -> C[end] */
    1392        1897 :     for (s = end; s > b; s--) gel(v,k++) = gel(index_backward,s);
    1393             :   } else {
    1394        1897 :     GEN index_forward = gel(W,13);
    1395        1897 :     for (s = a; s < end; s++) gel(v,k++) = gel(index_forward,s);
    1396        1897 :     gel(v,k++) = gel(index_forward,end); /* C[end] -> oo */
    1397        1897 :     gel(v,k++) = gel(index_oo,1); /* oo -> C[1] */
    1398        1897 :     for (s = 1; s < b; s++) gel(v,k++) = gel(index_forward,s);
    1399             :   }
    1400        3794 :   return v;
    1401             : }
    1402             : /* index of oo -> C[1], oo -> C[end] */
    1403             : static GEN
    1404        1666 : indices_oo(GEN W, GEN C)
    1405             : {
    1406        1666 :   long end = lg(C)-1;
    1407        1666 :   GEN w, v = cgetg(2+1, t_VEC), oo = cusp_infinity();
    1408        1666 :   w = mkpath(oo, gel(C,1)); /* oo -> C[1]=0 */
    1409        1666 :   gel(v,1) = path_Gamma0N_decompose(W, w);
    1410        1666 :   w = mkpath(oo, gel(C,end)); /* oo -> C[end]=1 */
    1411        1666 :   gel(v,2) = path_Gamma0N_decompose(W, w);
    1412        1666 :   return v;
    1413             : }
    1414             : 
    1415             : /* index of C[1]->C[2], C[2]->C[3], ... C[end-1]->C[end], C[end]->oo
    1416             :  * Recall that C[1] = 0, C[end] = 1 */
    1417             : static GEN
    1418        1666 : indices_forward(GEN W, GEN C)
    1419             : {
    1420        1666 :   long s, k = 1, end = lg(C)-1;
    1421        1666 :   GEN v = cgetg(end+1, t_VEC);
    1422       79366 :   for (s = 1; s <= end; s++)
    1423             :   {
    1424       77700 :     GEN w = mkpath(gel(C,s), s == end? cusp_infinity(): gel(C,s+1));
    1425       77700 :     gel(v,k++) = path_Gamma0N_decompose(W, w);
    1426             :   }
    1427        1666 :   return v;
    1428             : }
    1429             : /* index of C[1]->oo, C[2]->C[1], ... C[end]->C[end-1] */
    1430             : static GEN
    1431        1666 : indices_backward(GEN W, GEN C)
    1432             : {
    1433        1666 :   long s, k = 1, end = lg(C)-1;
    1434        1666 :   GEN v = cgetg(end+1, t_VEC);
    1435       79366 :   for (s = 1; s <= end; s++)
    1436             :   {
    1437       77700 :     GEN w = mkpath(gel(C,s), s == 1? cusp_infinity(): gel(C,s-1));
    1438       77700 :     gel(v,k++) = path_Gamma0N_decompose(W, w);
    1439             :   }
    1440        1666 :   return v;
    1441             : }
    1442             : 
    1443             : /*[0,-1;1,-1]*/
    1444             : static GEN
    1445        1708 : mkTAU()
    1446        1708 : { return mkmat22(gen_0,gen_m1, gen_1,gen_m1); }
    1447             : /* S */
    1448             : static GEN
    1449          42 : mkS()
    1450          42 : { return mkmat22(gen_0,gen_1, gen_m1,gen_0); }
    1451             : /* N = integer > 1. Returns data describing Delta_0 = Z[P^1(Q)]_0 seen as
    1452             :  * a Gamma_0(N) - module. */
    1453             : static GEN
    1454        1687 : msinit_N(ulong N)
    1455             : {
    1456             :   GEN p1N, C, vecF, vecT2, vecT31, TAU, W, W2, singlerel, annT2, annT31;
    1457             :   GEN F_index;
    1458             :   ulong r, s, width;
    1459             :   long nball, nbgen, nbp1N;
    1460             :   hashtable *F, *T2, *T31, *T32, *E1, *E2;
    1461             :   PS_sets_t S;
    1462             : 
    1463        1687 :   W = zerovec(16);
    1464        1687 :   gel(W,1) = p1N = create_p1mod(N);
    1465        1687 :   gel(W,16)= inithashcusps(p1N);
    1466        1687 :   TAU = mkTAU();
    1467        1687 :   if (N == 1)
    1468             :   {
    1469          21 :     gel(W,5) = mkvecsmall(1);
    1470             :     /* cheat because sets are not disjoint if N=1 */
    1471          21 :     gel(W,11) = mkvecsmall5(0, 0, 1, 1, 2);
    1472          21 :     gel(W,12) = mkvec(mat2(1,0,0,1));
    1473          21 :     gel(W,8) = mkvec( mkmat22(gen_1,gen_1, mkS(),gen_1) );
    1474          21 :     gel(W,9) = mkvec( mkmat2(mkcol3(gen_1,TAU,ZM_sqr(TAU)),
    1475             :                              mkcol3(gen_1,gen_1,gen_1)) );
    1476          21 :     return W;
    1477             :   }
    1478        1666 :   nbp1N = p1_size(p1N);
    1479        1666 :   form_E_F_T(N,p1N, &C, &S);
    1480        1666 :   E1  = S.E1;
    1481        1666 :   E2  = S.E2;
    1482        1666 :   T31 = S.T31;
    1483        1666 :   T32 = S.T32;
    1484        1666 :   F   = S.F;
    1485        1666 :   T2  = S.T2;
    1486        1666 :   nbgen = lg(C)-1;
    1487             : 
    1488             :  /* Put our paths in the order: F,E2,T32,E1,T2,T31
    1489             :   * W2[j] associates to the j-th element of this list its index in P1. */
    1490        1666 :   fill_W2_W12(W, &S);
    1491        1666 :   W2 = gel(W, 2);
    1492        1666 :   nball = lg(W2)-1;
    1493        1666 :   gel(W,3) = reverse_list(W2, nbp1N);
    1494        1666 :   gel(W,5) = vecslice(gel(W,2), F->nb + E2->nb + T32->nb + 1, nball);
    1495        1666 :   gel(W,4) = reverse_list(gel(W,5), nbp1N);
    1496        1666 :   gel(W,13) = indices_forward(W, C);
    1497        1666 :   gel(W,14) = indices_backward(W, C);
    1498        1666 :   gel(W,15) = indices_oo(W, C);
    1499        8330 :   gel(W,11) = mkvecsmall5(F->nb,
    1500        1666 :                           F->nb + E2->nb,
    1501        1666 :                           F->nb + E2->nb + T32->nb,
    1502        1666 :                           F->nb + E2->nb + T32->nb + E1->nb,
    1503        1666 :                           F->nb + E2->nb + T32->nb + E1->nb + T2->nb);
    1504             :   /* relations between T32 and T31 [not stored!]
    1505             :    * T32[i] = - T31[i] */
    1506             : 
    1507             :   /* relations of F */
    1508        1666 :   width = E1->nb + T2->nb + T31->nb;
    1509             :   /* F_index[r] = [index_1, ..., index_k], where index_i is the p1_index()
    1510             :    * of the elementary unimodular path between 2 consecutive cusps
    1511             :    * [in E1,E2,T2,T31 or T32] */
    1512        1666 :   F_index = cgetg(F->nb+1, t_VEC);
    1513        1666 :   vecF = hash_to_vec(F);
    1514      150402 :   for (r = 1; r <= F->nb; r++)
    1515             :   {
    1516      148736 :     GEN w = gel(gel(vecF,r), 2);
    1517      148736 :     long a = w[1], b = w[2], d = labs(b - a);
    1518             :     /* c1 = cusp_list[a],  c2 = cusp_list[b], ci != oo */
    1519      297472 :     gel(F_index,r) = (nbgen-d >= d-1)? F_indices(W, a,b)
    1520      148736 :                                      : F_indices_oo(W, lg(C)-1,a,b);
    1521             :   }
    1522             : 
    1523        1666 :   singlerel = cgetg(width+1, t_VEC);
    1524             :   /* form the single boundary relation */
    1525       40579 :   for (s = 1; s <= E2->nb; s++)
    1526             :   { /* reverse(E2[s]) = gamma * E1[c] */
    1527       38913 :     GEN T = gel(S.E2fromE1,s), gamma = E2fromE1_gamma(T);
    1528       38913 :     gel(singlerel, E2fromE1_c(T)) = mkmat22(gen_1,gen_1, gamma,gen_m1);
    1529             :   }
    1530        1666 :   for (r = E1->nb + 1; r <= width; r++) gel(singlerel, r) = gen_1;
    1531             : 
    1532             :   /* form the 2-torsion relations */
    1533        1666 :   annT2 = cgetg(T2->nb+1, t_VEC);
    1534        1666 :   vecT2 = hash_to_vec(T2);
    1535        2478 :   for (r = 1; r <= T2->nb; r++)
    1536             :   {
    1537         812 :     GEN w = gel(vecT2,r);
    1538         812 :     GEN gamma = gamma_equiv_matrix(vecreverse(w), w);
    1539         812 :     gel(annT2, r) = mkmat22(gen_1,gen_1, gamma,gen_1);
    1540             :   }
    1541             : 
    1542             :   /* form the 3-torsion relations */
    1543        1666 :   annT31 = cgetg(T31->nb+1, t_VEC);
    1544        1666 :   vecT31 = hash_to_vec(T31);
    1545        2394 :   for (r = 1; r <= T31->nb; r++)
    1546             :   {
    1547         728 :     GEN M = path_to_ZM( vecreverse(gel(vecT31,r)) );
    1548         728 :     GEN gamma = ZM_mul(ZM_mul(M, TAU), SL2_inv(M));
    1549         728 :     gel(annT31, r) = mkmat2(mkcol3(gen_1,gamma,ZM_sqr(gamma)),
    1550             :                             mkcol3(gen_1,gen_1,gen_1));
    1551             :   }
    1552        1666 :   gel(W,6) = F_index;
    1553        1666 :   gel(W,7) = S.E2fromE1;
    1554        1666 :   gel(W,8) = annT2;
    1555        1666 :   gel(W,9) = annT31;
    1556        1666 :   gel(W,10)= singlerel;
    1557        1666 :   return W;
    1558             : }
    1559             : static GEN
    1560         112 : cusp_to_P1Q(GEN c) { return c[2]? gdivgs(stoi(c[1]), c[2]): mkoo(); }
    1561             : static GEN
    1562          21 : mspathgens_i(GEN W)
    1563             : {
    1564             :   GEN R, r, g, section, gen, annT2, annT31;
    1565             :   long i, l;
    1566          21 :   checkms(W); W = get_msN(W);
    1567          21 :   section = msN_get_section(W);
    1568          21 :   gen = ms_get_genindex(W);
    1569          21 :   l = lg(gen);
    1570          21 :   g = cgetg(l,t_VEC);
    1571          77 :   for (i = 1; i < l; i++)
    1572             :   {
    1573          56 :     GEN p = gel(section,gen[i]);
    1574          56 :     gel(g,i) = mkvec2(cusp_to_P1Q(gel(p,1)), cusp_to_P1Q(gel(p,2)));
    1575             :   }
    1576          21 :   annT2 = msN_get_annT2(W);
    1577          21 :   annT31= msN_get_annT31(W);
    1578          21 :   if (ms_get_N(W) == 1)
    1579             :   {
    1580           7 :     R = cgetg(3, t_VEC);
    1581           7 :     gel(R,1) = mkvec( mkvec2(gel(annT2,1), gen_1) );
    1582           7 :     gel(R,2) = mkvec( mkvec2(gel(annT31,1), gen_1) );
    1583             :   }
    1584             :   else
    1585             :   {
    1586          14 :     GEN singlerel = msN_get_singlerel(W);
    1587          14 :     long j, nbT2 = lg(annT2)-1, nbT31 = lg(annT31)-1, nbE1 = ms_get_nbE1(W);
    1588          14 :     R = cgetg(nbT2+nbT31+2, t_VEC);
    1589          14 :     l = lg(singlerel);
    1590          14 :     r = cgetg(l, t_VEC);
    1591          42 :     for (i = 1; i <= nbE1; i++)
    1592          28 :       gel(r,i) = mkvec2(gel(singlerel, i), utoi(i));
    1593          35 :     for (; i < l; i++)
    1594          21 :       gel(r,i) = mkvec2(gen_1, utoi(i));
    1595          14 :     gel(R,1) = r; j = 2;
    1596          35 :     for (i = 1; i <= nbT2; i++,j++)
    1597          21 :       gel(R,j) = mkvec( mkvec2(gel(annT2,i), utoi(i + nbE1)) );
    1598          14 :     for (i = 1; i <= nbT31; i++,j++)
    1599           0 :       gel(R,j) = mkvec( mkvec2(gel(annT31,i), utoi(i + nbE1 + nbT2)) );
    1600             :   }
    1601          21 :   return mkvec2(g,R);
    1602             : }
    1603             : GEN
    1604          21 : mspathgens(GEN W)
    1605             : {
    1606          21 :   pari_sp av = avma;
    1607          21 :   return gerepilecopy(av, mspathgens_i(W));
    1608             : }
    1609             : /* Modular symbols in weight k: Hom_Gamma(Delta, Q[x,y]_{k-2}) */
    1610             : /* A symbol phi is represented by the {phi(g_i)}, {phi(g'_i)}, {phi(g''_i)}
    1611             :  * where the {g_i, g'_i, g''_i} are the Z[\Gamma]-generators of Delta,
    1612             :  * g_i corresponds to E1, g'_i to T2, g''_i to T31.
    1613             :  */
    1614             : 
    1615             : /* FIXME: export. T^1, ..., T^n */
    1616             : static GEN
    1617      702520 : RgX_powers(GEN T, long n)
    1618             : {
    1619      702520 :   GEN v = cgetg(n+1, t_VEC);
    1620             :   long i;
    1621      702520 :   gel(v, 1) = T;
    1622      702520 :   for (i = 1; i < n; i++) gel(v,i+1) = RgX_mul(gel(v,i), T);
    1623      702520 :   return v;
    1624             : }
    1625             : 
    1626             : /* g = [a,b;c,d] a mat2. Return (X^{k-2} | g)(X,Y)[X = 1]. */
    1627             : static GEN
    1628        2926 : voo_act_Gl2Q(GEN g, long k)
    1629             : {
    1630        2926 :   GEN mc = stoi(-coeff(g,2,1)), d = stoi(coeff(g,2,2));
    1631        2926 :   return RgX_to_RgC(gpowgs(deg1pol_shallow(mc, d, 0), k-2), k-1);
    1632             : }
    1633             : 
    1634             : struct m_act {
    1635             :   long dim, k, p;
    1636             :   GEN q;
    1637             :   GEN(*act)(struct m_act *,GEN);
    1638             : };
    1639             : 
    1640             : /* g = [a,b;c,d]. Return (P | g)(X,Y)[X = 1] = P(dX - cY, -b X + aY)[X = 1],
    1641             :  * for P = X^{k-2}, X^{k-3}Y, ..., Y^{k-2} */
    1642             : GEN
    1643      351260 : RgX_act_Gl2Q(GEN g, long k)
    1644             : {
    1645             :   GEN a,b,c,d, V1,V2,V;
    1646             :   long i;
    1647      351260 :   if (k == 2) return matid(1);
    1648      351260 :   a = gcoeff(g,1,1); b = gcoeff(g,1,2);
    1649      351260 :   c = gcoeff(g,2,1); d = gcoeff(g,2,2);
    1650      351260 :   V1 = RgX_powers(deg1pol_shallow(gneg(c), d, 0), k-2); /* d - c Y */
    1651      351260 :   V2 = RgX_powers(deg1pol_shallow(a, gneg(b), 0), k-2); /*-b + a Y */
    1652      351260 :   V = cgetg(k, t_MAT);
    1653      351260 :   gel(V,1)   = RgX_to_RgC(gel(V1, k-2), k-1);
    1654      820708 :   for (i = 1; i < k-2; i++)
    1655             :   {
    1656      469448 :     GEN v1 = gel(V1, k-2-i); /* (d-cY)^(k-2-i) */
    1657      469448 :     GEN v2 = gel(V2, i); /* (-b+aY)^i */
    1658      469448 :     gel(V,i+1) = RgX_to_RgC(RgX_mul(v1,v2), k-1);
    1659             :   }
    1660      351260 :   gel(V,k-1) = RgX_to_RgC(gel(V2, k-2), k-1);
    1661      351260 :   return V; /* V[i+1] = X^i | g */
    1662             : }
    1663             : /* z in Z[Gl2(Q)], return the matrix of z acting on V */
    1664             : static GEN
    1665      600824 : act_ZGl2Q(GEN z, struct m_act *T, hashtable *H)
    1666             : {
    1667      600824 :   GEN S = NULL, G, E;
    1668             :   pari_sp av;
    1669             :   long l, j;
    1670             :   /* paranoia: should not occur */
    1671      600824 :   if (typ(z) == t_INT) return scalarmat_shallow(z, T->dim);
    1672      600824 :   G = gel(z,1); l = lg(G);
    1673      600824 :   E = gel(z,2); av = avma;
    1674     1771210 :   for (j = 1; j < l; j++)
    1675             :   {
    1676     1170386 :     GEN M, g = gel(G,j), n = gel(E,j);
    1677     1170386 :     if (typ(g) == t_INT) /* = 1 */
    1678        3927 :       M = n; /* n*Id_dim */
    1679             :     else
    1680             :     { /*search in H succeeds because of preload*/
    1681     1166459 :       M = H? (GEN)hash_search(H,g)->val: T->act(T,g);
    1682     1166459 :       if (is_pm1(n))
    1683     1158934 :       { if (signe(n) < 0) M = RgM_neg(M); }
    1684             :       else
    1685        7525 :         M = RgM_Rg_mul(M, n);
    1686             :     }
    1687     1170386 :     if (!S) { S = M; continue; }
    1688      569562 :     S = gadd(S, M);
    1689      569562 :     if (gc_needed(av,1))
    1690             :     {
    1691           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"act_ZGl2Q, j = %ld",j);
    1692           0 :       S = gerepileupto(av, S);
    1693             :     }
    1694             :   }
    1695      600824 :   return gerepilecopy(av, S);
    1696             : }
    1697             : static GEN
    1698      351113 : _RgX_act_Gl2Q(struct m_act *S, GEN z) { return RgX_act_Gl2Q(z, S->k); }
    1699             : /* acting on (X^{k-2},...,Y^{k-2}) */
    1700             : GEN
    1701       60718 : RgX_act_ZGl2Q(GEN z, long k)
    1702             : {
    1703             :   struct m_act T;
    1704       60718 :   T.k = k;
    1705       60718 :   T.dim = k-1;
    1706       60718 :   T.act=&_RgX_act_Gl2Q;
    1707       60718 :   return act_ZGl2Q(z, &T, NULL);
    1708             : }
    1709             : 
    1710             : /* First pass, identify matrices in Sl_2 to convert to operators;
    1711             :  * insert operators in hashtable. This allows GC in act_ZGl2Q */
    1712             : static void
    1713     1070832 : hash_preload(GEN M, struct m_act *S, hashtable *H)
    1714             : {
    1715     1070832 :   if (typ(M) != t_INT)
    1716             :   {
    1717     1070832 :     ulong h = H->hash(M);
    1718     1070832 :     hashentry *e = hash_search2(H, M, h);
    1719     1070832 :     if (!e) hash_insert2(H, M, S->act(S,M), h);
    1720             :   }
    1721     1070832 : }
    1722             : /* z a sparse operator */
    1723             : static void
    1724      540092 : hash_vecpreload(GEN z, struct m_act *S, hashtable *H)
    1725             : {
    1726      540092 :   GEN G = gel(z,1);
    1727      540092 :   long i, l = lg(G);
    1728      540092 :   for (i = 1; i < l; i++) hash_preload(gel(G,i), S, H);
    1729      540092 : }
    1730             : static void
    1731       40817 : ZGl2QC_preload(struct m_act *S, GEN v, hashtable *H)
    1732             : {
    1733       40817 :   GEN val = gel(v,2);
    1734       40817 :   long i, l = lg(val);
    1735       40817 :   for (i = 1; i < l; i++) hash_vecpreload(gel(val,i), S, H);
    1736       40817 : }
    1737             : /* Given a sparse vector of elements in Z[G], convert it to a (sparse) vector
    1738             :  * of operators on V (given by t_MAT) */
    1739             : static void
    1740       40831 : ZGl2QC_to_act(struct m_act *S, GEN v, hashtable *H)
    1741             : {
    1742       40831 :   GEN val = gel(v,2);
    1743       40831 :   long i, l = lg(val);
    1744       40831 :   for (i = 1; i < l; i++) gel(val,i) = act_ZGl2Q(gel(val,i), S, H);
    1745       40831 : }
    1746             : 
    1747             : /* For all V[i] in Z[\Gamma], find the P such that  P . V[i]^* = 0;
    1748             :  * write P in basis X^{k-2}, ..., Y^{k-2} */
    1749             : static GEN
    1750        1246 : ZGV_tors(GEN V, long k)
    1751             : {
    1752        1246 :   long i, l = lg(V);
    1753        1246 :   GEN v = cgetg(l, t_VEC);
    1754        1750 :   for (i = 1; i < l; i++)
    1755             :   {
    1756         504 :     GEN a = ZSl2_star(gel(V,i));
    1757         504 :     gel(v,i) = ZM_ker(RgX_act_ZGl2Q(a,k));
    1758             :   }
    1759        1246 :   return v;
    1760             : }
    1761             : 
    1762             : static long
    1763    12397784 : set_from_index(GEN W11, long i)
    1764             : {
    1765    12397784 :   if (i <= W11[1]) return 1;
    1766    11194533 :   if (i <= W11[2]) return 2;
    1767     5813017 :   if (i <= W11[3]) return 3;
    1768     5804274 :   if (i <= W11[4]) return 4;
    1769       31367 :   if (i <= W11[5]) return 5;
    1770        8463 :   return 6;
    1771             : }
    1772             : 
    1773             : /* det M = 1 */
    1774             : static void
    1775     1536703 : treat_index(GEN W, GEN M, long index, GEN v)
    1776             : {
    1777     1536703 :   GEN W11 = gel(W,11);
    1778     1536703 :   long shift = W11[3]; /* #F + #E2 + T32 */
    1779     1536703 :   switch(set_from_index(W11, index))
    1780             :   {
    1781             :     case 1: /*F*/
    1782             :     {
    1783      251363 :       GEN F_index = gel(W,6), ind = gel(F_index, index);
    1784      251363 :       long j, l = lg(ind);
    1785     1289043 :       for (j = 1; j < l; j++)
    1786             :       {
    1787     1037680 :         GEN IND = gel(ind,j), M0 = gel(IND,2);
    1788     1037680 :         long index = mael(IND,1,1);
    1789     1037680 :         treat_index(W, ZM_mul(M,M0), index, v);
    1790             :       }
    1791      251363 :       break;
    1792             :     }
    1793             : 
    1794             :     case 2: /*E2, E2[r] + gamma * E1[s] = 0 */
    1795             :     {
    1796      580237 :       long r = index - W11[1];
    1797      580237 :       GEN z = gel(msN_get_E2fromE1(W), r);
    1798             : 
    1799      580237 :       index = E2fromE1_c(z);
    1800      580237 :       M = G_ZG_mul(M, E2fromE1_Zgamma(z)); /* M * (-gamma) */
    1801      580237 :       gel(v, index) = ZG_add(gel(v, index), M);
    1802      580237 :       break;
    1803             :     }
    1804             : 
    1805             :     case 3: /*T32, T32[i] = -T31[i] */
    1806             :     {
    1807        6006 :       long T3shift = W11[5] - W11[2]; /* #T32 + #E1 + #T2 */
    1808        6006 :       index += T3shift;
    1809        6006 :       index -= shift;
    1810        6006 :       gel(v, index) = ZG_add(gel(v, index), to_famat_shallow(M,gen_m1));
    1811        6006 :       break;
    1812             :     }
    1813             :     default: /*E1,T2,T31*/
    1814      699097 :       index -= shift;
    1815      699097 :       gel(v, index) = ZG_add(gel(v, index), to_famat_shallow(M,gen_1));
    1816      699097 :       break;
    1817             :   }
    1818     1536703 : }
    1819             : static void
    1820    10861081 : treat_index_trivial(GEN v, GEN W, long index)
    1821             : {
    1822    10861081 :   GEN W11 = gel(W,11);
    1823    10861081 :   long shift = W11[3]; /* #F + #E2 + T32 */
    1824    10861081 :   switch(set_from_index(W11, index))
    1825             :   {
    1826             :     case 1: /*F*/
    1827             :     {
    1828      951888 :       GEN F_index = gel(W,6), ind = gel(F_index, index);
    1829      951888 :       long j, l = lg(ind);
    1830    10324657 :       for (j = 1; j < l; j++)
    1831             :       {
    1832     9372769 :         GEN IND = gel(ind,j);
    1833     9372769 :         treat_index_trivial(v, W, mael(IND,1,1));
    1834             :       }
    1835      951888 :       break;
    1836             :     }
    1837             : 
    1838             :     case 2: /*E2, E2[r] + gamma * E1[s] = 0 */
    1839             :     {
    1840     4801279 :       long r = index - W11[1];
    1841     4801279 :       long s = E2fromE1_c(gel(msN_get_E2fromE1(W), r));
    1842     4801279 :       v[s]--;
    1843     4801279 :       break;
    1844             :     }
    1845             : 
    1846             :     case 3: case 5: case 6: /*T32,T2,T31*/
    1847       15463 :       break;
    1848             : 
    1849             :     case 4: /*E1*/
    1850     5092451 :       v[index-shift]++;
    1851     5092451 :       break;
    1852             :   }
    1853    10861081 : }
    1854             : 
    1855             : static GEN
    1856      178542 : M2_log(GEN W, GEN M)
    1857             : {
    1858      178542 :   GEN a = gcoeff(M,1,1), b = gcoeff(M,1,2);
    1859      178542 :   GEN c = gcoeff(M,2,1), d = gcoeff(M,2,2);
    1860             :   GEN  u, v, D, V;
    1861             :   long index, s;
    1862             : 
    1863      178542 :   W = get_msN(W);
    1864      178542 :   V = zerovec(ms_get_nbgen(W));
    1865             : 
    1866      178542 :   D = subii(mulii(a,d), mulii(b,c));
    1867      178542 :   s = signe(D);
    1868      178542 :   if (!s) return V;
    1869      177198 :   if (is_pm1(D))
    1870             :   { /* shortcut, no need to apply Manin's trick */
    1871       63532 :     if (s < 0) { b = negi(b); d = negi(d); }
    1872       63532 :     M = Gamma0N_decompose(W, mkmat22(a,b, c,d), &index);
    1873       63532 :     treat_index(W, M, index, V);
    1874             :   }
    1875             :   else
    1876             :   {
    1877             :     GEN U, B, P, Q, PQ, C1,C2;
    1878             :     long i, l;
    1879      113666 :     (void)bezout(a,c,&u,&v);
    1880      113666 :     B = addii(mulii(b,u), mulii(d,v));
    1881             :     /* [u,v;-c,a] [a,b; c,d] = [1,B; 0,D], i.e. M = U [1,B;0,D] */
    1882      113666 :     U = mkmat22(a,negi(v), c,u);
    1883             : 
    1884             :     /* {1/0 -> B/D} as \sum g_i, g_i unimodular paths */
    1885      113666 :     PQ = ZV_allpnqn( gboundcf(gdiv(B,D), 0) );
    1886      113666 :     P = gel(PQ,1); l = lg(P);
    1887      113666 :     Q = gel(PQ,2);
    1888      113666 :     C1 = gel(U,1);
    1889      549157 :     for (i = 1; i < l; i++, C1 = C2)
    1890             :     {
    1891             :       GEN M;
    1892      435491 :       C2 = ZM_ZC_mul(U, mkcol2(gel(P,i), gel(Q,i)));
    1893      435491 :       if (!odd(i)) C1 = ZC_neg(C1);
    1894      435491 :       M = Gamma0N_decompose(W, mkmat2(C1,C2), &index);
    1895      435491 :       treat_index(W, M, index, V);
    1896             :     }
    1897             :   }
    1898      177198 :   return V;
    1899             : }
    1900             : 
    1901             : /* express +oo->q=a/b in terms of the Z[G]-generators, trivial action */
    1902             : static void
    1903        7630 : Q_log_trivial(GEN v, GEN W, GEN q)
    1904             : {
    1905        7630 :   GEN Q, W3 = gel(W,3), p1N = msN_get_p1N(W);
    1906        7630 :   ulong c,d, N = p1N_get_N(p1N);
    1907             :   long i, lx;
    1908             : 
    1909        7630 :   Q = Q_log_init(N, q);
    1910        7630 :   lx = lg(Q);
    1911        7630 :   c = 0;
    1912       31906 :   for (i = 1; i < lx; i++, c = d)
    1913             :   {
    1914             :     long index;
    1915       24276 :     d = Q[i];
    1916       24276 :     if (c && !odd(i)) c = N - c;
    1917       24276 :     index = W3[ p1_index(c,d,p1N) ];
    1918       24276 :     treat_index_trivial(v, W, index);
    1919             :   }
    1920        7630 : }
    1921             : static void
    1922      587223 : M2_log_trivial(GEN V, GEN W, GEN M)
    1923             : {
    1924      587223 :   GEN p1N = gel(W,1), W3 = gel(W,3);
    1925      587223 :   ulong N = p1N_get_N(p1N);
    1926      587223 :   GEN a = gcoeff(M,1,1), b = gcoeff(M,1,2);
    1927      587223 :   GEN c = gcoeff(M,2,1), d = gcoeff(M,2,2);
    1928             :   GEN  u, v, D;
    1929             :   long index, s;
    1930             : 
    1931      587223 :   D = subii(mulii(a,d), mulii(b,c));
    1932      587223 :   s = signe(D);
    1933      592879 :   if (!s) return;
    1934      587216 :   if (is_pm1(D))
    1935             :   { /* shortcut, not need to apply Manin's trick */
    1936      214683 :     if (s < 0) d = negi(d);
    1937      214683 :     index = W3[ p1_index(umodiu(c,N),umodiu(d,N),p1N) ];
    1938      214683 :     treat_index_trivial(V, W, index);
    1939             :   }
    1940             :   else
    1941             :   {
    1942             :     GEN U, B, P, Q, PQ;
    1943             :     long i, l;
    1944      372533 :     if (!signe(c)) { Q_log_trivial(V,W,gdiv(b,d)); return; }
    1945      366884 :     (void)bezout(a,c,&u,&v);
    1946      366884 :     B = addii(mulii(b,u), mulii(d,v));
    1947             :     /* [u,v;-c,a] [a,b; c,d] = [1,B; 0,D], i.e. M = U [1,B;0,D] */
    1948      366884 :     U = mkvec2(c, u);
    1949             : 
    1950             :     /* {1/0 -> B/D} as \sum g_i, g_i unimodular paths */
    1951      366884 :     PQ = ZV_allpnqn( gboundcf(gdiv(B,D), 0) );
    1952      366884 :     P = gel(PQ,1); l = lg(P);
    1953      366884 :     Q = gel(PQ,2);
    1954     1616237 :     for (i = 1; i < l; i++, c = d)
    1955             :     {
    1956     1249353 :       d = addii(mulii(gel(U,1),gel(P,i)), mulii(gel(U,2),gel(Q,i)));
    1957     1249353 :       if (!odd(i)) c = negi(c);
    1958     1249353 :       index = W3[ p1_index(umodiu(c,N),umodiu(d,N),p1N) ];
    1959     1249353 :       treat_index_trivial(V, W, index);
    1960             :     }
    1961             :   }
    1962             : }
    1963             : 
    1964             : static GEN
    1965       16772 : cusp_to_ZC(GEN c)
    1966             : {
    1967       16772 :   switch(typ(c))
    1968             :   {
    1969             :     case t_INFINITY:
    1970          35 :       return mkcol2(gen_1,gen_0);
    1971             :     case t_INT:
    1972          91 :       return mkcol2(c,gen_1);
    1973             :     case t_FRAC:
    1974         140 :       return mkcol2(gel(c,1),gel(c,2));
    1975             :     case t_VECSMALL:
    1976       16506 :       return mkcol2(stoi(c[1]), stoi(c[2]));
    1977             :     default:
    1978           0 :       pari_err_TYPE("mspathlog",c);
    1979           0 :       return NULL;
    1980             :   }
    1981             : }
    1982             : static GEN
    1983        8386 : path2_to_M2(GEN p)
    1984        8386 : { return mkmat2(cusp_to_ZC(gel(p,1)), cusp_to_ZC(gel(p,2))); }
    1985             : static GEN
    1986       15876 : path_to_M2(GEN p)
    1987             : {
    1988       15876 :   if (lg(p) != 3) pari_err_TYPE("mspathlog",p);
    1989       15869 :   switch(typ(p))
    1990             :   {
    1991             :     case t_MAT:
    1992        9597 :       RgM_check_ZM(p,"mspathlog");
    1993        9597 :       break;
    1994             :     case t_VEC:
    1995        6272 :       p = path2_to_M2(p);
    1996        6272 :       break;
    1997           0 :     default: pari_err_TYPE("mspathlog",p);
    1998             :   }
    1999       15869 :   return p;
    2000             : }
    2001             : /* Expresses path p as \sum x_i g_i, where the g_i are our distinguished
    2002             :  * generators and x_i \in Z[\Gamma]. Returns [x_1,...,x_n] */
    2003             : GEN
    2004       12523 : mspathlog(GEN W, GEN p)
    2005             : {
    2006       12523 :   pari_sp av = avma;
    2007       12523 :   checkms(W);
    2008       12523 :   return gerepilecopy(av, M2_log(W, path_to_M2(p)));
    2009             : }
    2010             : 
    2011             : /** HECKE OPERATORS **/
    2012             : /* [a,b;c,d] * cusp */
    2013             : static GEN
    2014     1489684 : cusp_mul(long a, long b, long c, long d, GEN cusp)
    2015             : {
    2016     1489684 :   long x = cusp[1], y = cusp[2];
    2017     1489684 :   long A = a*x+b*y, B = c*x+d*y, u = cgcd(A,B);
    2018     1489684 :   if (u != 1) { A /= u; B /= u; }
    2019     1489684 :   return mkcol2s(A, B);
    2020             : }
    2021             : /* f in Gl2(Q), act on path (zm), return path_to_M2(f.path) */
    2022             : static GEN
    2023      744842 : Gl2Q_act_path(GEN f, GEN path)
    2024             : {
    2025      744842 :   long a = coeff(f,1,1), b = coeff(f,1,2);
    2026      744842 :   long c = coeff(f,2,1), d = coeff(f,2,2);
    2027      744842 :   GEN c1 = cusp_mul(a,b,c,d, gel(path,1));
    2028      744842 :   GEN c2 = cusp_mul(a,b,c,d, gel(path,2));
    2029      744842 :   return mkmat2(c1,c2);
    2030             : }
    2031             : 
    2032             : static GEN
    2033      151809 : init_act_trivial(GEN W) { return const_vecsmall(ms_get_nbE1(W), 0); }
    2034             : static GEN
    2035        3339 : mspathlog_trivial(GEN W, GEN p)
    2036             : {
    2037             :   GEN v;
    2038        3339 :   W = get_msN(W);
    2039        3339 :   v = init_act_trivial(W);
    2040        3339 :   M2_log_trivial(v, W, path_to_M2(p));
    2041        3332 :   return v;
    2042             : }
    2043             : 
    2044             : /* map from W1=Hom(Delta_0(N1),Q) -> W2=Hom(Delta_0(N2),Q), weight 2,
    2045             :  * trivial action. v a Gl2_Q or a t_VEC of Gl2_Q (\sum v[i] in Z[Gl2(Q)]).
    2046             :  * Return the matrix attached to the action of v. */
    2047             : static GEN
    2048        2702 : getMorphism_trivial(GEN WW1, GEN WW2, GEN v)
    2049             : {
    2050        2702 :   GEN T, section, gen, W1 = get_msN(WW1), W2 = get_msN(WW2);
    2051             :   long j, lv, d2;
    2052        2702 :   if (ms_get_N(W1) == 1) return cgetg(1,t_MAT);
    2053        2702 :   if (ms_get_N(W2) == 1) return zeromat(0, ms_get_nbE1(W1));
    2054        2702 :   section = msN_get_section(W2);
    2055        2702 :   gen = msN_get_genindex(W2);
    2056        2702 :   d2 = ms_get_nbE1(W2);
    2057        2702 :   T = cgetg(d2+1, t_MAT);
    2058        2702 :   lv = lg(v);
    2059      149191 :   for (j = 1; j <= d2; j++)
    2060             :   {
    2061      146489 :     GEN w = gel(section, gen[j]);
    2062      146489 :     GEN t = init_act_trivial(W1);
    2063      146489 :     pari_sp av = avma;
    2064             :     long l;
    2065      146489 :     for (l = 1; l < lv; l++) M2_log_trivial(t, W1, Gl2Q_act_path(gel(v,l), w));
    2066      146489 :     gel(T,j) = t; avma = av;
    2067             :   }
    2068        2702 :   return shallowtrans(zm_to_ZM(T));
    2069             : }
    2070             : 
    2071             : static GEN
    2072      166019 : RgV_sparse(GEN v, GEN *pind)
    2073             : {
    2074             :   long i, l, k;
    2075      166019 :   GEN w = cgetg_copy(v,&l), ind = cgetg(l, t_VECSMALL);
    2076    17145044 :   for (i = k = 1; i < l; i++)
    2077             :   {
    2078    16979025 :     GEN c = gel(v,i);
    2079    16979025 :     if (typ(c) == t_INT) continue;
    2080      785603 :     gel(w,k) = c; ind[k] = i; k++;
    2081             :   }
    2082      166019 :   setlg(w,k); setlg(ind,k);
    2083      166019 :   *pind = ind; return w;
    2084             : }
    2085             : 
    2086             : static int
    2087      163065 : mat2_isidentity(GEN M)
    2088             : {
    2089      163065 :   GEN A = gel(M,1), B = gel(M,2);
    2090      163065 :   return A[1] == 1 && A[2] == 0 && B[1] == 0 && B[2] == 1;
    2091             : }
    2092             : /* path a mat22/mat22s, return log(f.path)^* . f in sparse form */
    2093             : static GEN
    2094      166019 : M2_logf(GEN Wp, GEN path, GEN f)
    2095             : {
    2096      166019 :   pari_sp av = avma;
    2097             :   GEN ind, L;
    2098             :   long i, l;
    2099      166019 :   if (f)
    2100      160951 :     path = Gl2Q_act_path(f, path);
    2101        5068 :   else if (typ(gel(path,1)) == t_VECSMALL)
    2102        2114 :     path = path2_to_M2(path);
    2103      166019 :   L = M2_log(Wp, path);
    2104      166019 :   L = RgV_sparse(L,&ind); l = lg(L);
    2105      166019 :   for (i = 1; i < l; i++) gel(L,i) = ZSl2_star(gel(L,i));
    2106      166019 :   if (f) ZGC_G_mul_inplace(L, mat2_to_ZM(f));
    2107      166019 :   return gerepilecopy(av, mkvec2(ind,L));
    2108             : }
    2109             : 
    2110             : static hashtable *
    2111        3738 : Gl2act_cache(long dim) { return set_init(dim*10); }
    2112             : 
    2113             : /* f zm/ZM in Gl_2(Q), acts from the left on Delta, which is generated by
    2114             :  * (g_i) as Z[Gamma1]-module, and by (G_i) as Z[Gamma2]-module.
    2115             :  * We have f.G_j = \sum_i \lambda_{i,j} g_i,   \lambda_{i,j} in Z[Gamma1]
    2116             :  * For phi in Hom_Gamma1(D,V), g in D, phi | f is in Hom_Gamma2(D,V) and
    2117             :  *  (phi | f)(G_j) = phi(f.G_j) | f
    2118             :  *                 = phi( \sum_i \lambda_{i,j} g_i ) | f
    2119             :  *                 = \sum_i phi(g_i) | (\lambda_{i,j}^* f)
    2120             :  *                 = \sum_i phi(g_i) | \mu_{i,j}(f)
    2121             :  * More generally
    2122             :  *  (\sum_k (phi |v_k))(G_j) = \sum_i phi(g_i) | \Mu_{i,j}
    2123             :  * with \Mu_{i,j} = \sum_k \mu{i,j}(v_k)
    2124             :  * Return the \Mu_{i,j} matrix as vector of sparse columns of operators on V */
    2125             : static GEN
    2126        3262 : init_dual_act(GEN v, GEN W1, GEN W2, struct m_act *S)
    2127             : {
    2128        3262 :   GEN section = ms_get_section(W2), gen = ms_get_genindex(W2);
    2129             :   /* HACK: the actions we consider in dimension 1 are trivial and in
    2130             :    * characteristic != 2, 3 => torsion generators are 0
    2131             :    * [satisfy e.g. (1+gamma).g = 0 => \phi(g) | 1+gamma  = 0 => \phi(g) = 0 */
    2132        3262 :   long j, lv = lg(v), dim = S->dim == 1? ms_get_nbE1(W2): lg(gen)-1;
    2133        3262 :   GEN T = cgetg(dim+1, t_VEC);
    2134        3262 :   hashtable *H = Gl2act_cache(dim);
    2135       41139 :   for (j = 1; j <= dim; j++)
    2136             :   {
    2137       37877 :     pari_sp av = avma;
    2138       37877 :     GEN w = gel(section, gen[j]); /* path_to_zm( E1/T2/T3 element ) */
    2139       37877 :     GEN t = NULL;
    2140             :     long k;
    2141      200942 :     for (k = 1; k < lv; k++)
    2142             :     {
    2143      163065 :       GEN tk, f = gel(v,k);
    2144      163065 :       if (typ(gel(f,1)) != t_VECSMALL) f = ZM_to_zm(f);
    2145      163065 :       if (mat2_isidentity(f)) f = NULL;
    2146      163065 :       tk = M2_logf(W1, w, f); /* mu_{.,j}(v[k]) as sparse vector */
    2147      163065 :       t = t? ZGCs_add(t, tk): tk;
    2148             :     }
    2149       37877 :     gel(T,j) = gerepilecopy(av, t);
    2150             :   }
    2151       41139 :   for (j = 1; j <= dim; j++)
    2152             :   {
    2153       37877 :     ZGl2QC_preload(S, gel(T,j), H);
    2154       37877 :     ZGl2QC_to_act(S, gel(T,j), H);
    2155             :   }
    2156        3262 :   return T;
    2157             : }
    2158             : 
    2159             : /* modular symbol given by phi[j] = \phi(G_j)
    2160             :  * \sum L[i]*phi[i], L a sparse column of operators */
    2161             : static GEN
    2162      354354 : dense_act_col(GEN col, GEN phi)
    2163             : {
    2164      354354 :   GEN s = NULL, colind = gel(col,1), colval = gel(col,2);
    2165      354354 :   long i, l = lg(colind), lphi = lg(phi);
    2166     5630121 :   for (i = 1; i < l; i++)
    2167             :   {
    2168     5278490 :     long a = colind[i];
    2169             :     GEN t;
    2170     5278490 :     if (a >= lphi) break; /* happens if k=2: torsion generator t omitted */
    2171     5275767 :     t = gel(phi, a); /* phi(G_a) */
    2172     5275767 :     t = RgM_RgC_mul(gel(colval,i), t);
    2173     5275767 :     s = s? RgC_add(s, t): t;
    2174             :   }
    2175      354354 :   return s;
    2176             : }
    2177             : /* modular symbol given by \phi( G[ind[j]] ) = val[j]
    2178             :  * \sum L[i]*phi[i], L a sparse column of operators */
    2179             : static GEN
    2180      780283 : sparse_act_col(GEN col, GEN phi)
    2181             : {
    2182      780283 :   GEN s = NULL, colind = gel(col,1), colval = gel(col,2);
    2183      780283 :   GEN ind = gel(phi,2), val = gel(phi,3);
    2184      780283 :   long a, l = lg(ind);
    2185      780283 :   if (lg(gel(phi,1)) == 1) return RgM_RgC_mul(gel(colval,1), gel(val,1));
    2186     3036978 :   for (a = 1; a < l; a++)
    2187             :   {
    2188     2257094 :     GEN t = gel(val, a); /* phi(G_i) */
    2189     2257094 :     long i = zv_search(colind, ind[a]);
    2190     2257094 :     if (!i) continue;
    2191      543228 :     t = RgM_RgC_mul(gel(colval,i), t);
    2192      543228 :     s = s? RgC_add(s, t): t;
    2193             :   }
    2194      779884 :   return s;
    2195             : }
    2196             : static int
    2197       69545 : phi_sparse(GEN phi) { return typ(gel(phi,1)) == t_VECSMALL; }
    2198             : /* phi in Hom_Gamma1(Delta, V), return the matrix whose colums are the
    2199             :  *   \sum_i phi(g_i) | \mu_{i,j} = (phi|f)(G_j),
    2200             :  * see init_dual_act. */
    2201             : static GEN
    2202       69545 : dual_act(long dimV, GEN act, GEN phi)
    2203             : {
    2204       69545 :   long l = lg(act), j;
    2205       69545 :   GEN v = cgetg(l, t_MAT);
    2206       69545 :   GEN (*ACT)(GEN,GEN) = phi_sparse(phi)? sparse_act_col: dense_act_col;
    2207     1200850 :   for (j = 1; j < l; j++)
    2208             :   {
    2209     1131305 :     pari_sp av = avma;
    2210     1131305 :     GEN s = ACT(gel(act,j), phi);
    2211     1131305 :     gel(v,j) = s? gerepileupto(av,s): zerocol(dimV);
    2212             :   }
    2213       69545 :   return v;
    2214             : }
    2215             : 
    2216             : /* in level N > 1 */
    2217             : static void
    2218       59507 : msk_get_st(GEN W, long *s, long *t)
    2219       59507 : { GEN st = gmael(W,3,3); *s = st[1]; *t = st[2]; }
    2220             : static GEN
    2221       59507 : msk_get_link(GEN W) { return gmael(W,3,4); }
    2222             : static GEN
    2223       59948 : msk_get_inv(GEN W) { return gmael(W,3,5); }
    2224             : /* \phi in Hom(Delta, V), \phi(G_k) = phi[k]. Write \phi as
    2225             :  *   \sum_{i,j} mu_{i,j} phi_{i,j}, mu_{i,j} in Q */
    2226             : static GEN
    2227       59948 : getMorphism_basis(GEN W, GEN phi)
    2228             : {
    2229       59948 :   GEN R, Q, Ls, T0, T1, Ts, link, basis, inv = msk_get_inv(W);
    2230             :   long i, j, r, s, t, dim, lvecT;
    2231             : 
    2232       59948 :   if (ms_get_N(W) == 1) return ZC_apply_dinv(inv, gel(phi,1));
    2233       59507 :   lvecT = lg(phi);
    2234       59507 :   basis = msk_get_basis(W);
    2235       59507 :   dim = lg(basis)-1;
    2236       59507 :   R = zerocol(dim);
    2237       59507 :   msk_get_st(W, &s, &t);
    2238       59507 :   link = msk_get_link(W);
    2239      791350 :   for (r = 2; r < lvecT; r++)
    2240             :   {
    2241             :     GEN Tr, L;
    2242      731843 :     if (r == s) continue;
    2243      672336 :     Tr = gel(phi,r); /* Phi(G_r), r != 1,s */
    2244      672336 :     L = gel(link, r);
    2245      672336 :     Q = ZC_apply_dinv(gel(inv,r), Tr);
    2246             :     /* write Phi(G_r) as sum_{a,b} mu_{a,b} Phi_{a,b}(G_r) */
    2247      672336 :     for (j = 1; j < lg(L); j++) gel(R, L[j]) = gel(Q,j);
    2248             :   }
    2249       59507 :   Ls = gel(link, s);
    2250       59507 :   T1 = gel(phi,1); /* Phi(G_1) */
    2251       59507 :   gel(R, Ls[t]) = gel(T1, 1);
    2252             : 
    2253       59507 :   T0 = NULL;
    2254      791350 :   for (i = 2; i < lg(link); i++)
    2255             :   {
    2256             :     GEN L;
    2257      731843 :     if (i == s) continue;
    2258      672336 :     L = gel(link,i);
    2259     3514154 :     for (j =1 ; j < lg(L); j++)
    2260             :     {
    2261     2841818 :       long n = L[j]; /* phi_{i,j} = basis[n] */
    2262     2841818 :       GEN mu_ij = gel(R, n);
    2263     2841818 :       GEN phi_ij = gel(basis, n), pols = gel(phi_ij,3);
    2264     2841818 :       GEN z = RgC_Rg_mul(gel(pols, 3), mu_ij);
    2265     2841818 :       T0 = T0? RgC_add(T0, z): z; /* += mu_{i,j} Phi_{i,j} (G_s) */
    2266             :     }
    2267             :   }
    2268       59507 :   Ts = gel(phi,s); /* Phi(G_s) */
    2269       59507 :   if (T0) Ts = RgC_sub(Ts, T0);
    2270             :   /* solve \sum_{j!=t} mu_{s,j} Phi_{s,j}(G_s) = Ts */
    2271       59507 :   Q = ZC_apply_dinv(gel(inv,s), Ts);
    2272       59507 :   for (j = 1; j < t; j++) gel(R, Ls[j]) = gel(Q,j);
    2273             :   /* avoid mu_{s,t} */
    2274       59507 :   for (j = t; j < lg(Q); j++) gel(R, Ls[j+1]) = gel(Q,j);
    2275       59507 :   return R;
    2276             : }
    2277             : 
    2278             : /* a = s(g_i) for some modular symbol s; b in Z[G]
    2279             :  * return s(b.g_i) = b^* . s(g_i) */
    2280             : static GEN
    2281      115458 : ZGl2Q_act_s(GEN b, GEN a, long k)
    2282             : {
    2283      115458 :   if (typ(b) == t_INT)
    2284             :   {
    2285       58604 :     if (!signe(b)) return gen_0;
    2286          14 :     switch(typ(a))
    2287             :     {
    2288             :       case t_POL:
    2289          14 :         a = RgX_to_RgC(a, k-1); /*fall through*/
    2290             :       case t_COL:
    2291          14 :         a = RgC_Rg_mul(a,b);
    2292          14 :         break;
    2293           0 :       default: a = scalarcol_shallow(b,k-1);
    2294             :     }
    2295             :   }
    2296             :   else
    2297             :   {
    2298       56854 :     b = RgX_act_ZGl2Q(ZSl2_star(b), k);
    2299       56854 :     switch(typ(a))
    2300             :     {
    2301             :       case t_POL:
    2302          63 :         a = RgX_to_RgC(a, k-1); /*fall through*/
    2303             :       case t_COL:
    2304       45094 :         a = RgM_RgC_mul(b,a);
    2305       45094 :         break;
    2306       11760 :       default: a = RgC_Rg_mul(gel(b,1),a);
    2307             :     }
    2308             :   }
    2309       56868 :   return a;
    2310             : }
    2311             : 
    2312             : static int
    2313          21 : checksymbol(GEN W, GEN s)
    2314             : {
    2315             :   GEN t, annT2, annT31, singlerel;
    2316             :   long i, k, l, nbE1, nbT2, nbT31;
    2317          21 :   k = msk_get_weight(W);
    2318          21 :   W = get_msN(W);
    2319          21 :   nbE1 = ms_get_nbE1(W);
    2320          21 :   singlerel = gel(W,10);
    2321          21 :   l = lg(singlerel);
    2322          21 :   if (k == 2)
    2323             :   {
    2324           0 :     for (i = nbE1+1; i < l; i++)
    2325           0 :       if (!gequal0(gel(s,i))) return 0;
    2326           0 :     return 1;
    2327             :   }
    2328          21 :   annT2 = msN_get_annT2(W); nbT2 = lg(annT2)-1;
    2329          21 :   annT31 = msN_get_annT31(W);nbT31 = lg(annT31)-1;
    2330          21 :   t = NULL;
    2331          84 :   for (i = 1; i < l; i++)
    2332             :   {
    2333          63 :     GEN a = gel(s,i);
    2334          63 :     a = ZGl2Q_act_s(gel(singlerel,i), a, k);
    2335          63 :     t = t? gadd(t, a): a;
    2336             :   }
    2337          21 :   if (!gequal0(t)) return 0;
    2338          14 :   for (i = 1; i <= nbT2; i++)
    2339             :   {
    2340           0 :     GEN a = gel(s,i + nbE1);
    2341           0 :     a = ZGl2Q_act_s(gel(annT2,i), a, k);
    2342           0 :     if (!gequal0(a)) return 0;
    2343             :   }
    2344          28 :   for (i = 1; i <= nbT31; i++)
    2345             :   {
    2346          14 :     GEN a = gel(s,i + nbE1 + nbT2);
    2347          14 :     a = ZGl2Q_act_s(gel(annT31,i), a, k);
    2348          14 :     if (!gequal0(a)) return 0;
    2349             :   }
    2350          14 :   return 1;
    2351             : }
    2352             : GEN
    2353          56 : msissymbol(GEN W, GEN s)
    2354             : {
    2355             :   long k, nbgen;
    2356          56 :   checkms(W);
    2357          56 :   k = msk_get_weight(W);
    2358          56 :   nbgen = ms_get_nbgen(W);
    2359          56 :   switch(typ(s))
    2360             :   {
    2361             :     case t_VEC: /* values s(g_i) */
    2362          21 :       if (lg(s)-1 != nbgen) return gen_0;
    2363          21 :       break;
    2364             :     case t_COL:
    2365          28 :       if (msk_get_sign(W))
    2366             :       {
    2367           0 :         GEN star = gel(msk_get_starproj(W), 1);
    2368           0 :         if (lg(star) == lg(s)) return gen_1;
    2369             :       }
    2370          28 :       if (k == 2) /* on the dual basis of (g_i) */
    2371             :       {
    2372           0 :         if (lg(s)-1 != nbgen) return gen_0;
    2373             :       }
    2374             :       else
    2375             :       {
    2376          28 :         GEN basis = msk_get_basis(W);
    2377          28 :         return (lg(s) == lg(basis))? gen_1: gen_0;
    2378             :       }
    2379           0 :       break;
    2380             :     case t_MAT:
    2381             :     {
    2382           7 :       long i, l = lg(s);
    2383           7 :       GEN v = cgetg(l, t_VEC);
    2384           7 :       for (i = 1; i < l; i++) gel(v,i) = msissymbol(W,gel(s,i))? gen_1: gen_0;
    2385           7 :       return v;
    2386             :     }
    2387           0 :     default: return gen_0;
    2388             :   }
    2389          21 :   return checksymbol(W,s)? gen_1: gen_0;
    2390             : }
    2391             : #if DEBUG
    2392             : /* phi is a sparse symbol from msk_get_basis, return phi(G_j) */
    2393             : static GEN
    2394             : phi_Gj(GEN W, GEN phi, long j)
    2395             : {
    2396             :   GEN ind = gel(phi,2), pols = gel(phi,3);
    2397             :   long i = vecsmall_isin(ind,j);
    2398             :   return i? gel(pols,i): NULL;
    2399             : }
    2400             : /* check that \sum d_i phi_i(G_j)  = T_j for all j */
    2401             : static void
    2402             : checkdec(GEN W, GEN D, GEN T)
    2403             : {
    2404             :   GEN B = msk_get_basis(W);
    2405             :   long i, j;
    2406             :   if (!checksymbol(W,T)) pari_err_BUG("checkdec");
    2407             :   for (j = 1; j < lg(T); j++)
    2408             :   {
    2409             :     GEN S = gen_0;
    2410             :     for (i = 1; i < lg(D); i++)
    2411             :     {
    2412             :       GEN d = gel(D,i), v = phi_Gj(W, gel(B,i), j);
    2413             :       if (!v || gequal0(d)) continue;
    2414             :       S = gadd(S, gmul(d, v));
    2415             :     }
    2416             :     /* S = \sum_i d_i phi_i(G_j) */
    2417             :     if (!gequal(S, gel(T,j)))
    2418             :       pari_warn(warner, "checkdec j = %ld\n\tS = %Ps\n\tT = %Ps", j,S,gel(T,j));
    2419             :   }
    2420             : }
    2421             : #endif
    2422             : 
    2423             : /* map op: W1 = Hom(Delta_0(N1),V) -> W2 = Hom(Delta_0(N2),V), given by
    2424             :  * \sum v[i], v[i] in Gl2(Q) */
    2425             : static GEN
    2426        5481 : getMorphism(GEN W1, GEN W2, GEN v)
    2427             : {
    2428             :   struct m_act S;
    2429             :   GEN B1, M, act;
    2430        5481 :   long a, l, k = msk_get_weight(W1);
    2431        5481 :   if (k == 2) return getMorphism_trivial(W1,W2,v);
    2432        2779 :   S.k = k;
    2433        2779 :   S.dim = k-1;
    2434        2779 :   S.act = &_RgX_act_Gl2Q;
    2435        2779 :   act = init_dual_act(v,W1,W2,&S);
    2436        2779 :   B1 = msk_get_basis(W1);
    2437        2779 :   l = lg(B1); M = cgetg(l, t_MAT);
    2438       61656 :   for (a = 1; a < l; a++)
    2439             :   {
    2440       58877 :     pari_sp av = avma;
    2441       58877 :     GEN phi = dual_act(S.dim, act, gel(B1,a));
    2442       58877 :     GEN D = getMorphism_basis(W2, phi);
    2443             : #if DEBUG
    2444             :     checkdec(W2,D,T);
    2445             : #endif
    2446       58877 :     gel(M,a) = gerepilecopy(av, D);
    2447             :   }
    2448        2779 :   return M;
    2449             : }
    2450             : static GEN
    2451        4347 : msendo(GEN W, GEN v) { return getMorphism(W, W, v); }
    2452             : 
    2453             : static GEN
    2454        2583 : endo_project(GEN W, GEN e, GEN H)
    2455             : {
    2456        2583 :   if (msk_get_sign(W)) e = Qevproj_apply(e, msk_get_starproj(W));
    2457        2583 :   if (H) e = Qevproj_apply(e, Qevproj_init0(H));
    2458        2583 :   return e;
    2459             : }
    2460             : static GEN
    2461        2926 : mshecke_i(GEN W, ulong p)
    2462             : {
    2463        2926 :   GEN v = ms_get_N(W) % p? Tp_matrices(p): Up_matrices(p);
    2464        2926 :   return msendo(W,v);
    2465             : }
    2466             : GEN
    2467        2534 : mshecke(GEN W, long p, GEN H)
    2468             : {
    2469        2534 :   pari_sp av = avma;
    2470             :   GEN T;
    2471        2534 :   checkms(W);
    2472        2534 :   if (p <= 1) pari_err_PRIME("mshecke",stoi(p));
    2473        2534 :   T = mshecke_i(W,p);
    2474        2534 :   T = endo_project(W,T,H);
    2475        2534 :   return gerepilecopy(av, T);
    2476             : }
    2477             : 
    2478             : static GEN
    2479          42 : msatkinlehner_i(GEN W, long Q)
    2480             : {
    2481          42 :   long N = ms_get_N(W);
    2482             :   GEN v;
    2483          42 :   if (Q == 1) return matid(msk_get_dim(W));
    2484          28 :   if (Q == N) return msendo(W, mkvec(mat2(0,1,-N,0)));
    2485          21 :   if (N % Q) pari_err_DOMAIN("msatkinlehner","N % Q","!=",gen_0,stoi(Q));
    2486          14 :   v = WQ_matrix(N, Q);
    2487          14 :   if (!v) pari_err_DOMAIN("msatkinlehner","gcd(Q,N/Q)","!=",gen_1,stoi(Q));
    2488          14 :   return msendo(W,mkvec(v));
    2489             : }
    2490             : GEN
    2491          42 : msatkinlehner(GEN W, long Q, GEN H)
    2492             : {
    2493          42 :   pari_sp av = avma;
    2494             :   GEN w;
    2495             :   long k;
    2496          42 :   checkms(W);
    2497          42 :   k = msk_get_weight(W);
    2498          42 :   if (Q <= 0) pari_err_DOMAIN("msatkinlehner","Q","<=",gen_0,stoi(Q));
    2499          42 :   w = msatkinlehner_i(W,Q);
    2500          35 :   w = endo_project(W,w,H);
    2501          35 :   if (k > 2 && Q != 1) w = RgM_Rg_div(w, powuu(Q,(k-2)>>1));
    2502          35 :   return gerepilecopy(av, w);
    2503             : }
    2504             : 
    2505             : static GEN
    2506        1400 : msstar_i(GEN W) { return msendo(W, mkvec(mat2(-1,0,0,1))); }
    2507             : GEN
    2508          14 : msstar(GEN W, GEN H)
    2509             : {
    2510          14 :   pari_sp av = avma;
    2511             :   GEN s;
    2512          14 :   checkms(W);
    2513          14 :   s = msstar_i(W);
    2514          14 :   s = endo_project(W,s,H);
    2515          14 :   return gerepilecopy(av, s);
    2516             : }
    2517             : 
    2518             : #if 0
    2519             : /* is \Gamma_0(N) cusp1 = \Gamma_0(N) cusp2 ? */
    2520             : static int
    2521             : iscuspeq(ulong N, GEN cusp1, GEN cusp2)
    2522             : {
    2523             :   long p1, q1, p2, q2, s1, s2, d;
    2524             :   p1 = cusp1[1]; p2 = cusp2[1];
    2525             :   q1 = cusp1[2]; q2 = cusp2[2];
    2526             :   d = Fl_mul(smodss(q1,N),smodss(q2,N), N);
    2527             :   d = ugcd(d, N);
    2528             : 
    2529             :   s1 = q1 > 2? Fl_inv(smodss(p1,q1), q1): 1;
    2530             :   s2 = q2 > 2? Fl_inv(smodss(p2,q2), q2): 1;
    2531             :   return Fl_mul(s1,q2,d) == Fl_mul(s2,q1,d);
    2532             : }
    2533             : #endif
    2534             : 
    2535             : /* return E_c(r) */
    2536             : static GEN
    2537        2926 : get_Ec_r(GEN c, long k)
    2538             : {
    2539        2926 :   long p = c[1], q = c[2], u, v;
    2540             :   GEN gr;
    2541        2926 :   (void)cbezout(p, q, &u, &v);
    2542        2926 :   gr = mat2(p, -v, q, u); /* g . (1:0) = (p:q) */
    2543        2926 :   return voo_act_Gl2Q(sl2_inv(gr), k);
    2544             : }
    2545             : /* N > 1; returns the modular symbol attached to the cusp c := p/q via the rule
    2546             :  * E_c(path from a to b in Delta_0) := E_c(b) - E_c(a), where
    2547             :  * E_c(r) := 0 if r != c mod Gamma
    2548             :  *           v_oo | gamma_r^(-1)
    2549             :  * where v_oo is stable by T = [1,1;0,1] (i.e x^(k-2)) and
    2550             :  * gamma_r . (1:0) = r, for some gamma_r in SL_2(Z) * */
    2551             : static GEN
    2552         441 : msfromcusp_trivial(GEN W, GEN c)
    2553             : {
    2554         441 :   GEN section = ms_get_section(W), gen = ms_get_genindex(W);
    2555         441 :   GEN S = ms_get_hashcusps(W);
    2556         441 :   long j, ic = cusp_index(c, S), l = ms_get_nbE1(W)+1;
    2557         441 :   GEN phi = cgetg(l, t_COL);
    2558       90062 :   for (j = 1; j < l; j++)
    2559             :   {
    2560       89621 :     GEN vj, g = gel(section, gen[j]); /* path_to_zm(generator) */
    2561       89621 :     GEN c1 = gel(g,1), c2 = gel(g,2);
    2562       89621 :     long i1 = cusp_index(c1, S);
    2563       89621 :     long i2 = cusp_index(c2, S);
    2564       89621 :     if (i1 == ic)
    2565        3206 :       vj = (i2 == ic)?  gen_0: gen_1;
    2566             :     else
    2567       86415 :       vj = (i2 == ic)? gen_m1: gen_0;
    2568       89621 :     gel(phi, j) = vj;
    2569             :   }
    2570         441 :   return phi;
    2571             : }
    2572             : static GEN
    2573        1512 : msfromcusp_i(GEN W, GEN c)
    2574             : {
    2575             :   GEN section, gen, S, phi;
    2576        1512 :   long j, ic, l, k = msk_get_weight(W);
    2577        1512 :   if (k == 2)
    2578             :   {
    2579         441 :     long N = ms_get_N(W);
    2580         441 :     return N == 1? cgetg(1,t_COL): msfromcusp_trivial(W, c);
    2581             :   }
    2582        1071 :   k = msk_get_weight(W);
    2583        1071 :   section = ms_get_section(W);
    2584        1071 :   gen = ms_get_genindex(W);
    2585        1071 :   S = ms_get_hashcusps(W);
    2586        1071 :   ic = cusp_index(c, S);
    2587        1071 :   l = lg(gen);
    2588        1071 :   phi = cgetg(l, t_COL);
    2589       12579 :   for (j = 1; j < l; j++)
    2590             :   {
    2591       11508 :     GEN vj = NULL, g = gel(section, gen[j]); /* path_to_zm(generator) */
    2592       11508 :     GEN c1 = gel(g,1), c2 = gel(g,2);
    2593       11508 :     long i1 = cusp_index(c1, S);
    2594       11508 :     long i2 = cusp_index(c2, S);
    2595       11508 :     if (i1 == ic) vj = get_Ec_r(c1, k);
    2596       11508 :     if (i2 == ic)
    2597             :     {
    2598        1463 :       GEN s = get_Ec_r(c2, k);
    2599        1463 :       vj = vj? gsub(vj, s): gneg(s);
    2600             :     }
    2601       11508 :     if (!vj) vj = zerocol(k-1);
    2602       11508 :     gel(phi, j) = vj;
    2603             :   }
    2604        1071 :   return getMorphism_basis(W, phi);
    2605             : }
    2606             : GEN
    2607          28 : msfromcusp(GEN W, GEN c)
    2608             : {
    2609          28 :   pari_sp av = avma;
    2610             :   long N;
    2611          28 :   checkms(W);
    2612          28 :   N = ms_get_N(W);
    2613          28 :   switch(typ(c))
    2614             :   {
    2615             :     case t_INFINITY:
    2616           7 :       c = mkvecsmall2(1,0);
    2617           7 :       break;
    2618             :     case t_INT:
    2619          14 :       c = mkvecsmall2(smodis(c,N), 1);
    2620          14 :       break;
    2621             :     case t_FRAC:
    2622           7 :       c = mkvecsmall2(smodis(gel(c,1),N), smodis(gel(c,2),N));
    2623           7 :       break;
    2624             :     default:
    2625           0 :       pari_err_TYPE("msfromcusp",c);
    2626             :   }
    2627          28 :   return gerepilecopy(av, msfromcusp_i(W,c));
    2628             : }
    2629             : 
    2630             : static GEN
    2631         357 : mseisenstein_i(GEN W)
    2632             : {
    2633         357 :   GEN M, S = ms_get_hashcusps(W), cusps = gel(S,3);
    2634         357 :   long i, l = lg(cusps);
    2635         357 :   if (msk_get_weight(W)==2) l--;
    2636         357 :   M = cgetg(l, t_MAT);
    2637         357 :   for (i = 1; i < l; i++) gel(M,i) = msfromcusp_i(W, gel(cusps,i));
    2638         357 :   return Qevproj_star(W, QM_image(M));
    2639             : }
    2640             : GEN
    2641         357 : mseisenstein(GEN W)
    2642             : {
    2643         357 :   pari_sp av = avma;
    2644         357 :   checkms(W);
    2645         357 :   return gerepilecopy(av, Qevproj_init(mseisenstein_i(W)));
    2646             : }
    2647             : 
    2648             : /* upper bound for log_2 |charpoly(T_p|S)|, where S is a cuspidal subspace of
    2649             :  * dimension d, k is the weight */
    2650             : #if 0
    2651             : static long
    2652             : TpS_char_bound(ulong p, long k, long d)
    2653             : { /* |eigenvalue| <= 2 p^(k-1)/2 */
    2654             :   return d * (2 + (log2((double)p)*(k-1))/2);
    2655             : }
    2656             : #endif
    2657             : static long
    2658         336 : TpE_char_bound(ulong p, long k, long d)
    2659             : { /* |eigenvalue| <= 2 p^(k-1) */
    2660         336 :   return d * (2 + log2((double)p)*(k-1));
    2661             : }
    2662             : 
    2663             : GEN
    2664         336 : mscuspidal(GEN W, long flag)
    2665             : {
    2666         336 :   pari_sp av = avma;
    2667             :   GEN S, E, M, T, TE, chE;
    2668             :   long bit;
    2669             :   forprime_t F;
    2670             :   ulong p, N;
    2671             :   pari_timer ti;
    2672             : 
    2673         336 :   E = mseisenstein(W);
    2674         336 :   N = ms_get_N(W);
    2675         336 :   (void)u_forprime_init(&F, 2, ULONG_MAX);
    2676         336 :   while ((p = u_forprime_next(&F)))
    2677         469 :     if (N % p) break;
    2678         336 :   if (DEBUGLEVEL) timer_start(&ti);
    2679         336 :   T = mshecke(W, p, NULL);
    2680         336 :   if (DEBUGLEVEL) timer_printf(&ti,"Tp, p = %ld", p);
    2681         336 :   TE = Qevproj_apply(T, E); /* T_p | E */
    2682         336 :   if (DEBUGLEVEL) timer_printf(&ti,"Qevproj_init(E)");
    2683         336 :   bit = TpE_char_bound(p, msk_get_weight(W), lg(TE)-1);
    2684         336 :   chE = QM_charpoly_ZX_bound(TE, bit);
    2685         336 :   chE = ZX_radical(chE);
    2686         336 :   M = RgX_RgM_eval(chE, T);
    2687         336 :   S = Qevproj_init(QM_image(M));
    2688         336 :   return gerepilecopy(av, flag? mkvec2(S,E): S);
    2689             : }
    2690             : 
    2691             : /** INIT ELLSYM STRUCTURE **/
    2692             : /* V a vector of ZM. If all of them have 0 last row, return NULL.
    2693             :  * Otherwise return [m,i,j], where m = V[i][last,j] contains the value
    2694             :  * of smallest absolute value */
    2695             : static GEN
    2696         931 : RgMV_find_non_zero_last_row(long offset, GEN V)
    2697             : {
    2698         931 :   long i, lasti = 0, lastj = 0, lV = lg(V);
    2699         931 :   GEN m = NULL;
    2700        4074 :   for (i = 1; i < lV; i++)
    2701             :   {
    2702        3143 :     GEN M = gel(V,i);
    2703        3143 :     long j, n, l = lg(M);
    2704        3143 :     if (l == 1) continue;
    2705        2835 :     n = nbrows(M);
    2706       13804 :     for (j = 1; j < l; j++)
    2707             :     {
    2708       10969 :       GEN a = gcoeff(M, n, j);
    2709       10969 :       if (!gequal0(a) && (!m || abscmpii(a, m) < 0))
    2710             :       {
    2711        1582 :         m = a; lasti = i; lastj = j;
    2712        1582 :         if (is_pm1(m)) goto END;
    2713             :       }
    2714             :     }
    2715             :   }
    2716             : END:
    2717         931 :   if (!m) return NULL;
    2718         623 :   return mkvec2(m, mkvecsmall2(lasti+offset, lastj));
    2719             : }
    2720             : /* invert the d_oo := (\gamma_oo - 1) operator, acting on
    2721             :  * [x^(k-2), ..., y^(k-2)] */
    2722             : static GEN
    2723         623 : Delta_inv(GEN doo, long k)
    2724             : {
    2725         623 :   GEN M = RgX_act_ZGl2Q(doo, k);
    2726         623 :   M = RgM_minor(M, k-1, 1); /* 1st column and last row are 0 */
    2727         623 :   return ZM_inv_denom(M);
    2728             : }
    2729             : /* The ZX P = \sum a_i x^i y^{k-2-i} is given by the ZV [a_0, ..., a_k-2]~,
    2730             :  * return Q and d such that P = doo Q + d y^k-2, where d in Z and Q */
    2731             : static GEN
    2732       12831 : doo_decompose(GEN dinv, GEN P, GEN *pd)
    2733             : {
    2734       12831 :   long l = lg(P); *pd = gel(P, l-1);
    2735       12831 :   P = vecslice(P, 1, l-2);
    2736       12831 :   return shallowconcat(gen_0, ZC_apply_dinv(dinv, P));
    2737             : }
    2738             : 
    2739             : static GEN
    2740       12831 : get_phi_ij(long i,long j,long n, long s,long t,GEN P_st,GEN Q_st,GEN d_st,
    2741             :            GEN P_ij, GEN lP_ij, GEN dinv)
    2742             : {
    2743             :   GEN ind, pols;
    2744       12831 :   if (i == s && j == t)
    2745             :   {
    2746         623 :     ind = mkvecsmall(1);
    2747         623 :     pols = mkvec(scalarcol_shallow(gen_1, lg(P_st)-1)); /* x^{k-2} */
    2748             :   }
    2749             :   else
    2750             :   {
    2751       12208 :     GEN d_ij, Q_ij = doo_decompose(dinv, lP_ij, &d_ij);
    2752       12208 :     GEN a = ZC_Z_mul(P_ij, d_st);
    2753       12208 :     GEN b = ZC_Z_mul(P_st, negi(d_ij));
    2754       12208 :     GEN c = RgC_sub(RgC_Rg_mul(Q_ij, d_st), RgC_Rg_mul(Q_st, d_ij));
    2755       12208 :     if (i == s) { /* j != t */
    2756        1645 :       ind = mkvecsmall2(1, s);
    2757        1645 :       pols = mkvec2(c, ZC_add(a, b));
    2758             :     } else {
    2759       10563 :       ind = mkvecsmall3(1, i, s);
    2760       10563 :       pols = mkvec3(c, a, b); /* image of g_1, g_i, g_s */
    2761             :     }
    2762       12208 :     pols = Q_primpart(pols);
    2763             :   }
    2764       12831 :   return mkvec3(mkvecsmall3(i,j,n), ind, pols);
    2765             : }
    2766             : 
    2767             : static GEN
    2768         749 : mskinit_trivial(GEN WN)
    2769             : {
    2770         749 :   long dim = ms_get_nbE1(WN);
    2771         749 :   return mkvec3(WN, gen_0, mkvec2(gen_0,mkvecsmall2(2, dim)));
    2772             : }
    2773             : /* sum of #cols of the matrices contained in V */
    2774             : static long
    2775        1246 : RgMV_dim(GEN V)
    2776             : {
    2777        1246 :   long l = lg(V), d = 0, i;
    2778        1246 :   for (i = 1; i < l; i++) d += lg(gel(V,i)) - 1;
    2779        1246 :   return d;
    2780             : }
    2781             : static GEN
    2782         623 : mskinit_nontrivial(GEN WN, long k)
    2783             : {
    2784         623 :   GEN annT2 = gel(WN,8), annT31 = gel(WN,9), singlerel = gel(WN,10);
    2785             :   GEN link, basis, monomials, Inv;
    2786         623 :   long nbE1 = ms_get_nbE1(WN);
    2787         623 :   GEN dinv = Delta_inv(ZG_neg( ZSl2_star(gel(singlerel,1)) ), k);
    2788         623 :   GEN p1 = cgetg(nbE1+1, t_VEC), remove;
    2789         623 :   GEN p2 = ZGV_tors(annT2, k);
    2790         623 :   GEN p3 = ZGV_tors(annT31, k);
    2791         623 :   GEN gentor = shallowconcat(p2, p3);
    2792             :   GEN P_st, lP_st, Q_st, d_st;
    2793             :   long n, i, dim, s, t, u;
    2794         623 :   gel(p1, 1) = cgetg(1,t_MAT); /* dummy */
    2795        3360 :   for (i = 2; i <= nbE1; i++) /* skip 1st element = (\gamma_oo-1)g_oo */
    2796             :   {
    2797        2737 :     GEN z = gel(singlerel, i);
    2798        2737 :     gel(p1, i) = RgX_act_ZGl2Q(ZSl2_star(z), k);
    2799             :   }
    2800         623 :   remove = RgMV_find_non_zero_last_row(nbE1, gentor);
    2801         623 :   if (!remove) remove = RgMV_find_non_zero_last_row(0, p1);
    2802         623 :   if (!remove) pari_err_BUG("msinit [no y^k-2]");
    2803         623 :   remove = gel(remove,2); /* [s,t] */
    2804         623 :   s = remove[1];
    2805         623 :   t = remove[2];
    2806             :   /* +1 because of = x^(k-2), but -1 because of Manin relation */
    2807         623 :   dim = (k-1)*(nbE1-1) + RgMV_dim(p2) + RgMV_dim(p3);
    2808             :   /* Let (g_1,...,g_d) be the Gamma-generators of Delta, g_1 = g_oo.
    2809             :    * We describe modular symbols by the collection phi(g_1), ..., phi(g_d)
    2810             :    * \in V := Q[x,y]_{k-2}, with right Gamma action.
    2811             :    * For each i = 1, .., d, let V_i \subset V be the Q-vector space of
    2812             :    * allowed values for phi(g_i): with basis (P^{i,j}) given by the monomials
    2813             :    * x^(j-1) y^{k-2-(j-1)}, j = 1 .. k-1
    2814             :    * (g_i in E_1) or the solution of the torsion equations (1 + gamma)P = 0
    2815             :    * (g_i in T2) or (1 + gamma + gamma^2)P = 0 (g_i in T31). All such P
    2816             :    * are chosen in Z[x,y] with Q_content 1.
    2817             :    *
    2818             :    * The Manin relation (singlerel) is of the form \sum_i \lambda_i g_i = 0,
    2819             :    * where \lambda_i = 1 if g_i in T2 or T31, and \lambda_i = (1 - \gamma_i)
    2820             :    * for g_i in E1.
    2821             :    *
    2822             :    * If phi \in Hom_Gamma(Delta, V), it is defined by phi(g_i) := P_i in V
    2823             :    * with \sum_i P_i . \lambda_i^* = 0, where (\sum n_i g_i)^* :=
    2824             :    * \sum n_i \gamma_i^(-1).
    2825             :    *
    2826             :    * We single out gamma_1 / g_1 (g_oo in Pollack-Stevens paper) and
    2827             :    * write P_{i,j} \lambda_i^* =  Q_{i,j} (\gamma_1 - 1)^* + d_{i,j} y^{k-2}
    2828             :    * where d_{i,j} is a scalar and Q_{i,j} in V; we normalize Q_{i,j} to
    2829             :    * that the coefficient of x^{k-2} is 0.
    2830             :    *
    2831             :    * There exist (s,t) such that d_{s,t} != 0.
    2832             :    * A Q-basis of the (dual) space of modular symbols is given by the
    2833             :    * functions phi_{i,j}, 2 <= i <= d, 1 <= j <= k-1, mapping
    2834             :    *  g_1 -> d_{s,t} Q_{i,j} - d_{i,j} Q_{s,t} + [(i,j)=(s,t)] x^{k-2}
    2835             :    * If i != s
    2836             :    *   g_i -> d_{s,t} P_{i,j}
    2837             :    *   g_s -> - d_{i,j} P_{s,t}
    2838             :    * If i = s, j != t
    2839             :    *   g_i -> d_{s,t} P_{i,j} - d_{i,j} P_{s,t}
    2840             :    * And everything else to 0. Again we normalize the phi_{i,j} such that
    2841             :    * their image has content 1. */
    2842         623 :   monomials = matid(k-1); /* represent the monomials x^{k-2}, ... , y^{k-2} */
    2843         623 :   if (s <= nbE1) /* in E1 */
    2844             :   {
    2845         308 :     P_st = gel(monomials, t);
    2846         308 :     lP_st = gmael(p1, s, t); /* P_{s,t} lambda_s^* */
    2847             :   }
    2848             :   else /* in T2, T31 */
    2849             :   {
    2850         315 :     P_st = gmael(gentor, s - nbE1, t);
    2851         315 :     lP_st = P_st;
    2852             :   }
    2853         623 :   Q_st = doo_decompose(dinv, lP_st, &d_st);
    2854         623 :   basis = cgetg(dim+1, t_VEC);
    2855         623 :   link = cgetg(nbE1 + lg(gentor), t_VEC);
    2856         623 :   gel(link,1) = cgetg(1,t_VECSMALL); /* dummy */
    2857         623 :   n = 1;
    2858        3360 :   for (i = 2; i <= nbE1; i++)
    2859             :   {
    2860        2737 :     GEN L = cgetg(k, t_VECSMALL);
    2861             :     long j;
    2862             :     /* link[i][j] = n gives correspondance between phi_{i,j} and basis[n] */
    2863        2737 :     gel(link,i) = L;
    2864       14000 :     for (j = 1; j < k; j++)
    2865             :     {
    2866       11263 :       GEN lP_ij = gmael(p1, i, j); /* P_{i,j} lambda_i^* */
    2867       11263 :       GEN P_ij = gel(monomials,j);
    2868       11263 :       L[j] = n;
    2869       11263 :       gel(basis, n) = get_phi_ij(i,j,n, s,t, P_st, Q_st, d_st, P_ij, lP_ij, dinv);
    2870       11263 :       n++;
    2871             :     }
    2872             :   }
    2873        1127 :   for (u = 1; u < lg(gentor); u++,i++)
    2874             :   {
    2875         504 :     GEN V = gel(gentor,u);
    2876         504 :     long j, lV = lg(V);
    2877         504 :     GEN L = cgetg(lV, t_VECSMALL);
    2878         504 :     gel(link,i) = L;
    2879        2072 :     for (j = 1; j < lV; j++)
    2880             :     {
    2881        1568 :       GEN lP_ij = gel(V, j); /* P_{i,j} lambda_i^* = P_{i,j} */
    2882        1568 :       GEN P_ij = lP_ij;
    2883        1568 :       L[j] = n;
    2884        1568 :       gel(basis, n) = get_phi_ij(i,j,n, s,t, P_st, Q_st, d_st, P_ij, lP_ij, dinv);
    2885        1568 :       n++;
    2886             :     }
    2887             :   }
    2888         623 :   Inv = cgetg(lg(link), t_VEC);
    2889         623 :   gel(Inv,1) = cgetg(1, t_MAT); /* dummy */
    2890        3864 :   for (i = 2; i < lg(link); i++)
    2891             :   {
    2892        3241 :     GEN M, inv, B = gel(link,i);
    2893        3241 :     long j, lB = lg(B);
    2894        3241 :     if (i == s) { B = vecsplice(B, t); lB--; } /* remove phi_st */
    2895        3241 :     M = cgetg(lB, t_MAT);
    2896       15449 :     for (j = 1; j < lB; j++)
    2897             :     {
    2898       12208 :       GEN phi_ij = gel(basis, B[j]), pols = gel(phi_ij,3);
    2899       12208 :       gel(M, j) = gel(pols, 2); /* phi_ij(g_i) */
    2900             :     }
    2901        3241 :     if (i <= nbE1 && i != s) /* maximal rank k-1 */
    2902        2429 :       inv = ZM_inv_denom(M);
    2903             :     else /* i = s (rank k-2) or from torsion: rank k/3 or k/2 */
    2904         812 :       inv = Qevproj_init(M);
    2905        3241 :     gel(Inv,i) = inv;
    2906             :   }
    2907         623 :   return mkvec3(WN, gen_0, mkvec5(basis, mkvecsmall2(k, dim), mkvecsmall2(s,t),
    2908             :                                   link, Inv));
    2909             : }
    2910             : static GEN
    2911        1386 : add_star(GEN W, long sign)
    2912             : {
    2913        1386 :   GEN s = msstar_i(W);
    2914        1386 :   GEN K = sign? QM_ker_r(gsubgs(s, sign)): cgetg(1,t_MAT);
    2915        1386 :   gel(W,2) = mkvec3(stoi(sign), s, Qevproj_init(K));
    2916        1386 :   return W;
    2917             : }
    2918             : /* WN = msinit_N(N) */
    2919             : static GEN
    2920        1386 : mskinit(ulong N, long k, long sign)
    2921             : {
    2922        1386 :   GEN W, WN = msinit_N(N);
    2923        1386 :   if (N == 1)
    2924             :   {
    2925          14 :     GEN basis, M = RgXV_to_RgM(mfperiodpolbasis(k, 0), k-1);
    2926          14 :     GEN T = cgetg(1, t_VECSMALL), ind = mkvecsmall(1);
    2927          14 :     long i, l = lg(M);
    2928          14 :     basis = cgetg(l, t_VEC);
    2929          14 :     for (i = 1; i < l; i++) gel(basis,i) = mkvec3(T, ind, mkvec(gel(M,i)));
    2930          14 :     W = mkvec3(WN, gen_0, mkvec5(basis, mkvecsmall2(k, l-1), mkvecsmall2(0,0),
    2931             :                                  gen_0, Qevproj_init(M)));
    2932             :   }
    2933             :   else
    2934        1372 :     W = k == 2? mskinit_trivial(WN)
    2935        1372 :               : mskinit_nontrivial(WN, k);
    2936        1386 :   return add_star(W, sign);
    2937             : }
    2938             : GEN
    2939         462 : msinit(GEN N, GEN K, long sign)
    2940             : {
    2941         462 :   pari_sp av = avma;
    2942             :   GEN W;
    2943             :   long k;
    2944         462 :   if (typ(N) != t_INT) pari_err_TYPE("msinit", N);
    2945         462 :   if (typ(K) != t_INT) pari_err_TYPE("msinit", K);
    2946         462 :   k = itos(K);
    2947         462 :   if (k < 2) pari_err_DOMAIN("msinit","k", "<", gen_2,K);
    2948         462 :   if (odd(k)) pari_err_IMPL("msinit [odd weight]");
    2949         462 :   if (signe(N) <= 0) pari_err_DOMAIN("msinit","N", "<=", gen_0,N);
    2950         462 :   W = mskinit(itou(N), k, sign);
    2951         462 :   return gerepilecopy(av, W);
    2952             : }
    2953             : 
    2954             : /* W = msinit, xpm integral modular symbol of weight 2, c t_FRAC
    2955             :  * Return image of <oo->c> */
    2956             : static GEN
    2957        1981 : Q_xpm(GEN W, GEN xpm, GEN c)
    2958             : {
    2959        1981 :   pari_sp av = avma;
    2960             :   GEN v;
    2961        1981 :   W = get_msN(W);
    2962        1981 :   v = init_act_trivial(W);
    2963        1981 :   Q_log_trivial(v, W, c); /* oo -> (a:b), c = a/b */
    2964        1981 :   return gerepileuptoint(av, ZV_zc_mul(xpm, v));
    2965             : }
    2966             : 
    2967             : static GEN
    2968       20146 : eval_single(GEN s, long k, GEN B, long v)
    2969             : {
    2970             :   long i, l;
    2971       20146 :   GEN A = cgetg_copy(s,&l);
    2972       20146 :   for (i=1; i<l; i++) gel(A,i) = ZGl2Q_act_s(gel(B,i), gel(s,i), k);
    2973       20146 :   A = RgV_sum(A);
    2974       20146 :   if (is_vec_t(typ(A))) A = RgV_to_RgX(A, v);
    2975       20146 :   return A;
    2976             : }
    2977             : /* Evaluate symbol s on mspathlog B (= sum p_i g_i, p_i in Z[G]). Allow
    2978             :  * s = t_MAT [ collection of symbols, return a vector ]*/
    2979             : static GEN
    2980       15799 : mseval_by_values(GEN W, GEN s, GEN p, long v)
    2981             : {
    2982       15799 :   long i, l, k = msk_get_weight(W);
    2983             :   GEN A;
    2984       15799 :   if (k == 2)
    2985             :   { /* trivial represention: don't bother with Z[G] */
    2986        3339 :     GEN B = mspathlog_trivial(W,p);
    2987        3332 :     if (typ(s) != t_MAT) return RgV_zc_mul(s,B);
    2988        3262 :     l = lg(s); A = cgetg(l, t_VEC);
    2989        3262 :     for (i = 1; i < l; i++) gel(A,i) = RgV_zc_mul(gel(s,i), B);
    2990             :   }
    2991             :   else
    2992             :   {
    2993       12460 :     GEN B = mspathlog(W,p);
    2994       12460 :     if (typ(s) != t_MAT) return eval_single(s, k, B, v);
    2995         812 :     l = lg(s); A = cgetg(l, t_VEC);
    2996         812 :     for (i = 1; i < l; i++) gel(A,i) = eval_single(gel(s,i), k, B, v);
    2997             :   }
    2998        4074 :   return A;
    2999             : }
    3000             : 
    3001             : /* express symbol on the basis phi_{i,j} */
    3002             : static GEN
    3003       20524 : symtophi(GEN W, GEN s)
    3004             : {
    3005       20524 :   GEN e, basis = msk_get_basis(W);
    3006       20524 :   long i, l = lg(basis);
    3007       20524 :   if (lg(s) != l) pari_err_TYPE("mseval",s);
    3008       20524 :   e = zerovec(ms_get_nbgen(W));
    3009      312914 :   for (i=1; i<l; i++)
    3010             :   {
    3011      292390 :     GEN phi, ind, pols, c = gel(s,i);
    3012             :     long j, m;
    3013      292390 :     if (gequal0(c)) continue;
    3014      122528 :     phi = gel(basis,i);
    3015      122528 :     ind = gel(phi,2); m = lg(ind);
    3016      122528 :     pols = gel(phi,3);
    3017      470470 :     for (j=1; j<m; j++)
    3018             :     {
    3019      347942 :       long t = ind[j];
    3020      347942 :       gel(e,t) = gadd(gel(e,t), gmul(c, gel(pols,j)));
    3021             :     }
    3022             :   }
    3023       20524 :   return e;
    3024             : }
    3025             : /* evaluate symbol s on path p */
    3026             : GEN
    3027       16758 : mseval(GEN W, GEN s, GEN p)
    3028             : {
    3029       16758 :   pari_sp av = avma;
    3030       16758 :   long i, k, l, v = 0;
    3031       16758 :   checkms(W);
    3032       16758 :   k = msk_get_weight(W);
    3033       16758 :   switch(typ(s))
    3034             :   {
    3035             :     case t_VEC: /* values s(g_i) */
    3036           7 :       if (lg(s)-1 != ms_get_nbgen(W)) pari_err_TYPE("mseval",s);
    3037           7 :       if (!p) return gcopy(s);
    3038           0 :       v = gvar(s);
    3039           0 :       break;
    3040             :     case t_COL:
    3041       12663 :       if (msk_get_sign(W))
    3042             :       {
    3043         357 :         GEN star = gel(msk_get_starproj(W), 1);
    3044         357 :         if (lg(star) == lg(s)) s = RgM_RgC_mul(star, s);
    3045             :       }
    3046       12663 :       if (k == 2) /* on the dual basis of (g_i) */
    3047             :       {
    3048         637 :         if (lg(s)-1 != ms_get_nbE1(W)) pari_err_TYPE("mseval",s);
    3049         637 :         if (!p) return gtrans(s);
    3050             :       }
    3051             :       else
    3052       12026 :         s = symtophi(W,s);
    3053       12103 :       break;
    3054             :     case t_MAT:
    3055        4088 :       l = lg(s);
    3056        4088 :       if (!p)
    3057             :       {
    3058           7 :         GEN v = cgetg(l, t_VEC);
    3059           7 :         for (i = 1; i < l; i++) gel(v,i) = mseval(W, gel(s,i), NULL);
    3060           7 :         return v;
    3061             :       }
    3062        4081 :       if (l == 1) return cgetg(1, t_VEC);
    3063        4074 :       if (msk_get_sign(W))
    3064             :       {
    3065          84 :         GEN star = gel(msk_get_starproj(W), 1);
    3066          84 :         if (lg(star) == lgcols(s)) s = RgM_mul(star, s);
    3067             :       }
    3068        4074 :       if (k == 2)
    3069        3262 :       { if (nbrows(s) != ms_get_nbE1(W)) pari_err_TYPE("mseval",s); }
    3070             :       else
    3071             :       {
    3072         812 :         GEN t = cgetg(l, t_MAT);
    3073         812 :         for (i = 1; i < l; i++) gel(t,i) = symtophi(W,gel(s,i));
    3074         812 :         s = t;
    3075             :       }
    3076        4074 :       break;
    3077           0 :     default: pari_err_TYPE("mseval",s);
    3078             :   }
    3079       16177 :   if (p)
    3080       15799 :     s = mseval_by_values(W, s, p, v);
    3081             :   else
    3082             :   {
    3083         378 :     l = lg(s);
    3084        3675 :     for (i = 1; i < l; i++)
    3085             :     {
    3086        3297 :       GEN c = gel(s,i);
    3087        3297 :       if (!isintzero(c)) gel(s,i) = RgV_to_RgX(gel(s,i), v);
    3088             :     }
    3089             :   }
    3090       16170 :   return gerepilecopy(av, s);
    3091             : }
    3092             : 
    3093             : static GEN
    3094         938 : allxpm(GEN W, GEN xpm, long f)
    3095             : {
    3096         938 :   GEN v, L = coprimes_zv(f);
    3097         938 :   long a, nonzero = 0;
    3098         938 :   v = const_vec(f, NULL);
    3099        3647 :   for (a = 1; a <= f; a++)
    3100             :   {
    3101             :     GEN c;
    3102        2709 :     if (!L[a]) continue;
    3103        1981 :     c = Q_xpm(W, xpm, sstoQ(a, f));
    3104        1981 :     if (!gequal0(c)) { gel(v,a) = c; nonzero = 1; }
    3105             :   }
    3106         938 :   return nonzero? v: NULL;
    3107             : }
    3108             : /* \sum_{a mod f} chi(a) x(a/f) */
    3109             : static GEN
    3110         441 : seval(GEN G, GEN chi, GEN vx)
    3111             : {
    3112         441 :   GEN vZ, T, s = gen_0, go = zncharorder(G,chi);
    3113         441 :   long i, l = lg(vx), o = itou(go);
    3114         441 :   T = polcyclo(o,0);
    3115         441 :   vZ = mkvec2(RgXQ_powers(RgX_rem(pol_x(0), T), o-1, T), go);
    3116        1813 :   for (i = 1; i < l; i++)
    3117             :   {
    3118        1372 :     GEN x = gel(vx,i);
    3119        1372 :     if (x) s = gadd(s, gmul(x, znchareval(G, chi, utoi(i), vZ)));
    3120             :   }
    3121         441 :   return gequal0(s)? NULL: poleval(s, rootsof1u_cx(o, DEFAULTPREC));
    3122             : }
    3123             : 
    3124             : static long
    3125         189 : nb_components(GEN E) { return signe(ell_get_disc(E)) > 0? 2: 1; }
    3126             : /* E minimal */
    3127             : static GEN
    3128         609 : ellperiod(GEN E, long s)
    3129             : {
    3130         609 :   GEN w = ellR_omega(E,DEFAULTPREC);
    3131         609 :   if (s == 1)
    3132         420 :     w = gel(w,1);
    3133         189 :   else if (nb_components(E) == 2)
    3134          77 :     w = gneg(gel(w,2));
    3135             :   else
    3136         112 :     w = mkcomplex(gen_0, gneg(gmul2n(imag_i(gel(w,2)), 1)));
    3137         609 :   return w;
    3138             : }
    3139             : 
    3140             : /* Let W = msinit(conductor(E), 2), xpm an integral modular symbol with the same
    3141             :  * eigenvalues as L_E. There exist a unique C such that
    3142             :  *   C*L(E,(D/.),1)_{xpm} = L(E,(D/.),1) / w1(E_D) != 0, for all D fundamental,
    3143             :  * sign(D) = s, and such that E_D has rank 0. Return C * ellperiod(E,s) */
    3144             : static GEN
    3145         441 : ell_get_Cw(GEN LE, GEN W, GEN xpm, long s)
    3146             : {
    3147         441 :   long f, NE = ms_get_N(W);
    3148         441 :   const long bit = 64;
    3149             : 
    3150        1372 :   for (f = 1;; f++)
    3151             :   { /* look for chi with conductor f coprime to N(E) and parity s
    3152             :      * such that L(E,chi,1) != 0 */
    3153        1372 :     pari_sp av = avma;
    3154             :     GEN vchi, vx, G;
    3155             :     long l, i;
    3156        1372 :     if ((f & 3) == 2 || cgcd(NE,f) != 1) continue;
    3157         938 :     vx = allxpm(W, xpm, f); if (!vx) continue;
    3158         441 :     G = znstar0(utoipos(f),1);
    3159         441 :     vchi = chargalois(G,NULL); l = lg(vchi);
    3160         763 :     for (i = 1; i < l; i++)
    3161             :     {
    3162         763 :       pari_sp av2 = avma;
    3163         763 :       GEN tau, z, S, L, chi = gel(vchi,i);
    3164         763 :       long o = zncharisodd(G,chi);
    3165         763 :       if ((s > 0 && o) || (s < 0 && !o)
    3166         553 :           || itos(zncharconductor(G, chi)) != f) continue;
    3167         441 :       S = seval(G, chi, vx);
    3168         441 :       if (!S) { avma = av2; continue; }
    3169             : 
    3170         441 :       L = lfuntwist(LE, mkvec2(G, zncharconj(G,chi)));
    3171         441 :       z = lfun(L, gen_1, bit);
    3172         441 :       tau = znchargauss(G, chi, gen_1, bit);
    3173         882 :       return gdiv(gmul(z, tau), S); /* C * w */
    3174             :     }
    3175           0 :     avma = av;
    3176         931 :   }
    3177             : }
    3178             : static GEN
    3179         357 : ell_get_scale(GEN LE, GEN W, long sign, GEN x)
    3180             : {
    3181         357 :   if (sign)
    3182         273 :     return ell_get_Cw(LE, W, gel(x,1), sign);
    3183             :   else
    3184             :   {
    3185          84 :     GEN Cwp = ell_get_Cw(LE, W, gel(x,1), 1);
    3186          84 :     GEN Cwm = ell_get_Cw(LE, W, gel(x,2),-1);
    3187          84 :     return mkvec2(Cwp, Cwm);
    3188             :   }
    3189             : }
    3190             : /* E minimal */
    3191             : static GEN
    3192         441 : msfromell_scale(GEN x, GEN Cw, GEN E, long s)
    3193             : {
    3194         441 :   GEN B = int2n(32);
    3195         441 :   if (s)
    3196             :   {
    3197         273 :     GEN C = gdiv(Cw, ellperiod(E,s));
    3198         273 :     return ZC_Q_mul(gel(x,1), bestappr(C,B));
    3199             :   }
    3200             :   else
    3201             :   {
    3202         168 :     GEN xp = gel(x,1), Cp = gdiv(gel(Cw,1), ellperiod(E, 1)), L;
    3203         168 :     GEN xm = gel(x,2), Cm = gdiv(gel(Cw,2), ellperiod(E,-1));
    3204         168 :     xp = ZC_Q_mul(xp, bestappr(Cp,B));
    3205         168 :     xm = ZC_Q_mul(xm, bestappr(Cm,B));
    3206         168 :     if (signe(ell_get_disc(E)) > 0)
    3207          70 :       L = mkmat2(xp, xm); /* E(R) has 2 connected components */
    3208             :     else
    3209          98 :       L = mkmat2(gsub(xp,xm), gmul2n(xm,1));
    3210         168 :     return mkvec3(xp, xm, L);
    3211             :   }
    3212             : }
    3213             : /* v != 0 */
    3214             : static GEN
    3215         443 : Flc_normalize(GEN v, ulong p)
    3216             : {
    3217         443 :   long i, l = lg(v);
    3218         782 :   for (i = 1; i < l; i++)
    3219         782 :     if (v[i])
    3220             :     {
    3221         443 :       if (v[i] != 1) v = Flv_Fl_div(v, v[i], p);
    3222         443 :       return v;
    3223             :     }
    3224           0 :   return NULL;
    3225             : }
    3226             : /* K \cap Ker M  [F_l vector spaces]. K = NULL means full space */
    3227             : static GEN
    3228         394 : msfromell_ker(GEN K, GEN M, ulong l)
    3229             : {
    3230         394 :   GEN B, Ml = ZM_to_Flm(M, l);
    3231         394 :   if (K) Ml = Flm_mul(Ml, K, l);
    3232         394 :   B = Flm_ker(Ml, l);
    3233         394 :   if (!K) K = B;
    3234          36 :   else if (lg(B) < lg(K))
    3235          36 :     K = Flm_mul(K, B, l);
    3236         394 :   return K;
    3237             : }
    3238             : /* K = \cap_p Ker(T_p - a_p), 2-dimensional. Set *xl to the 1-dimensional
    3239             :  * Fl-basis  such that star . xl = sign . xl if sign != 0 and
    3240             :  * star * xl[1] = xl[1]; star * xl[2] = -xl[2] if sign = 0 */
    3241             : static void
    3242         358 : msfromell_l(GEN *pxl, GEN K, GEN star, long sign, ulong l)
    3243             : {
    3244         358 :   GEN s = ZM_to_Flm(star, l);
    3245         358 :   GEN a = gel(K,1), Sa = Flm_Flc_mul(s,a,l);
    3246         358 :   GEN b = gel(K,2);
    3247         358 :   GEN t = Flv_add(a,Sa,l), xp, xm;
    3248         358 :   if (zv_equal0(t))
    3249             :   {
    3250          14 :     xm = a;
    3251          14 :     xp = Flv_add(b,Flm_Flc_mul(s,b,l), l);
    3252             :   }
    3253             :   else
    3254             :   {
    3255         344 :     xp = t; t = Flv_sub(a, Sa, l);
    3256         344 :     xm = zv_equal0(t)? Flv_sub(b, Flm_Flc_mul(s,b,l), l): t;
    3257             :   }
    3258             :   /* xp = 0 on Im(S - 1), xm = 0 on Im(S + 1) */
    3259         358 :   if (sign > 0)
    3260         252 :     *pxl = mkmat(Flc_normalize(xp, l));
    3261         106 :   else if (sign < 0)
    3262          21 :     *pxl = mkmat(Flc_normalize(xm, l));
    3263             :   else
    3264          85 :     *pxl = mkmat2(Flc_normalize(xp, l), Flc_normalize(xm, l));
    3265         358 : }
    3266             : /* return a primitive symbol */
    3267             : static GEN
    3268         358 : msfromell_ratlift(GEN x, GEN q)
    3269             : {
    3270         358 :   GEN B = sqrti(shifti(q,-1));
    3271         358 :   GEN r = FpM_ratlift(x, q, B, B, NULL);
    3272         358 :   if (r) r = Q_primpart(r);
    3273         358 :   return r;
    3274             : }
    3275             : static int
    3276         358 : msfromell_check(GEN x, GEN vT, GEN star, long sign)
    3277             : {
    3278             :   long i, l;
    3279             :   GEN sx;
    3280         358 :   if (!x) return 0;
    3281         357 :   l = lg(vT);
    3282         749 :   for (i = 1; i < l; i++)
    3283             :   {
    3284         392 :     GEN T = gel(vT,i);
    3285         392 :     if (!gequal0(ZM_mul(T, x))) return 0; /* fail */
    3286             :   }
    3287         357 :   sx = ZM_mul(star,x);
    3288         357 :   if (sign)
    3289         273 :     return ZV_equal(gel(sx,1), sign > 0? gel(x,1): ZC_neg(gel(x,1)));
    3290             :   else
    3291          84 :     return ZV_equal(gel(sx,1),gel(x,1)) && ZV_equal(gel(sx,2),ZC_neg(gel(x,2)));
    3292             : }
    3293             : GEN
    3294         357 : msfromell(GEN E0, long sign)
    3295             : {
    3296         357 :   pari_sp av = avma, av2;
    3297         357 :   GEN T, Cw, E, NE, star, q, vT, xl, xr, W, x = NULL, K = NULL;
    3298             :   long lE, single;
    3299             :   ulong p, l, N;
    3300             :   forprime_t S, Sl;
    3301             : 
    3302         357 :   if (typ(E0) != t_VEC) pari_err_TYPE("msfromell",E0);
    3303         357 :   lE = lg(E0);
    3304         357 :   if (lE == 1) return cgetg(1,t_VEC);
    3305         357 :   single = (typ(gel(E0,1)) != t_VEC);
    3306         357 :   E = single ? E0: gel(E0,1);
    3307         357 :   NE = ellQ_get_N(E);
    3308             :   /* must make it integral for ellap; we have minimal model at hand */
    3309         357 :   T = obj_check(E, Q_MINIMALMODEL); if (lg(T) != 2) E = gel(T,3);
    3310         357 :   N = itou(NE); av2 = avma;
    3311         357 :   W = gerepilecopy(av2, mskinit(N,2,0));
    3312         357 :   star = msk_get_star(W);
    3313         357 :   init_modular_small(&Sl);
    3314             :   /* loop for p <= count_Manin_symbols(N) / 6 would be enough */
    3315         357 :   (void)u_forprime_init(&S, 2, ULONG_MAX);
    3316         357 :   vT = cgetg(1, t_VEC);
    3317         357 :   l = u_forprime_next(&Sl);
    3318         357 :   while( (p = u_forprime_next(&S)) )
    3319             :   {
    3320             :     GEN M;
    3321         476 :     if (N % p == 0) continue;
    3322         392 :     av2 = avma;
    3323         392 :     M = RgM_Rg_sub_shallow(mshecke_i(W, p), ellap(E, utoipos(p)));
    3324         392 :     M = gerepilecopy(av2, M);
    3325         392 :     vT = shallowconcat(vT, mkvec(M)); /* for certification at the end */
    3326         392 :     K = msfromell_ker(K, M, l);
    3327         392 :     if (lg(K) == 3) break;
    3328             :   }
    3329         357 :   if (!p) pari_err_BUG("msfromell: ran out of primes");
    3330             : 
    3331             :   /* mod one l should be enough */
    3332         357 :   msfromell_l(&xl, K, star, sign, l);
    3333         357 :   x = ZM_init_CRT(xl, l);
    3334         357 :   q = utoipos(l);
    3335         357 :   xr = msfromell_ratlift(x, q);
    3336             :   /* paranoia */
    3337         715 :   while (!msfromell_check(xr, vT, star, sign) && (l = u_forprime_next(&Sl)) )
    3338             :   {
    3339           1 :     GEN K = NULL;
    3340           1 :     long i, lvT = lg(vT);
    3341           2 :     for (i = 1; i < lvT; i++)
    3342             :     {
    3343           2 :       K = msfromell_ker(K, gel(vT,i), l);
    3344           2 :       if (lg(K) == 3) break;
    3345             :     }
    3346           1 :     if (i >= lvT) { x = NULL; continue; }
    3347           1 :     msfromell_l(&xl, K, star, sign, l);
    3348           1 :     ZM_incremental_CRT(&x, xl, &q, l);
    3349           1 :     xr = msfromell_ratlift(x, q);
    3350             :   }
    3351             :   /* linear form = 0 on all Im(Tp - ap) and Im(S - sign) if sign != 0 */
    3352         357 :   Cw = ell_get_scale(lfuncreate(E), W, sign, xr);
    3353         357 :   if (single)
    3354         329 :     x = msfromell_scale(xr, Cw, E, sign);
    3355             :   else
    3356             :   { /* assume all E0[i] isogenous, given by minimal models */
    3357          28 :     GEN v = cgetg(lE, t_VEC);
    3358             :     long i;
    3359          28 :     for (i=1; i<lE; i++) gel(v,i) = msfromell_scale(xr, Cw, gel(E0,i), sign);
    3360          28 :     x = v;
    3361             :   }
    3362         357 :   return gerepilecopy(av, mkvec2(W, x));
    3363             : }
    3364             : 
    3365             : GEN
    3366          21 : msfromhecke(GEN W, GEN v, GEN H)
    3367             : {
    3368          21 :   pari_sp av = avma;
    3369          21 :   long i, l = lg(v);
    3370          21 :   GEN K = NULL;
    3371          21 :   checkms(W);
    3372          21 :   if (typ(v) != t_VEC) pari_err_TYPE("msfromhecke",v);
    3373          49 :   for (i = 1; i < l; i++)
    3374             :   {
    3375          28 :     GEN K2, T, p, P, c = gel(v,i);
    3376          28 :     if (typ(c) != t_VEC || lg(c) != 3) pari_err_TYPE("msfromhecke",v);
    3377          28 :     p = gel(c,1);
    3378          28 :     if (typ(p) != t_INT) pari_err_TYPE("msfromhecke",v);
    3379          28 :     P = gel(c,2);
    3380          28 :     switch(typ(P))
    3381             :     {
    3382             :       case t_INT:
    3383          21 :         P = deg1pol_shallow(gen_1, negi(P), 0);
    3384          21 :         break;
    3385             :       case t_POL:
    3386           7 :         if (RgX_is_ZX(P)) break;
    3387             :       default:
    3388           0 :         pari_err_TYPE("msfromhecke",v);
    3389             :     };
    3390          28 :     T = mshecke(W, itos(p), H);
    3391          28 :     T = Q_primpart(RgX_RgM_eval(P, T));
    3392          28 :     if (K) T = ZM_mul(T,K);
    3393          28 :     K2 = ZM_ker(T);
    3394          28 :     if (!K) K = K2;
    3395           7 :     else if (lg(K2) < lg(K)) K = ZM_mul(K,K2);
    3396             :   }
    3397          21 :   return gerepilecopy(av, K);
    3398             : }
    3399             : 
    3400             : /* OVERCONVERGENT MODULAR SYMBOLS */
    3401             : 
    3402             : static GEN
    3403        2933 : mspadic_get_Wp(GEN W) { return gel(W,1); }
    3404             : static GEN
    3405         483 : mspadic_get_Tp(GEN W) { return gel(W,2); }
    3406             : static GEN
    3407         483 : mspadic_get_bin(GEN W) { return gel(W,3); }
    3408             : static GEN
    3409         476 : mspadic_get_actUp(GEN W) { return gel(W,4); }
    3410             : static GEN
    3411         476 : mspadic_get_q(GEN W) { return gel(W,5); }
    3412             : static long
    3413        1456 : mspadic_get_p(GEN W) { return gel(W,6)[1]; }
    3414             : static long
    3415        1211 : mspadic_get_n(GEN W) { return gel(W,6)[2]; }
    3416             : static long
    3417         161 : mspadic_get_flag(GEN W) { return gel(W,6)[3]; }
    3418             : static GEN
    3419         483 : mspadic_get_M(GEN W) { return gel(W,7); }
    3420             : static GEN
    3421         483 : mspadic_get_C(GEN W) { return gel(W,8); }
    3422             : static long
    3423         973 : mspadic_get_weight(GEN W) { return msk_get_weight(mspadic_get_Wp(W)); }
    3424             : 
    3425             : void
    3426         980 : checkmspadic(GEN W)
    3427             : {
    3428         980 :   if (typ(W) != t_VEC || lg(W) != 9) pari_err_TYPE("checkmspadic",W);
    3429         980 :   checkms(mspadic_get_Wp(W));
    3430         980 : }
    3431             : 
    3432             : /* f in M_2(Z) \cap GL_2(Q), p \nmid a [ and for the result to mean anything
    3433             :  * p | c, but not needed here]. Return the matrix M in M_D(Z), D = M+k-1
    3434             :  * such that, if v = \int x^i d mu, i < D, is a vector of D moments of mu,
    3435             :  * then M * v is the vector of moments of mu | f  mod p^D */
    3436             : static GEN
    3437      276073 : moments_act_i(struct m_act *S, GEN f)
    3438             : {
    3439      276073 :   long j, k = S->k, D = S->dim;
    3440      276073 :   GEN a = gcoeff(f,1,1), b = gcoeff(f,1,2);
    3441      276073 :   GEN c = gcoeff(f,2,1), d = gcoeff(f,2,2);
    3442      276073 :   GEN u,z,C, q = S->q, mat = cgetg(D+1, t_MAT);
    3443             : 
    3444      276073 :   a = modii(a,q);
    3445      276073 :   z = FpX_powu(deg1pol(c,a,0), k-2, q); /* (a+cx)^(k-2) */
    3446             :   /* u := (b+dx) / (a+cx) mod (q,x^D) = (b/a +d/a*x) / (1 - (-c/a)*x) */
    3447      276073 :   if (!equali1(a))
    3448             :   {
    3449      271229 :     GEN ai = Fp_inv(a,q);
    3450      271229 :     b = Fp_mul(b,ai,q);
    3451      271229 :     c = Fp_mul(c,ai,q);
    3452      271229 :     d = Fp_mul(d,ai,q);
    3453             :   }
    3454      276073 :   u = cgetg(D+2,t_POL); u[1] = evalsigne(1)|evalvarn(0);
    3455      276073 :   gel(u, 2) = gen_1;
    3456      276073 :   gel(u, 3) = C = Fp_neg(c,q);
    3457      276073 :   for (j = 4; j < D+2; j++) gel(u,j) = Fp_mul(gel(u,j-1), C, q);
    3458      276073 :   u = FpX_red(RgXn_mul(deg1pol(d,b,0), u, D), q);
    3459     2369024 :   for (j = 1; j <= D; j++)
    3460             :   {
    3461     2092951 :     gel(mat,j) = RgX_to_RgC(z, D); /* (a+cx)^(k-2) * ((b+dx)/(a+cx))^(j-1) */
    3462     2092951 :     if (j != D) z = FpX_red(RgXn_mul(z, u, D), q);
    3463             :   }
    3464      276073 :   return shallowtrans(mat);
    3465             : }
    3466             : static GEN
    3467      275611 : moments_act(struct m_act *S, GEN f)
    3468      275611 : { pari_sp av = avma; return gerepilecopy(av, moments_act_i(S,f)); }
    3469             : static GEN
    3470         483 : init_moments_act(GEN W, long p, long n, GEN q, GEN v)
    3471             : {
    3472             :   struct m_act S;
    3473         483 :   long k = msk_get_weight(W);
    3474         483 :   S.p = p;
    3475         483 :   S.k = k;
    3476         483 :   S.q = q;
    3477         483 :   S.dim = n+k-1;
    3478         483 :   S.act = &moments_act; return init_dual_act(v,W,W,&S);
    3479             : }
    3480             : 
    3481             : static void
    3482        6762 : clean_tail(GEN phi, long c, GEN q)
    3483             : {
    3484        6762 :   long a, l = lg(phi);
    3485      214438 :   for (a = 1; a < l; a++)
    3486             :   {
    3487      207676 :     GEN P = FpV_red(gel(phi,a), q); /* phi(G_a) = vector of moments */
    3488      207676 :     long j, lP = lg(P);
    3489      207676 :     for (j = c; j < lP; j++) gel(P,j) = gen_0; /* reset garbage to 0 */
    3490      207676 :     gel(phi,a) = P;
    3491             :   }
    3492        6762 : }
    3493             : /* concat z to all phi[i] */
    3494             : static GEN
    3495         630 : concat2(GEN phi, GEN z)
    3496             : {
    3497             :   long i, l;
    3498         630 :   GEN v = cgetg_copy(phi,&l);
    3499         630 :   for (i = 1; i < l; i++) gel(v,i) = shallowconcat(gel(phi,i), z);
    3500         630 :   return v;
    3501             : }
    3502             : static GEN
    3503         630 : red_mod_FilM(GEN phi, ulong p, long k, long flag)
    3504             : {
    3505             :   long a, l;
    3506         630 :   GEN den = gen_1, v = cgetg_copy(phi, &l);
    3507         630 :   if (flag)
    3508             :   {
    3509         343 :     phi = Q_remove_denom(phi, &den);
    3510         343 :     if (!den) { den = gen_1; flag = 0; }
    3511             :   }
    3512       29386 :   for (a = 1; a < l; a++)
    3513             :   {
    3514       28756 :     GEN P = gel(phi,a), q = den;
    3515             :     long j;
    3516      207676 :     for (j = lg(P)-1; j >= k+1; j--)
    3517             :     {
    3518      178920 :       q = muliu(q,p);
    3519      178920 :       gel(P,j) = modii(gel(P,j),q);
    3520             :     }
    3521       28756 :     q = muliu(q,p);
    3522       93380 :     for (     ; j >= 1; j--)
    3523       64624 :       gel(P,j) = modii(gel(P,j),q);
    3524       28756 :     gel(v,a) = P;
    3525             :   }
    3526         630 :   if (flag) v = gdiv(v, den);
    3527         630 :   return v;
    3528             : }
    3529             : 
    3530             : /* denom(C) | p^(2(k-1) - v_p(ap)) */
    3531             : static GEN
    3532         154 : oms_dim2(GEN W, GEN phi, GEN C, GEN ap)
    3533             : {
    3534         154 :   long t, i, k = mspadic_get_weight(W);
    3535         154 :   long p = mspadic_get_p(W), n = mspadic_get_n(W);
    3536         154 :   GEN phi1 = gel(phi,1), phi2 = gel(phi,2);
    3537         154 :   GEN v, q = mspadic_get_q(W);
    3538         154 :   GEN act = mspadic_get_actUp(W);
    3539             : 
    3540         154 :   t = signe(ap)? Z_lval(ap,p) : k-1;
    3541         154 :   phi1 = concat2(phi1, zerovec(n));
    3542         154 :   phi2 = concat2(phi2, zerovec(n));
    3543        2107 :   for (i = 1; i <= n; i++)
    3544             :   {
    3545        1953 :     phi1 = dual_act(k-1, act, phi1);
    3546        1953 :     phi1 = dual_act(k-1, act, phi1);
    3547        1953 :     clean_tail(phi1, k + i*t, q);
    3548             : 
    3549        1953 :     phi2 = dual_act(k-1, act, phi2);
    3550        1953 :     phi2 = dual_act(k-1, act, phi2);
    3551        1953 :     clean_tail(phi2, k + i*t, q);
    3552             :   }
    3553         154 :   C = gpowgs(C,n);
    3554         154 :   v = RgM_RgC_mul(C, mkcol2(phi1,phi2));
    3555         154 :   phi1 = red_mod_FilM(gel(v,1), p, k, 1);
    3556         154 :   phi2 = red_mod_FilM(gel(v,2), p, k, 1);
    3557         154 :   return mkvec2(phi1,phi2);
    3558             : }
    3559             : 
    3560             : /* flag = 0 iff alpha is a p-unit */
    3561             : static GEN
    3562         322 : oms_dim1(GEN W, GEN phi, GEN alpha, long flag)
    3563             : {
    3564         322 :   long i, k = mspadic_get_weight(W);
    3565         322 :   long p = mspadic_get_p(W), n = mspadic_get_n(W);
    3566         322 :   GEN q = mspadic_get_q(W);
    3567         322 :   GEN act = mspadic_get_actUp(W);
    3568         322 :   phi = concat2(phi, zerovec(n));
    3569        3178 :   for (i = 1; i <= n; i++)
    3570             :   {
    3571        2856 :     phi = dual_act(k-1, act, phi);
    3572        2856 :     clean_tail(phi, k + i, q);
    3573             :   }
    3574         322 :   phi = gmul(lift_shallow(gpowgs(alpha,n)), phi);
    3575         322 :   phi = red_mod_FilM(phi, p, k, flag);
    3576         322 :   return mkvec(phi);
    3577             : }
    3578             : 
    3579             : /* lift polynomial P in RgX[X,Y]_{k-2} to a distribution \mu such that
    3580             :  * \int (Y - X z)^(k-2) d\mu(z) = P(X,Y)
    3581             :  * Return the t_VEC of k-1 first moments of \mu: \int z^i d\mu(z), 0<= i < k-1.
    3582             :  *   \sum_j (-1)^(k-2-j) binomial(k-2,j) Y^j \int z^(k-2-j) d\mu(z) = P(1,Y)
    3583             :  * Input is P(1,Y), bin = vecbinomial(k-2): bin[j] = binomial(k-2,j-1) */
    3584             : static GEN
    3585       38626 : RgX_to_moments(GEN P, GEN bin)
    3586             : {
    3587       38626 :   long j, k = lg(bin);
    3588             :   GEN Pd, Bd;
    3589       38626 :   if (typ(P) != t_POL) P = scalarpol(P,0);
    3590       38626 :   P = RgX_to_RgC(P, k-1); /* deg <= k-2 */
    3591       38626 :   settyp(P, t_VEC);
    3592       38626 :   Pd = P+1;  /* Pd[i] = coeff(P,i) */
    3593       38626 :   Bd = bin+1;/* Bd[i] = binomial(k-2,i) */
    3594       46249 :   for (j = 1; j < k-2; j++)
    3595             :   {
    3596        7623 :     GEN c = gel(Pd,j);
    3597        7623 :     if (odd(j)) c = gneg(c);
    3598        7623 :     gel(Pd,j) = gdiv(c, gel(Bd,j));
    3599             :   }
    3600       38626 :   return vecreverse(P);
    3601             : }
    3602             : static GEN
    3603         882 : RgXC_to_moments(GEN v, GEN bin)
    3604             : {
    3605             :   long i, l;
    3606         882 :   GEN w = cgetg_copy(v,&l);
    3607         882 :   for (i=1; i<l; i++) gel(w,i) = RgX_to_moments(gel(v,i),bin);
    3608         882 :   return w;
    3609             : }
    3610             : 
    3611             : /* W an mspadic, assume O[2] is integral, den is the cancelled denominator
    3612             :  * or NULL, L = log(path)^* in sparse form */
    3613             : static GEN
    3614        2954 : omseval_int(struct m_act *S, GEN PHI, GEN L, hashtable *H)
    3615             : {
    3616             :   long i, l;
    3617        2954 :   GEN v = cgetg_copy(PHI, &l);
    3618        2954 :   ZGl2QC_to_act(S, L, H); /* as operators on V */
    3619        6286 :   for (i = 1; i < l; i++)
    3620             :   {
    3621        3332 :     GEN T = dense_act_col(L, gel(PHI,i));
    3622        3332 :     gel(v,i) = T? FpC_red(T,S->q): zerocol(S->dim);
    3623             :   }
    3624        2954 :   return v;
    3625             : }
    3626             : 
    3627             : GEN
    3628          14 : msomseval(GEN W, GEN phi, GEN path)
    3629             : {
    3630             :   struct m_act S;
    3631          14 :   pari_sp av = avma;
    3632             :   GEN v, Wp;
    3633             :   long n, vden;
    3634          14 :   checkmspadic(W);
    3635          14 :   if (typ(phi) != t_COL || lg(phi) != 4)  pari_err_TYPE("msomseval",phi);
    3636          14 :   vden = itos(gel(phi,2));
    3637          14 :   phi = gel(phi,1);
    3638          14 :   n = mspadic_get_n(W);
    3639          14 :   Wp= mspadic_get_Wp(W);
    3640          14 :   S.k = mspadic_get_weight(W);
    3641          14 :   S.p = mspadic_get_p(W);
    3642          14 :   S.q = powuu(S.p, n+vden);
    3643          14 :   S.dim = n + S.k - 1;
    3644          14 :   S.act = &moments_act;
    3645          14 :   path = path_to_M2(path);
    3646          14 :   v = omseval_int(&S, phi, M2_logf(Wp,path,NULL), NULL);
    3647          14 :   return gerepilecopy(av, v);
    3648             : }
    3649             : /* W = msinit(N,k,...); if flag < 0 or flag >= k-1, allow all symbols;
    3650             :  * else commit to v_p(a_p) <= flag (ordinary if flag = 0)*/
    3651             : GEN
    3652         490 : mspadicinit(GEN W, long p, long n, long flag)
    3653             : {
    3654         490 :   pari_sp av = avma;
    3655             :   long a, N, k;
    3656             :   GEN P, C, M, bin, Wp, Tp, q, pn, actUp, teich, pas;
    3657             : 
    3658         490 :   checkms(W);
    3659         490 :   N = ms_get_N(W);
    3660         490 :   k = msk_get_weight(W);
    3661         490 :   if (flag < 0) flag = 1; /* worst case */
    3662         357 :   else if (flag >= k) flag = k-1;
    3663             : 
    3664         490 :   bin = vecbinomial(k-2);
    3665         490 :   Tp = mshecke(W, p, NULL);
    3666         490 :   if (N % p == 0)
    3667             :   {
    3668          91 :     if ((N/p) % p == 0) pari_err_IMPL("mspadicinit when p^2 | N");
    3669             :     /* a_p != 0 */
    3670          84 :     Wp = W;
    3671          84 :     M = gen_0;
    3672          84 :     flag = (k-2) / 2; /* exact valuation */
    3673             :     /* will multiply by matrix with denominator p^(k-2)/2 in mspadicint.
    3674             :      * Except if p = 2 (multiply by alpha^2) */
    3675          84 :     if (p == 2) n += k-2; else n += (k-2)/2;
    3676          84 :     pn = powuu(p,n);
    3677             :     /* For accuracy mod p^n, oms_dim1 require p^(k/2*n) */
    3678          84 :     q = powiu(pn, k/2);
    3679             :   }
    3680             :   else
    3681             :   { /* p-stabilize */
    3682         399 :     long s = msk_get_sign(W);
    3683             :     GEN M1, M2;
    3684             : 
    3685         399 :     Wp = mskinit(N*p, k, s);
    3686         399 :     M1 = getMorphism(W, Wp, mkvec(mat2(1,0,0,1)));
    3687         399 :     M2 = getMorphism(W, Wp, mkvec(mat2(p,0,0,1)));
    3688         399 :     if (s)
    3689             :     {
    3690         147 :       GEN SW = msk_get_starproj(W), SWp = msk_get_starproj(Wp);
    3691         147 :       M1 = Qevproj_apply2(M1, SW, SWp);
    3692         147 :       M2 = Qevproj_apply2(M2, SW, SWp);
    3693             :     }
    3694         399 :     M = mkvec2(M1,M2);
    3695         399 :     n += Z_lval(Q_denom(M), p); /*den. introduced by p-stabilization*/
    3696             :     /* in supersingular case: will multiply by matrix with denominator p^k
    3697             :      * in mspadicint. Except if p = 2 (multiply by alpha^2) */
    3698         399 :     if (flag) { if (p == 2) n += 2*k-2; else n += k; }
    3699         399 :     pn = powuu(p,n);
    3700             :     /* For accuracy mod p^n, supersingular require p^((2k-1-v_p(a_p))*n) */
    3701         399 :     if (flag) /* k-1 also takes care of a_p = 0. Worst case v_p(a_p) = flag */
    3702         231 :       q = powiu(pn, 2*k-1 - flag);
    3703             :     else
    3704         168 :       q = pn;
    3705             :   }
    3706         483 :   actUp = init_moments_act(Wp, p, n, q, Up_matrices(p));
    3707             : 
    3708         483 :   if (p == 2) C = gen_0;
    3709             :   else
    3710             :   {
    3711         427 :     pas = matpascal(n);
    3712         427 :     teich = teichmullerinit(p, n+1);
    3713         427 :     P = gpowers(utoipos(p), n);
    3714         427 :     C = cgetg(p, t_VEC);
    3715        2317 :     for (a = 1; a < p; a++)
    3716             :     { /* powb[j+1] = ((a - w(a)) / p)^j mod p^n */
    3717        1890 :       GEN powb = Fp_powers(diviuexact(subui(a, gel(teich,a)), p), n, pn);
    3718        1890 :       GEN Ca = cgetg(n+2, t_VEC);
    3719        1890 :       long j, r, ai = Fl_inv(a, p); /* a^(-1) */
    3720        1890 :       gel(C,a) = Ca;
    3721       22134 :       for (j = 0; j <= n; j++)
    3722             :       {
    3723       20244 :         GEN Caj = cgetg(j+2, t_VEC);
    3724       20244 :         GEN atij = gel(teich, Fl_powu(ai,j,p));/* w(a)^(-j) = w(a^(-j) mod p) */
    3725       20244 :         gel(Ca,j+1) = Caj;
    3726      158200 :         for (r = 0; r <= j; r++)
    3727             :         {
    3728      137956 :           GEN c = Fp_mul(gcoeff(pas,j+1,r+1), gel(powb, j-r+1), pn);
    3729      137956 :           c = Fp_mul(c,atij,pn); /* binomial(j,r)*b^(j-r)*w(a)^(-j) mod p^n */
    3730      137956 :           gel(Caj,r+1) = mulii(c, gel(P,j+1)); /* p^j * c mod p^(n+j) */
    3731             :         }
    3732             :       }
    3733             :     }
    3734             :   }
    3735         483 :   return gerepilecopy(av, mkvecn(8, Wp,Tp, bin, actUp, q,
    3736             :                                  mkvecsmall3(p,n,flag), M, C));
    3737             : }
    3738             : 
    3739             : #if 0
    3740             : /* assume phi an ordinary OMS */
    3741             : static GEN
    3742             : omsactgl2(GEN W, GEN phi, GEN M)
    3743             : {
    3744             :   GEN q, Wp, act;
    3745             :   long p, k, n;
    3746             :   checkmspadic(W);
    3747             :   Wp = mspadic_get_Wp(W);
    3748             :   p = mspadic_get_p(W);
    3749             :   k = mspadic_get_weight(W);
    3750             :   n = mspadic_get_n(W);
    3751             :   q = mspadic_get_q(W);
    3752             :   act = init_moments_act(Wp, p, n, q, M);
    3753             :   phi = gel(phi,1);
    3754             :   return dual_act(k-1, act, gel(phi,1));
    3755             : }
    3756             : #endif
    3757             : 
    3758             : static GEN
    3759         483 : eigenvalue(GEN T, GEN x)
    3760             : {
    3761         483 :   long i, l = lg(x);
    3762         637 :   for (i = 1; i < l; i++)
    3763         637 :     if (!isintzero(gel(x,i))) break;
    3764         483 :   if (i == l) pari_err_DOMAIN("mstooms", "phi", "=", gen_0, x);
    3765         483 :   return gdiv(RgMrow_RgC_mul(T,x,i), gel(x,i));
    3766             : }
    3767             : 
    3768             : /* p coprime to ap, return unit root of x^2 - ap*x + p^(k-1), accuracy p^n */
    3769             : GEN
    3770         266 : mspadic_unit_eigenvalue(GEN ap, long k, GEN p, long n)
    3771             : {
    3772         266 :   GEN sqrtD, D = subii(sqri(ap), shifti(powiu(p,k-1),2));
    3773         266 :   if (absequaliu(p,2))
    3774             :   {
    3775           7 :     n++; sqrtD = Zp_sqrt(D, p, n);
    3776           7 :     if (mod4(sqrtD) != mod4(ap)) sqrtD = negi(sqrtD);
    3777             :   }
    3778             :   else
    3779         259 :     sqrtD = Zp_sqrtlift(D, ap, p, n);
    3780             :   /* sqrtD = ap (mod p) */
    3781         266 :   return gmul2n(gadd(ap, cvtop(sqrtD,p,n)), -1);
    3782             : }
    3783             : 
    3784             : /* W = msinit(N,k,...); phi = T_p/U_p - eigensymbol */
    3785             : GEN
    3786         483 : mstooms(GEN W, GEN phi)
    3787             : {
    3788         483 :   pari_sp av = avma;
    3789             :   GEN Wp, bin, Tp, c, alpha, ap, phi0, M;
    3790             :   long k, p, vden;
    3791             : 
    3792         483 :   checkmspadic(W);
    3793         483 :   if (typ(phi) != t_COL)
    3794             :   {
    3795         161 :     if (!is_Qevproj(phi)) pari_err_TYPE("mstooms",phi);
    3796         161 :     phi = gel(phi,1);
    3797         161 :     if (lg(phi) != 2) pari_err_TYPE("mstooms [dim_Q (eigenspace) > 1]",phi);
    3798         161 :     phi = gel(phi,1);
    3799             :   }
    3800             : 
    3801         483 :   Wp = mspadic_get_Wp(W);
    3802         483 :   Tp = mspadic_get_Tp(W);
    3803         483 :   bin = mspadic_get_bin(W);
    3804         483 :   k = msk_get_weight(Wp);
    3805         483 :   p = mspadic_get_p(W);
    3806         483 :   M = mspadic_get_M(W);
    3807             : 
    3808         483 :   phi = Q_remove_denom(phi, &c);
    3809         483 :   ap = eigenvalue(Tp, phi);
    3810         483 :   vden = c? Z_lvalrem(c, p, &c): 0;
    3811             : 
    3812         483 :   if (typ(M) == t_INT)
    3813             :   { /* p | N */
    3814             :     GEN c1;
    3815          84 :     alpha = ap;
    3816          84 :     alpha = ginv(alpha);
    3817          84 :     phi0 = mseval(Wp, phi, NULL);
    3818          84 :     phi0 = RgXC_to_moments(phi0, bin);
    3819          84 :     phi0 = Q_remove_denom(phi0, &c1);
    3820          84 :     if (c1) { vden += Z_lvalrem(c1, p, &c1); c = mul_denom(c,c1); }
    3821          84 :     if (umodiu(ap,p)) /* p \nmid a_p */
    3822          49 :       phi = oms_dim1(W, phi0, alpha, 0);
    3823             :     else
    3824             :     {
    3825          35 :       phi = oms_dim1(W, phi0, alpha, 1);
    3826          35 :       phi = Q_remove_denom(phi, &c1);
    3827          35 :       if (c1) { vden += Z_lvalrem(c1, p, &c1); c = mul_denom(c,c1); }
    3828             :     }
    3829             :   }
    3830             :   else
    3831             :   { /* p-stabilize */
    3832             :     GEN M1, M2, phi1, phi2, c1;
    3833         399 :     if (typ(M) != t_VEC || lg(M) != 3) pari_err_TYPE("mstooms",W);
    3834         399 :     M1 = gel(M,1);
    3835         399 :     M2 = gel(M,2);
    3836             : 
    3837         399 :     phi1 = RgM_RgC_mul(M1, phi);
    3838         399 :     phi2 = RgM_RgC_mul(M2, phi);
    3839         399 :     phi1 = mseval(Wp, phi1, NULL);
    3840         399 :     phi2 = mseval(Wp, phi2, NULL);
    3841             : 
    3842         399 :     phi1 = RgXC_to_moments(phi1, bin);
    3843         399 :     phi2 = RgXC_to_moments(phi2, bin);
    3844         399 :     phi = Q_remove_denom(mkvec2(phi1,phi2), &c1);
    3845         399 :     phi1 = gel(phi,1);
    3846         399 :     phi2 = gel(phi,2);
    3847         399 :     if (c1) { vden += Z_lvalrem(c1, p, &c1); c = mul_denom(c,c1); }
    3848             :     /* all polynomials multiplied by c p^vden */
    3849         399 :     if (umodiu(ap, p))
    3850             :     {
    3851         238 :       alpha = mspadic_unit_eigenvalue(ap, k, utoipos(p), mspadic_get_n(W));
    3852         238 :       alpha = ginv(alpha);
    3853         238 :       phi0 = gsub(phi1, gmul(lift_shallow(alpha),phi2));
    3854         238 :       phi = oms_dim1(W, phi0, alpha, 0);
    3855             :     }
    3856             :     else
    3857             :     { /* p | ap, alpha = [a_p, -1; p^(k-1), 0] */
    3858         161 :       long flag = mspadic_get_flag(W);
    3859         161 :       if (!flag || (signe(ap) && Z_lval(ap,p) < flag))
    3860           7 :         pari_err_TYPE("mstooms [v_p(ap) > mspadicinit flag]", phi);
    3861         154 :       alpha = mkmat22(ap,gen_m1, powuu(p, k-1),gen_0);
    3862         154 :       alpha = ginv(alpha);
    3863         154 :       phi = oms_dim2(W, mkvec2(phi1,phi2), gsqr(alpha), ap);
    3864         154 :       phi = Q_remove_denom(phi, &c1);
    3865         154 :       if (c1) { vden += Z_lvalrem(c1, p, &c1); c = mul_denom(c,c1); }
    3866             :     }
    3867             :   }
    3868         476 :   if (vden) c = mul_denom(c, powuu(p,vden));
    3869         476 :   if (p == 2) alpha = gsqr(alpha);
    3870         476 :   if (c) alpha = gdiv(alpha,c);
    3871         476 :   if (typ(alpha) == t_MAT)
    3872             :   { /* express in basis (omega,-p phi(omega)) */
    3873         154 :     gcoeff(alpha,2,1) = gdivgs(gcoeff(alpha,2,1), -p);
    3874         154 :     gcoeff(alpha,2,2) = gdivgs(gcoeff(alpha,2,2), -p);
    3875             :     /* at the end of mspadicint we shall multiply result by [1,0;0,-1/p]*alpha
    3876             :      * vden + k is the denominator of this matrix */
    3877             :   }
    3878             :   /* phi is integral-valued */
    3879         476 :   return gerepilecopy(av, mkcol3(phi, stoi(vden), alpha));
    3880             : }
    3881             : 
    3882             : /* HACK: the v[j] have different lengths */
    3883             : static GEN
    3884        2156 : FpVV_dotproduct(GEN v, GEN w, GEN p)
    3885             : {
    3886        2156 :   long j, l = lg(v);
    3887        2156 :   GEN T = cgetg(l, t_VEC);
    3888        2156 :   for (j = 1; j < l; j++) gel(T,j) = FpV_dotproduct(gel(v,j),w,p);
    3889        2156 :   return T;
    3890             : }
    3891             : 
    3892             : /* \int (-4z)^j given \int z^j */
    3893             : static GEN
    3894          98 : twistmoment_m4(GEN v)
    3895             : {
    3896             :   long i, l;
    3897          98 :   GEN w = cgetg_copy(v, &l);
    3898        2009 :   for (i = 1; i < l; i++)
    3899             :   {
    3900        1911 :     GEN c = gel(v,i);
    3901        1911 :     if (i > 1) c = gmul2n(c, (i-1)<<1);
    3902        1911 :     gel(w,i) = odd(i)? c: gneg(c);
    3903             :   }
    3904          98 :   return w;
    3905             : }
    3906             : /* \int (4z)^j given \int z^j */
    3907             : static GEN
    3908          98 : twistmoment_4(GEN v)
    3909             : {
    3910             :   long i, l;
    3911          98 :   GEN w = cgetg_copy(v, &l);
    3912        2009 :   for (i = 1; i < l; i++)
    3913             :   {
    3914        1911 :     GEN c = gel(v,i);
    3915        1911 :     if (i > 1) c = gmul2n(c, (i-1)<<1);
    3916        1911 :     gel(w,i) = c;
    3917             :   }
    3918          98 :   return w;
    3919             : }
    3920             : /* W an mspadic, phi eigensymbol, p \nmid D. Return C(x) mod FilM */
    3921             : GEN
    3922         483 : mspadicmoments(GEN W, GEN PHI, long D)
    3923             : {
    3924         483 :   pari_sp av = avma;
    3925         483 :   long na, ia, b, lphi, aD = labs(D), pp, p, k, n, vden;
    3926             :   GEN Wp, Dact, vL, v, C, pn, phi;
    3927             :   struct m_act S;
    3928             :   hashtable *H;
    3929             : 
    3930         483 :   checkmspadic(W);
    3931         483 :   Wp = mspadic_get_Wp(W);
    3932         483 :   p = mspadic_get_p(W);
    3933         483 :   k = mspadic_get_weight(W);
    3934         483 :   n = mspadic_get_n(W);
    3935         483 :   C = mspadic_get_C(W);
    3936         483 :   if (typ(PHI) != t_COL || lg(PHI) != 4 || typ(gel(PHI,1)) != t_VEC)
    3937         476 :     PHI = mstooms(W, PHI);
    3938         476 :   vden = itos( gel(PHI,2) );
    3939         476 :   phi = gel(PHI,1); lphi = lg(phi);
    3940         476 :   if (p == 2) { na = 2; pp = 4; }
    3941         420 :   else        { na = p-1; pp = p; }
    3942         476 :   pn = powuu(p, n + vden);
    3943             : 
    3944         476 :   S.p = p;
    3945         476 :   S.k = k;
    3946         476 :   S.q = pn;
    3947         476 :   S.dim = n+k-1;
    3948         476 :   S.act = &moments_act;
    3949         476 :   H = Gl2act_cache(ms_get_nbgen(Wp));
    3950         476 :   if (D == 1) Dact = NULL;
    3951             :   else
    3952             :   {
    3953          63 :     GEN gaD = utoi(aD), Dk = Fp_pows(stoi(D), 2-k, pn);
    3954          63 :     if (!sisfundamental(D)) pari_err_TYPE("mspadicmoments", stoi(D));
    3955          63 :     if (D % p == 0) pari_err_DOMAIN("mspadicmoments","p","|", stoi(D), utoi(p));
    3956          63 :     Dact = cgetg(aD, t_VEC);
    3957         532 :     for (b = 1; b < aD; b++)
    3958             :     {
    3959         469 :       GEN z = NULL;
    3960         469 :       long s = kross(D, b);
    3961         469 :       if (s)
    3962             :       {
    3963         462 :         pari_sp av2 = avma;
    3964             :         GEN d;
    3965         462 :         z = moments_act_i(&S, mkmat22(gaD,utoipos(b), gen_0,gaD));
    3966         462 :         d = s > 0? Dk: Fp_neg(Dk, pn);
    3967         924 :         z = equali1(d)? gerepilecopy(av2, z)
    3968         462 :                       : gerepileupto(av2, FpM_Fp_mul(z, d, pn));
    3969             :       }
    3970         469 :       gel(Dact,b) = z;
    3971             :     }
    3972             :   }
    3973         476 :   vL = cgetg(na+1,t_VEC);
    3974             :   /* first pass to precompute log(paths), preload matrices and allow GC later */
    3975        2464 :   for (ia = 1; ia <= na; ia++)
    3976             :   {
    3977             :     GEN path, La;
    3978        1988 :     long a = (p == 2 && ia == 2)? -1: ia;
    3979        1988 :     if (Dact)
    3980             :     { /* twist by D */
    3981         224 :       La = cgetg(aD, t_VEC);
    3982        1442 :       for (b = 1; b < aD; b++)
    3983             :       {
    3984        1218 :         GEN Actb = gel(Dact,b);
    3985        1218 :         if (!Actb) continue;
    3986             :         /* oo -> a/pp + b/|D|*/
    3987        1176 :         path = mkmat22(gen_1, addii(mulss(a, aD), muluu(pp, b)),
    3988             :                        gen_0, muluu(pp, aD));
    3989        1176 :         gel(La,b) = M2_logf(Wp,path,NULL);
    3990        1176 :         ZGl2QC_preload(&S, gel(La,b), H);
    3991             :       }
    3992             :     }
    3993             :     else
    3994             :     {
    3995        1764 :       path = mkmat22(gen_1,stoi(a), gen_0, utoipos(pp));
    3996        1764 :       La = M2_logf(Wp,path,NULL);
    3997        1764 :       ZGl2QC_preload(&S, La, H);
    3998             :     }
    3999        1988 :     gel(vL,ia) = La;
    4000             :   }
    4001         476 :   v = cgetg(na+1,t_VEC);
    4002             :   /* second pass, with GC */
    4003        2464 :   for (ia = 1; ia <= na; ia++)
    4004             :   {
    4005        1988 :     pari_sp av2 = avma;
    4006        1988 :     GEN vca, Ca = gel(C,ia), La = gel(vL,ia), va = cgetg(lphi, t_VEC);
    4007             :     long i;
    4008        1988 :     if (!Dact) vca = omseval_int(&S, phi, La, H);
    4009             :     else
    4010             :     { /* twist by D */
    4011         224 :       vca = cgetg(lphi,t_VEC);
    4012        1442 :       for (b = 1; b < aD; b++)
    4013             :       {
    4014        1218 :         GEN T, Actb = gel(Dact,b);
    4015        1218 :         if (!Actb) continue;
    4016        1176 :         T = omseval_int(&S, phi, gel(La,b), H);
    4017        2352 :         for (i = 1; i < lphi; i++)
    4018             :         {
    4019        1176 :           GEN z = FpM_FpC_mul(Actb, gel(T,i), pn);
    4020        1176 :           gel(vca,i) = b==1? z: ZC_add(gel(vca,i), z);
    4021             :         }
    4022             :       }
    4023             :     }
    4024        1988 :     if (p != 2)
    4025        1876 :     { for (i=1; i<lphi; i++) gel(va,i) = FpVV_dotproduct(Ca,gel(vca,i),pn); }
    4026         112 :     else if (ia == 1) /* \tilde{a} = 1 */
    4027          56 :     { for (i=1; i<lphi; i++) gel(va,i) = twistmoment_4(gel(vca,i)); }
    4028             :     else /* \tilde{a} = -1 */
    4029          56 :     { for (i=1; i<lphi; i++) gel(va,i) = twistmoment_m4(gel(vca,i)); }
    4030        1988 :     gel(v,ia) = gerepilecopy(av2, va);
    4031             :   }
    4032         476 :   return gerepilecopy(av, mkvec3(v, gel(PHI,3), mkvecsmall4(p,n+vden,n,D)));
    4033             : }
    4034             : static void
    4035        1918 : checkoms(GEN v)
    4036             : {
    4037        1918 :   if (typ(v) != t_VEC || lg(v) != 4 || typ(gel(v,1)) != t_VEC
    4038        1918 :       || typ(gel(v,3))!=t_VECSMALL)
    4039           0 :     pari_err_TYPE("checkoms [apply mspadicmoments]", v);
    4040        1918 : }
    4041             : static long
    4042        4284 : oms_get_p(GEN oms) { return gel(oms,3)[1]; }
    4043             : static long
    4044        4186 : oms_get_n(GEN oms) { return gel(oms,3)[2]; }
    4045             : static long
    4046        2464 : oms_get_n0(GEN oms) { return gel(oms,3)[3]; }
    4047             : static long
    4048        1918 : oms_get_D(GEN oms) { return gel(oms,3)[4]; }
    4049             : static int
    4050          98 : oms_is_supersingular(GEN oms) { GEN v = gel(oms,1); return lg(gel(v,1)) == 3; }
    4051             : 
    4052             : /* sum(j = 1, n, (-1)^(j+1)/j * x^j) */
    4053             : static GEN
    4054         784 : log1x(long n)
    4055             : {
    4056         784 :   long i, l = n+3;
    4057         784 :   GEN v = cgetg(l, t_POL);
    4058         784 :   v[1] = evalvarn(0)|evalsigne(1); gel(v,2) = gen_0;
    4059         784 :   for (i = 3; i < l; i++) gel(v,i) = ginv(stoi(odd(i)? i-2: 2-i));
    4060         784 :   return v;
    4061             : }
    4062             : 
    4063             : /* S = (1+x)^zk log(1+x)^logj (mod x^(n+1)) */
    4064             : static GEN
    4065        1820 : xlog1x(long n, long zk, long logj, long *pteich)
    4066             : {
    4067        1820 :   GEN S = logj? RgXn_powu_i(log1x(n), logj, n+1): NULL;
    4068        1820 :   if (zk)
    4069             :   {
    4070        1183 :     GEN L = deg1pol_shallow(gen_1, gen_1, 0); /* x+1 */
    4071        1183 :     *pteich += zk;
    4072        1183 :     if (zk < 0) { L = RgXn_inv(L,n+1); zk = -zk; }
    4073        1183 :     if (zk != 1) L = RgXn_powu_i(L, zk, n+1);
    4074        1183 :     S = S? RgXn_mul(S, L, n+1): L;
    4075             :   }
    4076        1820 :   return S;
    4077             : }
    4078             : 
    4079             : /* oms from mspadicmoments; integrate teichmuller^i * S(x) [S = NULL: 1]*/
    4080             : static GEN
    4081        2366 : mspadicint(GEN oms, long teichi, GEN S)
    4082             : {
    4083        2366 :   pari_sp av = avma;
    4084        2366 :   long p = oms_get_p(oms), n = oms_get_n(oms), n0 = oms_get_n0(oms);
    4085        2366 :   GEN vT = gel(oms,1), alpha = gel(oms,2), gp = utoipos(p);
    4086        2366 :   long loss = S? Z_lval(Q_denom(S), p): 0;
    4087        2366 :   long nfinal = minss(n-loss, n0);
    4088        2366 :   long i, la, l = lg(gel(vT,1));
    4089        2366 :   GEN res = cgetg(l, t_COL), teich = NULL;
    4090             : 
    4091        2366 :   if (S) S = RgX_to_RgC(S,lg(gmael(vT,1,1))-1);
    4092        2366 :   if (p == 2)
    4093             :   {
    4094         448 :     la = 3; /* corresponds to [1,-1] */
    4095         448 :     teichi &= 1;
    4096             :   }
    4097             :   else
    4098             :   {
    4099        1918 :     la = p; /* corresponds to [1,2,...,p-1] */
    4100        1918 :     teichi = smodss(teichi, p-1);
    4101        1918 :     if (teichi) teich = teichmullerinit(p, n);
    4102             :   }
    4103        5446 :   for (i=1; i<l; i++)
    4104             :   {
    4105        3080 :     pari_sp av2 = avma;
    4106        3080 :     GEN s = gen_0;
    4107             :     long ia;
    4108       14756 :     for (ia = 1; ia < la; ia++)
    4109             :     { /* Ta[j+1] correct mod p^n */
    4110       11676 :       GEN Ta = gmael(vT,ia,i), v = S? RgV_dotproduct(Ta, S): gel(Ta,1);
    4111       11676 :       if (teichi && ia != 1)
    4112             :       {
    4113        3843 :         if (p != 2)
    4114        3626 :           v = gmul(v, gel(teich, Fl_powu(ia,teichi,p)));
    4115             :         else
    4116         217 :           if (teichi) v = gneg(v);
    4117             :       }
    4118       11676 :       s = gadd(s, v);
    4119             :     }
    4120        3080 :     s = gadd(s, zeropadic(gp,nfinal));
    4121        3080 :     gel(res,i) = gerepileupto(av2, s);
    4122             :   }
    4123        2366 :   return gerepileupto(av, gmul(alpha, res));
    4124             : }
    4125             : /* integrate P = polynomial in log(x); vlog[j+1] = mspadicint(0,log(1+x)^j) */
    4126             : static GEN
    4127         539 : mspadicint_RgXlog(GEN P, GEN vlog)
    4128             : {
    4129         539 :   long i, d = degpol(P);
    4130         539 :   GEN s = gmul(gel(P,2), gel(vlog,1));
    4131         539 :   for (i = 1; i <= d; i++) s = gadd(s, gmul(gel(P,i+2), gel(vlog,i+1)));
    4132         539 :   return s;
    4133             : };
    4134             : 
    4135             : /* oms from mspadicmoments */
    4136             : GEN
    4137          98 : mspadicseries(GEN oms, long teichi)
    4138             : {
    4139          98 :   pari_sp av = avma;
    4140             :   GEN S, L, X, vlog, s, s2, u, logu, bin;
    4141             :   long j, p, m, n, step, stop;
    4142          98 :   checkoms(oms);
    4143          98 :   n = oms_get_n0(oms);
    4144          98 :   if (n < 1)
    4145             :   {
    4146           0 :     s = zeroser(0,0);
    4147           0 :     if (oms_is_supersingular(oms)) s = mkvec2(s,s);
    4148           0 :     return gerepilecopy(av, s);
    4149             :   }
    4150          98 :   p = oms_get_p(oms);
    4151          98 :   vlog = cgetg(n+1, t_VEC);
    4152          98 :   step = p == 2? 2: 1;
    4153          98 :   stop = 0;
    4154          98 :   S = NULL;
    4155          98 :   L = log1x(n);
    4156         644 :   for (j = 0; j < n; j++)
    4157             :   {
    4158         616 :     if (j) stop += step + u_lval(j,p); /* = step*j + v_p(j!) */
    4159         616 :     if (stop >= n) break;
    4160             :     /* S = log(1+x)^j */
    4161         546 :     gel(vlog,j+1) = mspadicint(oms,teichi,S);
    4162         546 :     S = S? RgXn_mul(S, L, n+1): L;
    4163             :   }
    4164          98 :   m = j;
    4165          98 :   u = utoipos(p == 2? 5: 1+p);
    4166          98 :   logu = glog(cvtop(u, utoipos(p), 4*m), 0);
    4167          98 :   X = gdiv(pol_x(0), logu);
    4168          98 :   s = cgetg(m+1, t_VEC);
    4169          98 :   s2 = oms_is_supersingular(oms)? cgetg(m+1, t_VEC): NULL;
    4170          98 :   bin = pol_1(0);
    4171         539 :   for (j = 0; j < m; j++)
    4172             :   { /* bin = binomial(x/log(1+p+O(p^(4*n))), j) mod x^m */
    4173         539 :     GEN a, v = mspadicint_RgXlog(bin, vlog);
    4174         539 :     int done = 1;
    4175         539 :     gel(s,j+1) = a = gel(v,1);
    4176         539 :     if (!gequal0(a) || valp(a) > 0) done = 0; else setlg(s,j+1);
    4177         539 :     if (s2)
    4178             :     {
    4179         119 :       gel(s2,j+1) = a = gel(v,2);
    4180         119 :       if (!gequal0(a) || valp(a) > 0) done = 0; else setlg(s2,j+1);
    4181             :     }
    4182         539 :     if (done || j == m-1) break;
    4183         441 :     bin = RgXn_mul(bin, gdivgs(gsubgs(X, j), j+1), m);
    4184             :   }
    4185          98 :   s = RgV_to_ser(s,0,lg(s)+1);
    4186          98 :   if (s2) { s2 = RgV_to_ser(s2,0,lg(s2)+1); s = mkvec2(s, s2); }
    4187          98 :   if (kross(oms_get_D(oms), p) >= 0) return gerepilecopy(av, s);
    4188           7 :   return gerepileupto(av, gneg(s));
    4189             : }
    4190             : void
    4191        1911 : mspadic_parse_chi(GEN s, GEN *s1, GEN *s2)
    4192             : {
    4193        1911 :   if (!s) *s1 = *s2 = gen_0;
    4194        1778 :   else switch(typ(s))
    4195             :   {
    4196        1274 :     case t_INT: *s1 = *s2 = s; break;
    4197             :     case t_VEC:
    4198         504 :       if (lg(s) == 3)
    4199             :       {
    4200         504 :         *s1 = gel(s,1);
    4201         504 :         *s2 = gel(s,2);
    4202         504 :         if (typ(*s1) == t_INT && typ(*s2) == t_INT) break;
    4203             :       }
    4204           0 :     default: pari_err_TYPE("mspadicL",s);
    4205           0 :              *s1 = *s2 = NULL;
    4206             :   }
    4207        1911 : }
    4208             : /* oms from mspadicmoments
    4209             :  * r-th derivative of L(f,chi^s,psi) in direction <chi>
    4210             :    - s \in Z_p \times \Z/(p-1)\Z, s-> chi^s=<\chi>^s_1 omega^s_2)
    4211             :    - Z -> Z_p \times \Z/(p-1)\Z par s-> (s, s mod p-1).
    4212             :  */
    4213             : GEN
    4214        1820 : mspadicL(GEN oms, GEN s, long r)
    4215             : {
    4216        1820 :   pari_sp av = avma;
    4217             :   GEN s1, s2, z, S;
    4218             :   long p, n, teich;
    4219        1820 :   checkoms(oms);
    4220        1820 :   p = oms_get_p(oms);
    4221        1820 :   n = oms_get_n(oms);
    4222        1820 :   mspadic_parse_chi(s, &s1,&s2);
    4223        1820 :   teich = umodiu(subii(s2,s1), p==2? 2: p-1);
    4224        1820 :   S = xlog1x(n, itos(s1), r, &teich);
    4225        1820 :   z = mspadicint(oms, teich, S);
    4226        1820 :   if (lg(z) == 2) z = gel(z,1);
    4227        1820 :   if (kross(oms_get_D(oms), p) < 0) z = gneg(z);
    4228        1820 :   return gerepilecopy(av, z);
    4229             : }
    4230             : 
    4231             : /****************************************************************************/
    4232             : 
    4233             : struct siegel
    4234             : {
    4235             :   GEN V, Ast;
    4236             :   long N; /* level */
    4237             :   long oo; /* index of the [oo,0] path */
    4238             :   long k1, k2; /* two distinguished indices */
    4239             :   long n; /* #W, W = initial segment [in siegelstepC] already normalized */
    4240             : };
    4241             : 
    4242             : static void
    4243         357 : siegel_init(struct siegel *C, GEN M)
    4244             : {
    4245             :   GEN CPI, CP, MM, V, W, Ast;
    4246         357 :   GEN m = gel(M,11), M2 = gel(M,2), S = msN_get_section(M);
    4247         357 :   GEN E2fromE1 = msN_get_E2fromE1(M);
    4248         357 :   long m0 = lg(M2)-1;
    4249         357 :   GEN E2  = vecslice(M2, m[1]+1, m[2]);/* E2 */
    4250         357 :   GEN E1T = vecslice(M2, m[3]+1, m0); /* E1,T2,T31 */
    4251         357 :   GEN L = shallowconcat(E1T, E2);
    4252         357 :   long i, l = lg(L), n = lg(E1T)-1, lE = lg(E2);
    4253             : 
    4254         357 :   Ast = cgetg(l, t_VECSMALL);
    4255        6069 :   for (i = 1; i < lE; ++i)
    4256             :   {
    4257        5712 :     long j = E2fromE1_c(gel(E2fromE1,i));
    4258        5712 :     Ast[n+i] = j;
    4259        5712 :     Ast[j] = n+i;
    4260             :   }
    4261         357 :   for (; i<=n; ++i) Ast[i] = i;
    4262         357 :   MM = cgetg (l,t_VEC);
    4263             : 
    4264       12089 :   for (i = 1; i < l; i++)
    4265             :   {
    4266       11732 :     GEN c = gel(S, L[i]);
    4267       11732 :     long c12, c22, c21 = ucoeff(c,2,1);
    4268       11732 :     if (!c21) { gel(MM,i) = gen_0; continue; }
    4269       11375 :     c22 = ucoeff(c,2,2);
    4270       11375 :     if (!c22) { gel(MM,i) = gen_m1; continue; }
    4271       11018 :     c12 = ucoeff(c,1,2);
    4272       11018 :     gel(MM,i) = gdivgs(stoi(c12), c22); /* right extremity > 0 */
    4273             :   }
    4274         357 :   CP = indexsort(MM);
    4275         357 :   CPI = cgetg(l, t_VECSMALL);
    4276         357 :   V = cgetg(l, t_VEC);
    4277         357 :   W = cgetg(l, t_VECSMALL);
    4278       12089 :   for (i = 1; i < l; ++i)
    4279             :   {
    4280       11732 :     gel(V,i) = zm_to_ZM(gel(S, L[CP[i]]));
    4281       11732 :     CPI[CP[i]] = i;
    4282             :   }
    4283         357 :   for (i = 1; i < l; ++i) W[CPI[i]] = CPI[Ast[i]];
    4284         357 :   C->V = V;
    4285         357 :   C->Ast = W;
    4286         357 :   C->n = 0;
    4287         357 :   C->oo = 2;
    4288         357 :   C->N = ms_get_N(M);
    4289         357 : }
    4290             : 
    4291             : static double
    4292           0 : ZMV_size(GEN v)
    4293             : {
    4294           0 :   long i, l = lg(v);
    4295           0 :   GEN z = cgetg(l, t_VECSMALL);
    4296           0 :   for (i = 1; i < l; i++) z[i] = gexpo(gel(v,i));
    4297           0 :   return ((double)zv_sum(z)) / (4*(l-1));
    4298             : }
    4299             : 
    4300             : /* apply permutation perm to struct S. Don't follow k1,k2 */
    4301             : static void
    4302        5558 : siegel_perm0(struct siegel *S, GEN perm)
    4303             : {
    4304        5558 :   long i, l = lg(S->V);
    4305        5558 :   GEN V2 = cgetg(l, t_VEC), Ast2 = cgetg(l, t_VECSMALL);
    4306        5558 :   GEN V = S->V, Ast = S->Ast;
    4307             : 
    4308        5558 :   for (i = 1; i < l; i++) gel(V2,perm[i]) = gel(V,i);
    4309        5558 :   for (i = 1; i < l; i++) Ast2[perm[i]] = perm[Ast[i]];
    4310        5558 :   S->oo = perm[S->oo];
    4311        5558 :   S->Ast = Ast2;
    4312        5558 :   S->V = V2;
    4313        5558 : }
    4314             : /* apply permutation perm to full struct S */
    4315             : static void
    4316        5194 : siegel_perm(struct siegel *S, GEN perm)
    4317             : {
    4318        5194 :   siegel_perm0(S, perm);
    4319        5194 :   S->k1 = perm[S->k1];
    4320        5194 :   S->k2 = perm[S->k2];
    4321        5194 : }
    4322             : /* cyclic permutation of lg = l-1 moving a -> 1, a+1 -> 2, etc.  */
    4323             : static GEN
    4324        2884 : rotate_perm(long l, long a)
    4325             : {
    4326        2884 :   GEN p = cgetg(l, t_VECSMALL);
    4327        2884 :   long i, j = 1;
    4328        2884 :   for (i = a; i < l; i++) p[i] = j++;
    4329        2884 :   for (i = 1; i < a; i++) p[i] = j++;
    4330        2884 :   return p;
    4331             : }
    4332             : 
    4333             : /* a1 < c1 <= a2 < c2*/
    4334             : static GEN
    4335        2520 : basic_op_perm(long l, long a1, long a2, long c1, long c2)
    4336             : {
    4337        2520 :   GEN p = cgetg(l, t_VECSMALL);
    4338        2520 :   long i, j = 1;
    4339        2520 :   p[a1] = j++;
    4340        2520 :   for (i = c1; i < a2; i++)   p[i] = j++;
    4341        2520 :   for (i = a1+1; i < c1; i++) p[i] = j++;
    4342        2520 :   p[a2] = j++;
    4343        2520 :   for (i = c2; i < l; i++)    p[i] = j++;
    4344        2520 :   for (i = 1; i < a1; i++)    p[i] = j++;
    4345        2520 :   for (i = a2+1; i < c2; i++) p[i] = j++;
    4346        2520 :   return p;
    4347             : }
    4348             : static GEN
    4349         154 : basic_op_perm_elliptic(long l, long a1)
    4350             : {
    4351         154 :   GEN p = cgetg(l, t_VECSMALL);
    4352         154 :   long i, j = 1;
    4353         154 :   p[a1] = j++;
    4354         154 :   for (i = 1; i < a1; i++)   p[i] = j++;
    4355         154 :   for (i = a1+1; i < l; i++) p[i] = j++;
    4356         154 :   return p;
    4357             : }
    4358             : static GEN
    4359       14406 : ZM2_rev(GEN T) { return mkmat2(gel(T,2), ZC_neg(gel(T,1))); }
    4360             : 
    4361             : /* In place, V = vector of consecutive paths, between x <= y.
    4362             :  * V[x..y-1] <- g*V[x..y-1] */
    4363             : static void
    4364        5733 : path_vec_mul(GEN V, long x, long y, GEN g)
    4365             : {
    4366             :   long j;
    4367             :   GEN M;
    4368       11466 :   if (x == y) return;
    4369        3360 :   M = gel(V,x); gel(V,x) = ZM_mul(g,M);
    4370       37709 :   for (j = x+1; j < y; j++) /* V[j] <- g*V[j], optimized */
    4371             :   {
    4372       34349 :     GEN Mnext = gel(V,j); /* Mnext[,1] = M[,2] */
    4373       34349 :     GEN gM = gel(V,j-1), u = gel(gM,2);
    4374       34349 :     if (!ZV_equal(gel(M,2), gel(Mnext,1))) u = ZC_neg(u);
    4375       34349 :     gel(V,j) = mkmat2(u, ZM_ZC_mul(g,gel(Mnext,2)));
    4376       34349 :     M = Mnext;
    4377             :   }
    4378             : }
    4379             : 
    4380        4830 : static long prev(GEN V, long i) { return (i == 1)? lg(V)-1: i-1; }
    4381        4830 : static long next(GEN V, long i) { return (i == lg(V)-1)? 1: i+1; }
    4382             : static GEN
    4383       19600 : ZM_det2(GEN u, GEN v)
    4384             : {
    4385       19600 :   GEN a = gel(u,1), c = gel(u,2);
    4386       19600 :   GEN b = gel(v,1), d = gel(v,2); return subii(mulii(a,d), mulii(b,c));
    4387             : }
    4388             : static GEN
    4389       14406 : ZM2_det(GEN T) { return ZM_det2(gel(T,1),gel(T,2)); }
    4390             : static void
    4391        4466 : fill1(GEN V, long a)
    4392             : {
    4393        4466 :   long p = prev(V,a), n = next(V,a);
    4394        4466 :   GEN u = gmael(V,p,2), v = gmael(V,n,1);
    4395        4466 :   if (signe(ZM_det2(u,v)) < 0) v = ZC_neg(v);
    4396        4466 :   gel(V,a) = mkmat2(u, v);
    4397        4466 : }
    4398             : /* a1 < a2 */
    4399             : static void
    4400        2520 : fill2(GEN V, long a1, long a2)
    4401             : {
    4402        2520 :   if (a2 != a1+1) { fill1(V,a1); fill1(V,a2); } /* non adjacent, reconnect */
    4403             :   else /* parabolic */
    4404             :   {
    4405         364 :     long p = prev(V,a1), n = next(V,a2);
    4406         364 :     GEN u, v, C = gmael(V,a1,2), mC = ZC_neg(C); /* = \pm V[a2][1] */
    4407         364 :     u = gmael(V,p,2); v = (signe(ZM_det2(u,C)) < 0)? mC: C;
    4408         364 :     gel(V,a1) = mkmat2(u,v);
    4409         364 :     v = gmael(V,n,1); u = (signe(ZM_det2(C,v)) < 0)? mC: C;
    4410         364 :     gel(V,a2) = mkmat2(u,v);
    4411             :   }
    4412        2520 : }
    4413             : 
    4414             : /* DU = det(U), return g = T*U^(-1) or NULL if not in Gamma0(N); if N = 0,
    4415             :  * only test whether g is integral */
    4416             : static GEN
    4417       14693 : ZM2_div(GEN T, GEN U, GEN DU, long N)
    4418             : {
    4419       14693 :   GEN a=gcoeff(U,1,1), b=gcoeff(U,1,2), c=gcoeff(U,2,1), d=gcoeff(U,2,2);
    4420       14693 :   GEN e=gcoeff(T,1,1), f=gcoeff(T,1,2), g=gcoeff(T,2,1), h=gcoeff(T,2,2);
    4421             :   GEN A, B, C, D, r;
    4422             : 
    4423       14693 :   C = dvmdii(subii(mulii(d,g), mulii(c,h)), DU, &r);
    4424       14693 :   if (r != gen_0 || (N && smodis(C,N))) return NULL;
    4425       14406 :   A = dvmdii(subii(mulii(d,e), mulii(c,f)), DU, &r);
    4426       14406 :   if (r != gen_0) return NULL;
    4427       14406 :   B = dvmdii(subii(mulii(a,f), mulii(b,e)), DU, &r);
    4428       14406 :   if (r != gen_0) return NULL;
    4429       14406 :   D = dvmdii(subii(mulii(a,h), mulii(g,b)), DU, &r);
    4430       14406 :   if (r != gen_0) return NULL;
    4431       14406 :   return mkmat22(A,B,C,D);
    4432             : }
    4433             : 
    4434             : static GEN
    4435       14406 : get_g(struct siegel *S, long a1)
    4436             : {
    4437       14406 :   long a2 = S->Ast[a1];
    4438       14406 :   GEN a = gel(S->V,a1), ar = ZM2_rev(gel(S->V,a2)), Dar = ZM2_det(ar);
    4439       14406 :   GEN g = ZM2_div(a, ar, Dar, S->N);
    4440       14406 :   if (!g)
    4441             :   {
    4442         287 :     GEN tau = mkmat22(gen_0,gen_m1, gen_1,gen_m1); /*[0,-1;1,-1]*/
    4443         287 :     g = ZM2_div(ZM_mul(ar, tau), ar, Dar, 0);
    4444             :   }
    4445       14406 :   return g;
    4446             : }
    4447             : /* input V = (X1 a X2 | X3 a^* X4) + Ast
    4448             :  * a1 = index of a
    4449             :  * a2 = index of a^*, inferred from a1. We must have a != a^*
    4450             :  * c1 = first cut [ index of first path in X3 ]
    4451             :  * c2 = second cut [ either in X4 or X1, index of first path ]
    4452             :  * Assume a < a^* (cf Paranoia below): c1 or c2 must be in
    4453             :  *    ]a,a^*], and the other in the "complement" ]a^*,a] */
    4454             : static void
    4455        2520 : basic_op(struct siegel *S, long a1, long c1, long c2)
    4456             : {
    4457        2520 :   long l = lg(S->V), a2 = S->Ast[a1];
    4458             :   GEN g;
    4459             : 
    4460        2520 :   if (a1 == a2)
    4461             :   { /* a = a^* */
    4462           0 :     g = get_g(S, a1);
    4463           0 :     path_vec_mul(S->V, a1+1, l, g);
    4464           0 :     siegel_perm(S, basic_op_perm_elliptic(l, a1));
    4465             :     /* fill the hole left at a1, reconnect the path */
    4466        2520 :     fill1(S->V, a1); return;
    4467             :   }
    4468             : 
    4469             :   /* Paranoia: (a,a^*) conjugate, call 'a' the first one */
    4470        2520 :   if (a2 < a1) lswap(a1,a2);
    4471             :   /* Now a1 < a2 */
    4472        2520 :   if (c1 <= a1 || c1 > a2) lswap(c1,c2); /* ensure a1 < c1 <= a2 */
    4473        2520 :   if (c2 < a1)
    4474             :   { /* if cut c2 is in X1 = X11|X12, rotate to obtain
    4475             :        (a X2 | X3 a^* X4 X11|X12): then a1 = 1 */
    4476        2520 :     GEN p = rotate_perm(l, a1);
    4477        2520 :     siegel_perm(S, p);
    4478        2520 :     a1 = 1; /* = p[a1] */
    4479        2520 :     a2 = S->Ast[1]; /* > a1 */
    4480        2520 :     c1 = p[c1];
    4481        2520 :     c2 = p[c2];
    4482             :   }
    4483             :   /* Now a1 < c1 <= a2 < c2; a != a^* */
    4484        2520 :   g = get_g(S, a1);
    4485        2520 :   if (S->oo >= c1 && S->oo < c2) /* W inside [c1..c2[ */
    4486         539 :   { /* c2 -> c1 excluding a1 */
    4487         539 :     GEN gi = SL2_inv(g); /* g a^* = a; gi a = a^* */
    4488         539 :     path_vec_mul(S->V, 1, a1, gi);
    4489         539 :     path_vec_mul(S->V, a1+1, c1, gi);
    4490         539 :     path_vec_mul(S->V, c2, l, gi);
    4491             :   }
    4492             :   else
    4493             :   { /* c1 -> c2 excluding a2 */
    4494        1981 :     path_vec_mul(S->V, c1, a2, g);
    4495        1981 :     path_vec_mul(S->V, a2+1, c2, g);
    4496             :   }
    4497        2520 :   siegel_perm(S, basic_op_perm(l, a1,a2, c1,c2));
    4498             :   /* fill the holes left at a1,a2, reconnect the path */
    4499        2520 :   fill2(S->V, a1, a2);
    4500             : }
    4501             : /* a = a^* (elliptic case) */
    4502             : static void
    4503         154 : basic_op_elliptic(struct siegel *S, long a1)
    4504             : {
    4505         154 :   long l = lg(S->V);
    4506         154 :   GEN g = get_g(S, a1);
    4507         154 :   path_vec_mul(S->V, a1+1, l, g);
    4508         154 :   siegel_perm(S, basic_op_perm_elliptic(l, a1));
    4509             :   /* fill the hole left at a1 (now at 1), reconnect the path */
    4510         154 :   fill1(S->V, 1);
    4511         154 : }
    4512             : 
    4513             : /* input V = W X a b Y a^* Z b^* T, W already normalized
    4514             :  * X = [n+1, k1-1], Y = [k2+1, Ast[k1]-1],
    4515             :  * Z = [Ast[k1]+1, Ast[k2]-1], T = [Ast[k2]+1, oo].
    4516             :  * Assume that X doesn't start by c c^* or a b a^* b^*. */
    4517             : static void
    4518        1057 : siegelstep(struct siegel *S)
    4519             : {
    4520        1057 :   if (S->Ast[S->k1] == S->k1)
    4521             :   {
    4522         154 :     basic_op_elliptic(S, S->k1);
    4523         154 :     S->n++;
    4524             :   }
    4525         903 :   else if (S->Ast[S->k1] == S->k1+1)
    4526             :   {
    4527         364 :     basic_op(S, S->k1, S->Ast[S->k1], 1); /* 1: W starts there */
    4528         364 :     S->n += 2;
    4529             :   }
    4530             :   else
    4531             :   {
    4532         539 :     basic_op(S, S->k2, S->Ast[S->k1], 1); /* 1: W starts there */
    4533         539 :     basic_op(S, S->k1, S->k2, S->Ast[S->k2]);
    4534         539 :     basic_op(S, S->Ast[S->k2], S->k2, S->Ast[S->k1]);
    4535         539 :     basic_op(S, S->k1, S->Ast[S->k1], S->Ast[S->k2]);
    4536         539 :     S->n += 4;
    4537             :   }
    4538        1057 : }
    4539             : 
    4540             : /* normalize hyperbolic polygon */
    4541             : static void
    4542         301 : mssiegel(struct siegel *S)
    4543             : {
    4544         301 :   pari_sp av = avma;
    4545             :   long k, t, nv;
    4546             : #ifdef COUNT
    4547             :   long countset[16];
    4548             :   for (k = 0; k < 16; k++) countset[k] = 0;
    4549             : #endif
    4550             : 
    4551         301 :   nv = lg(S->V)-1;
    4552         301 :   if (DEBUGLEVEL>1) err_printf("nv = %ld, expo = %.2f\n", nv,ZMV_size(S->V));
    4553         301 :   t = 0;
    4554        2506 :   while (S->n < nv)
    4555             :   {
    4556        1904 :     if (S->Ast[S->n+1] == S->n+1) { S->n++; continue; }
    4557        1778 :     if (S->Ast[S->n+1] == S->n+2) { S->n += 2; continue; }
    4558        1134 :     if (S->Ast[S->n+1] == S->n+3 && S->Ast[S->n+2] == S->n+4) { S->n += 4; continue; }
    4559        1057 :     k = nv;
    4560        2184 :     while (k > S->n)
    4561             :     {
    4562        1127 :       if (S->Ast[k] == k) { k--; continue; }
    4563        1099 :       if (S->Ast[k] == k-1) { k -= 2; continue; }
    4564        1057 :       if (S->Ast[k] == k-2 && S->Ast[k-1] == k-3) { k -= 4; continue; }
    4565        1057 :       break;
    4566             :     }
    4567        1057 :     if (k != nv) { siegel_perm0(S, rotate_perm(nv+1, k+1)); S->n += nv-k; }
    4568             : 
    4569        6223 :     for (k = S->n+1; k <= nv; k++)
    4570        6223 :       if (S->Ast[k] <= k) { t = S->Ast[k]; break; }
    4571        1057 :     S->k1 = t;
    4572        1057 :     S->k2 = t+1;
    4573             : #ifdef COUNT
    4574             :     countset[ ((S->k1-1 == S->n)
    4575             :               | ((S->k2 == S->Ast[S->k1]-1) << 1)
    4576             :               | ((S->Ast[S->k1] == S->Ast[S->k2]-1) << 2)
    4577             :               | ((S->Ast[S->k2] == nv) << 3)) ]++;
    4578             : #endif
    4579        1057 :     siegelstep(S);
    4580        1057 :     if (gc_needed(av,2))
    4581             :     {
    4582           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"mspolygon, n = %ld",S->n);
    4583           0 :       gerepileall(av, 2, &S->V, &S->Ast);
    4584             :     }
    4585             :   }
    4586         301 :   if (DEBUGLEVEL>1) err_printf("expo = %.2f\n", ZMV_size(S->V));
    4587             : #ifdef COUNT
    4588             :   for (k = 0; k < 16; k++)
    4589             :     err_printf("%3ld: %6ld\n", k, countset[k]);
    4590             : #endif
    4591         301 : }
    4592             : GEN
    4593         378 : mspolygon(GEN M, long flag)
    4594             : {
    4595         378 :   pari_sp av = avma;
    4596             :   struct siegel S;
    4597             :   long i, l;
    4598             :   GEN G, msN;
    4599         378 :   if (typ(M) == t_INT)
    4600             :   {
    4601         301 :     long N = itos(M);
    4602         301 :     if (N <= 0) pari_err_DOMAIN("msinit","N", "<=", gen_0,M);
    4603         301 :     msN = msinit_N(N);
    4604             :   }
    4605          77 :   else { checkms(M); msN = get_msN(M); }
    4606         378 :   if (flag < 0 || flag > 1) pari_err_FLAG("mspolygon");
    4607         378 :   if (ms_get_N(msN) == 1)
    4608             :   {
    4609          21 :     GEN S = mkS(), V = mkvec2(matid(2), S);
    4610          21 :     return gerepilecopy(av, mkvec3(V, mkvecsmall2(1,2), mkvec2(S, mkTAU())));
    4611             :   }
    4612         357 :   siegel_init(&S, msN);
    4613         357 :   l = lg(S.V);
    4614         357 :   if (flag)
    4615             :   {
    4616         301 :     long oo2 = 0;
    4617         301 :     mssiegel(&S);
    4618        3451 :     for (i = 1; i < l; i++)
    4619             :     {
    4620        3451 :       GEN c = gel(S.V, i);
    4621        3451 :       GEN c22 = gcoeff(c,2,2); if (!signe(c22)) { oo2 = i; break; }
    4622             :     }
    4623         301 :     if (!oo2) pari_err_BUG("mspolygon");
    4624         301 :     siegel_perm0(&S, rotate_perm(l, oo2));
    4625             :   }
    4626         357 :   G = cgetg(l, t_VEC);
    4627         357 :   for (i = 1; i < l; i++) gel(G,i) = get_g(&S, i);
    4628         357 :   return gerepilecopy(av, mkvec3(S.V, S.Ast, G));
    4629             : }
    4630             : 
    4631             : #if 0
    4632             : static int
    4633             : iselliptic(GEN Ast, long i) { return i == Ast[i]; }
    4634             : static int
    4635             : isparabolic(GEN Ast, long i)
    4636             : { long i2 = Ast[i]; return (i2 == i+1 || i2 == i-1); }
    4637             : #endif
    4638             : 
    4639             : /* M from msinit, F QM maximal rank */
    4640             : GEN
    4641          56 : mslattice(GEN M, GEN F)
    4642             : {
    4643          56 :   pari_sp av = avma;
    4644             :   long i, ivB, j, k, l, lF;
    4645             :   GEN D, U, G, A, vB, m, d;
    4646             : 
    4647          56 :   checkms(M);
    4648          56 :   if (!F) F = gel(mscuspidal(M, 0), 1);
    4649             :   else
    4650             :   {
    4651          28 :     if (is_Qevproj(F)) F = gel(F,1);
    4652          28 :     if (typ(F) != t_MAT) pari_err_TYPE("mslattice",F);
    4653             :   }
    4654          56 :   lF = lg(F); if (lF == 1) return cgetg(1, t_MAT);
    4655          56 :   D = mspolygon(M,0);
    4656          56 :   k = msk_get_weight(M);
    4657          56 :   F = vec_Q_primpart(F);
    4658          56 :   if (typ(F)!=t_MAT || !RgM_is_ZM(F)) pari_err_TYPE("mslattice",F);
    4659          56 :   G = gel(D,3); l = lg(G);
    4660          56 :   A = gel(D,2);
    4661          56 :   vB = cgetg(l, t_COL);
    4662          56 :   m = mkmat2(mkcol2(gen_0,gen_1), NULL);
    4663        6860 :   for (i = ivB = 1; i < l; i++)
    4664             :   {
    4665        6804 :     GEN B, vb, g = gel(G,i);
    4666        6804 :     if (A[i] < i) continue;
    4667        3409 :     gel(m,2) = SL2_inv2(g);
    4668        3409 :     vb = mseval(M, F, m);
    4669        3409 :     if (k == 2) B = vb;
    4670             :     else
    4671             :     {
    4672             :       long lB;
    4673         147 :       B = RgXV_to_RgM(vb, k-1);
    4674             :       /* add coboundaries */
    4675         147 :       B = shallowconcat(B, RgM_Rg_sub(RgX_act_Gl2Q(g, k), gen_1));
    4676             :       /* beware: the basis for RgX_act_Gl2Q is (X^(k-2),...,Y^(k-2)) */
    4677         147 :       lB = lg(B);
    4678         147 :       for (j = 1; j < lB; j++) gel(B,j) = vecreverse(gel(B,j));
    4679             :     }
    4680        3409 :     gel(vB, ivB++) = B;
    4681             :   }
    4682          56 :   setlg(vB, ivB);
    4683          56 :   vB = shallowmatconcat(vB);
    4684          56 :   if (ZM_equal0(vB)) return gerepilecopy(av, F);
    4685             : 
    4686          56 :   (void)QM_ImQ_hnfall(vB, &U, 0);
    4687          56 :   if (k > 2) U = rowslice(U, 1, lgcols(U)-k); /* remove coboundary part */
    4688          56 :   U = Q_remove_denom(U, &d);
    4689          56 :   F = ZM_hnf(ZM_mul(F, U));
    4690          56 :   if (d) F = RgM_Rg_div(F, d);
    4691          56 :   return gerepileupto(av, F);
    4692             : }
    4693             : 
    4694             : /**** Petersson scalar product ****/
    4695             : 
    4696             : /* oo -> g^(-1) oo */
    4697             : static GEN
    4698        6181 : cocycle(GEN g)
    4699        6181 : { return mkmat22(gen_1, gcoeff(g,2,2), gen_0, negi(gcoeff(g,2,1))); }
    4700             : 
    4701             : /* C = vecbinome(k-2) */
    4702             : static GEN
    4703       17983 : bil(GEN P, GEN Q, GEN C)
    4704             : {
    4705       17983 :   long i, n = lg(C)-2; /* k - 2 */
    4706             :   GEN s;
    4707       17983 :   if (!n) return gmul(P,Q);
    4708       17962 :   if (typ(P) != t_POL) P = scalarpol_shallow(P,0);
    4709       17962 :   if (typ(Q) != t_POL) Q = scalarpol_shallow(Q,0);
    4710       17962 :   s = gen_0;
    4711       71960 :   for (i = 0; i <= n; i++)
    4712             :   {
    4713       53998 :     GEN t = gdiv(gmul(RgX_coeff(P,i), RgX_coeff(Q, n-i)), gel(C,i+1));
    4714       53998 :     s = odd(i)? gsub(s, t): gadd(s, t);
    4715             :   }
    4716       17962 :   return s;
    4717             : }
    4718             : 
    4719             : static void
    4720        1351 : mspetersson_i(GEN W, GEN F, GEN G, GEN *pvf, GEN *pvg, GEN *pC)
    4721             : {
    4722        1351 :   GEN WN = get_msN(W), annT2, annT31, section, c, vf, vg;
    4723             :   long i, n1, n2, n3;
    4724             : 
    4725        1351 :   annT2 = msN_get_annT2(WN);
    4726        1351 :   annT31 = msN_get_annT31(WN);
    4727        1351 :   section = msN_get_section(WN);
    4728             : 
    4729        1351 :   if (ms_get_N(WN) == 1)
    4730             :   {
    4731           7 :     vf = cgetg(3, t_VEC);
    4732           7 :     vg = cgetg(3, t_VEC);
    4733           7 :     gel(vf,1) = mseval(W, F, gel(section,1));
    4734           7 :     gel(vf,2) = gneg(gel(vf,1));
    4735           7 :     n1 = 0;
    4736             :   }
    4737             :   else
    4738             :   {
    4739        1344 :     GEN singlerel = msN_get_singlerel(WN);
    4740        1344 :     GEN gen = msN_get_genindex(WN);
    4741        1344 :     long l = lg(gen);
    4742        1344 :     vf = cgetg(l, t_VEC);
    4743        1344 :     vg = cgetg(l, t_VEC); /* generators of Delta ordered as E1,T2,T31 */
    4744        1344 :     for (i = 1; i < l; i++) gel(vf, i) = mseval(W, F, gel(section,gen[i]));
    4745        1344 :     n1 = ms_get_nbE1(WN); /* E1 */
    4746        7420 :     for (i = 1; i <= n1; i++)
    4747             :     {
    4748        6076 :       c = cocycle(gcoeff(gel(singlerel,i),2,1));
    4749        6076 :       gel(vg, i) = mseval(W, G, c);
    4750             :     }
    4751             :   }
    4752        1351 :   n2 = lg(annT2)-1; /* T2 */
    4753        1386 :   for (i = 1; i <= n2; i++)
    4754             :   {
    4755          35 :     c = cocycle(gcoeff(gel(annT2,i), 2,1));
    4756          35 :     gel(vg, i+n1) = gmul2n(mseval(W, G, c), -1);
    4757             :   }
    4758        1351 :   n3 = lg(annT31)-1; /* T31 */
    4759        1386 :   for (i = 1; i <= n3; i++)
    4760             :   {
    4761             :     GEN f;
    4762          35 :     c = cocycle(gcoeff(gel(annT31,i), 2,1));
    4763          35 :     f = mseval(W, G, c);
    4764          35 :     c = cocycle(gcoeff(gel(annT31,i), 3,1));
    4765          35 :     gel(vg, i+n1+n2) = gdivgs(gadd(f, mseval(W, G, c)), 3);
    4766             :   }
    4767        1351 :   *pC = vecbinome(msk_get_weight(W) - 2);
    4768        1351 :   *pvf = vf;
    4769        1351 :   *pvg = vg;
    4770        1351 : }
    4771             : 
    4772             : /* Petersson product on Hom_G(Delta_0, V_k) */
    4773             : GEN
    4774        1351 : mspetersson(GEN W, GEN F, GEN G)
    4775             : {
    4776        1351 :   pari_sp av = avma;
    4777             :   GEN vf, vg, C;
    4778             :   long k, l, tG, tF;
    4779        1351 :   checkms(W);
    4780        1351 :   if (!F) F = matid(msdim(W));
    4781        1351 :   if (!G) G = F;
    4782        1351 :   tF = typ(F);
    4783        1351 :   tG = typ(G);
    4784        1351 :   if (tF == t_MAT && tG != t_MAT) pari_err_TYPE("mspetersson",G);
    4785        1351 :   if (tG == t_MAT && tF != t_MAT) pari_err_TYPE("mspetersson",F);
    4786        1351 :   mspetersson_i(W, F, G, &vf, &vg, &C);
    4787        1351 :   l = lg(vf);
    4788        1351 :   if (tF != t_MAT)
    4789             :   { /* <F,G>, two symbols */
    4790        1274 :     GEN s = gen_0;
    4791        1274 :     for (k = 1; k < l; k++) s = gadd(s, bil(gel(vf,k), gel(vg,k), C));
    4792        1274 :     return gerepileupto(av, s);
    4793             :   }
    4794          77 :   else if (F != G)
    4795             :   { /* <(f_1,...,f_m), (g_1,...,g_n)> */
    4796           0 :     long iF, iG, lF = lg(F), lG = lg(G);
    4797           0 :     GEN M = cgetg(lG, t_MAT);
    4798           0 :     for (iG = 1; iG < lG; iG++)
    4799             :     {
    4800           0 :       GEN c = cgetg(lF, t_COL);
    4801           0 :       gel(M,iG) = c;
    4802           0 :       for (iF = 1; iF < lF; iF++)
    4803             :       {
    4804           0 :         GEN s = gen_0;
    4805           0 :         for (k = 1; k < l; k++)
    4806           0 :           s = gadd(s, bil(gmael(vf,k,iF), gmael(vg,k,iG), C));
    4807           0 :         gel(c,iF) = s; /* M[iF,iG] = <F[iF], G[iG] > */
    4808             :       }
    4809             :     }
    4810           0 :     return gerepilecopy(av, M);
    4811             :   }
    4812             :   else
    4813             :   { /* <(f_1,...,f_n), (f_1,...,f_n)> */
    4814          77 :     long iF, iG, n = lg(F)-1;
    4815          77 :     GEN M = zeromatcopy(n,n);
    4816         693 :     for (iG = 1; iG <= n; iG++)
    4817        3192 :       for (iF = iG+1; iF <= n; iF++)
    4818             :       {
    4819        2576 :         GEN s = gen_0;
    4820       14728 :         for (k = 1; k < l; k++)
    4821       12152 :           s = gadd(s, bil(gmael(vf,k,iF), gmael(vg,k,iG), C));
    4822        2576 :         gcoeff(M,iF,iG) = s; /* <F[iF], F[iG] > */
    4823        2576 :         gcoeff(M,iG,iF) = gneg(s);
    4824             :       }
    4825          77 :     return gerepilecopy(av, M);
    4826             :   }
    4827             : }

Generated by: LCOV version 1.11