Code coverage tests

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

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

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

LCOV - code coverage report
Current view: top level - basemath - pclgp.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.18.0 lcov report (development 29818-b3e15d99d2) Lines: 1838 2381 77.2 %
Date: 2024-12-27 09:09:37 Functions: 153 180 85.0 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2000, 2012  The PARI group.
       2             : 
       3             : This file is part of the PARI/GP package.
       4             : 
       5             : PARI/GP is free software; you can redistribute it and/or modify it under the
       6             : terms of the GNU General Public License as published by the Free Software
       7             : Foundation. It is distributed in the hope that it will be useful, but WITHOUT
       8             : ANY WARRANTY WHATSOEVER.
       9             : 
      10             : Check the License for details. You should have received a copy of it, along
      11             : with the package; see the file 'COPYING'. If not, write to the Free Software
      12             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      13             : 
      14             : #include "pari.h"
      15             : #include "paripriv.h"
      16             : 
      17             : #define DEBUGLEVEL DEBUGLEVEL_subcyclo
      18             : 
      19             : /* written by Takashi Fukuda */
      20             : 
      21             : #define onevec(x) const_vec(x,gen_1)
      22             : #define nullvec() cgetg(1, t_VEC)
      23             : #define order_f_x(f, x) Fl_order(x%f, eulerphiu(f), f)
      24             : 
      25             : #define USE_MLL         (1L<<0)
      26             : #define NO_PLUS_PART    (1L<<1)
      27             : #define NO_MINUS_PART   (1L<<2)
      28             : #define SKIP_PROPER     (1L<<3)
      29             : #define SAVE_MEMORY     (1L<<4)
      30             : #define USE_FULL_EL     (1L<<5)
      31             : #define USE_BASIS       (1L<<6)
      32             : #define USE_FACTOR      (1L<<7)
      33             : #define USE_GALOIS_POL  (1L<<8)
      34             : #define USE_F           (1L<<9)
      35             : 
      36             : static ulong
      37      178002 : _get_d(GEN H) { return umael(H, 2, 1);}
      38             : static ulong
      39      119820 : _get_f(GEN H) { return umael(H, 2, 2);}
      40             : static ulong
      41       30843 : _get_h(GEN H) { return umael(H, 2, 3);}
      42             : static long
      43       28392 : _get_s(GEN H) { return umael(H, 2, 4);}
      44             : static long
      45       53957 : _get_g(GEN H) { return umael(H, 2, 5);}
      46             : static GEN
      47       30787 : _get_H(GEN H) { return gel(H, 3);}
      48             : static ulong
      49      112755 : K_get_d(GEN K) { return _get_d(gel(K,1)); }
      50             : static ulong
      51       82965 : K_get_f(GEN K) { return _get_f(gel(K,1)); }
      52             : static ulong
      53       17424 : K_get_h(GEN K) { return _get_h(gel(K,1)); }
      54             : static long
      55           0 : K_get_s(GEN K) { return _get_s(gel(K,1)); }
      56             : static ulong
      57       17102 : K_get_g(GEN K) { return _get_g(gel(K,1)); }
      58             : static GEN
      59       17368 : K_get_H(GEN K) { return _get_H(gel(K,1)); }
      60             : static ulong
      61       59314 : K_get_dchi(GEN K) { return gel(K,6)[1]; }
      62             : static ulong
      63       23660 : K_get_nconj(GEN K) { return gel(K,6)[2]; }
      64             : 
      65             : /* G=<s> is a cyclic group of order n and t=s^(-1).
      66             :  *  convert sum_i a_i*s^i to sum_i b_i*t^i */
      67             : static GEN
      68          14 : Flx_recip1_inplace(GEN x, long pn)
      69             : {
      70          14 :   long i, lx = lg(x);
      71          14 :   if(lx-2 != pn) /* This case scarcely occurs */
      72             :   {
      73           0 :     long ly = pn+2;
      74           0 :     GEN y = const_vecsmall(ly, 0);
      75           0 :     y[1] = x[1];y[2] = x[2];
      76           0 :     for(i=3;i<lx;i++) y[ly+2-i] = x[i];
      77           0 :     return Flx_renormalize(y, ly);
      78             :   }
      79             :   else /* almost all cases */
      80             :   {
      81          14 :     long t, mid = (lx+1)>>1;
      82        7168 :     for(i=3;i<=mid;i++)
      83             :     {
      84        7154 :       t = x[i];x[i] = x[lx+2-i];x[lx+2-i] = t;
      85             :     }
      86          14 :     return Flx_renormalize(x, lx);
      87             :   }
      88             : }
      89             : 
      90             : /* Return h^degpol(P) P(x / h) */
      91             : static GEN
      92          14 : Flx_rescale_inplace(GEN P, ulong h, ulong p)
      93             : {
      94          14 :   long i, l = lg(P);
      95          14 :   ulong hi = h;
      96       14322 :   for (i=l-2; i>=2; i--)
      97             :   {
      98       14322 :     P[i] = Fl_mul(P[i], hi, p);
      99       14322 :     if (i == 2) break;
     100       14308 :     hi = Fl_mul(hi,h, p);
     101             :   }
     102          14 :   return P;
     103             : }
     104             : 
     105             : static GEN
     106          14 : zx_to_Flx_inplace(GEN x, ulong p)
     107             : {
     108          14 :   long i, lx = lg(x);
     109       14350 :   for (i=2; i<lx; i++) uel(x,i) = umodsu(x[i], p);
     110          14 :   return Flx_renormalize(x, lx);
     111             : }
     112             : 
     113             : /* zero pol of n components (i.e. deg=n-1). need to pass to ZX_renormalize */
     114             : INLINE GEN
     115       53263 : pol_zero(long n)
     116             : {
     117             :   long i;
     118       53263 :   GEN p = cgetg(n+2, t_POL);
     119       53263 :   p[1] = evalsigne(1) | evalvarn(0);
     120     1799812 :   for (i = 2; i < n+2; i++) gel(p, i) = gen_0;
     121       53263 :   return p;
     122             : }
     123             : 
     124             : /* e[i+1] = L*i + K for i >= n; determine K,L and reduce n if possible */
     125             : static GEN
     126          28 : vecsmall2vec2(GEN e, long n)
     127             : {
     128          28 :   long L = e[n+1] - e[n], K = e[n+1] - L*n;
     129          56 :   n--; while (n >= 0 && e[n+1] - L*n == K) n--;
     130          28 :   if (n < 0) e = nullvec(); else { setlg(e, n+2); e = zv_to_ZV(e); }
     131          28 :   return mkvec3(utoi(L), stoi(K), e); /* L >= 0 */
     132             : }
     133             : 
     134             : /* z=zeta_{p^n}; return k s.t. (z-1)^k || f(z) assuming deg(f)<phi(p^n) */
     135             : static long
     136          49 : zx_p_val(GEN f, ulong p, ulong n)
     137             : {
     138          49 :   pari_sp av = avma;
     139          49 :   ulong x = zx_lval(f, p);
     140          49 :   if (x) { f = zx_z_divexact(f, upowuu(p, x)); x *= (p-1)*upowuu(p, n-1); }
     141          49 :   x += Flx_val(Flx_translate1(zx_to_Flx(f, p), p));
     142          49 :   return gc_long(av, x);
     143             : }
     144             : 
     145             : static long
     146         315 : ZX_p_val(GEN f, ulong p, ulong n)
     147             : {
     148         315 :   pari_sp av = avma;
     149         315 :   ulong x = ZX_lval(f, p);
     150         315 :   if (x) { f = ZX_Z_divexact(f, powuu(p, x)); x *= (p-1)*upowuu(p, n-1); }
     151         315 :   x += Flx_val(Flx_translate1(ZX_to_Flx(f, p), p));
     152         315 :   return gc_long(av, x);
     153             : }
     154             : 
     155             : static GEN
     156          35 : set_A(GEN B, int *chi)
     157             : {
     158          35 :   long a, i, j, B1 = B[1], l = lg(B);
     159          35 :   GEN A = cgetg(l, t_VECSMALL);
     160     1687014 :   for (a = 0, j = 1; j < B1; j++) a += chi[j];
     161          35 :   A[1] = a;
     162         714 :   for (i = 2; i < l; i++)
     163             :   {
     164         679 :     long Bi = B[i];
     165    28159243 :     for (a = A[i-1], j = B[i-1]; j < Bi; j++) a += chi[j];
     166         679 :     A[i] = a;
     167             :   }
     168          35 :   return A;
     169             : }
     170             : 
     171             : /* g_n(a)=g_n(b) <==> a^2=b^2 mod 2^(n+2) <==> a=b,-b mod 2^(n+2)
     172             :  * g_n(a)=g_n(1+q0)^k <==> a=x(1+q0)^k x=1,-1
     173             :  * gam[1+a]=k, k<0  ==> g_n(a)=0
     174             :  *             k>=0 ==> g_n(a)^(-1)=gamma^k, gamma=g_n(1+q0) */
     175             : static GEN
     176          14 : set_gam2(long q01, long n)
     177             : {
     178             :   long i, x, x1, pn, pn2;
     179             :   GEN gam;
     180          14 :   pn = (1L<<n);
     181          14 :   pn2 = (pn<<2);
     182          14 :   gam = const_vecsmall(pn2, -1);
     183          14 :   x=Fl_inv(q01, pn2); x1=1;
     184       14350 :   for (i=0; i<pn; i++)
     185             :   {
     186       14336 :     gam[1+x1] = gam[1+Fl_neg(x1, pn2)] = i;
     187       14336 :     x1 = Fl_mul(x1, x, pn2);
     188             :   }
     189          14 :   return gam;
     190             : }
     191             : 
     192             : /* g_n(a)=g_n(b) <==> a^(p-1)=b^(p-1) mod p^(n+1) <==> a=xb x=<g^(p^n)>
     193             :  * g_n(a)=g_n(1+q0)^k <==> a=x(1+q0)^k x=<g^(p^n)>
     194             :  * gam[1+a]=k, k<0  ==> g_n(a)=0
     195             :  *             k>=0 ==> g_n(a)^(-1)=gamma^k, gamma=g_n(1+q0) */
     196             : static GEN
     197         476 : set_gam(long q01, long p, long n)
     198             : {
     199             :   long i, j, g, g1, x, x1, p1, pn, pn1;
     200             :   GEN A, gam;
     201         476 :   p1 = p-1; pn = upowuu(p, n); pn1 = p*pn;
     202         476 :   gam = const_vecsmall(pn1, -1);
     203         476 :   g = pgener_Zl(p); g1 = Fl_powu(g, pn, pn1);
     204         476 :   A = Fl_powers(g1, p1-1, pn1);  /* A[1+i]=g^(i*p^n) mod p^(n+1), 0<=i<=p-2 */
     205         476 :   x = Fl_inv(q01, pn1); x1 = 1;
     206      568694 :   for (i=0; i<pn; i++)
     207             :   {
     208     2086126 :     for (j=1; j<=p1; j++) gam[1+Fl_mul(x1, A[j], pn1)] = i;
     209      568218 :     x1 = Fl_mul(x1, x, pn1);
     210             :   }
     211         476 :   return gam;
     212             : }
     213             : 
     214             : /* k=Q(sqrt(m)), A_n=p-class gr. of k_n, |A_n|=p^(e_n)
     215             :  * return e_n-e_(n-1)
     216             :  * essential assumption : m is not divisible by p
     217             :  * Gold, Acta Arith. XXVI (1974), p.22 formula (3) */
     218             : static long
     219          35 : ediff(ulong p, long m, ulong n, int *chi)
     220             : {
     221          35 :   pari_sp av = avma;
     222             :   long j, lx, *px;
     223             :   ulong i, d, s, y, g, p1, pn, pn1, pn_1, phipn, phipn1;
     224             :   GEN A, B, x, gs, cs;
     225             : 
     226          35 :   d=((m-1)%4==0)?labs(m):4*labs(m);
     227          35 :   p1=p-1; pn_1=upowuu(p, n-1); pn=p*pn_1; pn1=p*pn; phipn=p1*pn_1; phipn1=p1*pn;
     228          35 :   lx=2*p1*phipn;
     229          35 :   y=Fl_inv(pn1%d, d); g=pgener_Zl(p);  /* pn1 may > d */
     230          35 :   cs = cgetg(2+phipn, t_VECSMALL); cs[1] = evalvarn(0);
     231          35 :   x = cgetg(1+lx, t_VECSMALL);
     232          35 :   gs = Fl_powers(g, phipn1-1, pn1); /* gs[1+i]=g^i(mod p^(n+1)), 0<=i<p^(n+1) */
     233             : 
     234         105 :   for (px=x,i=0; i<p1; i++)
     235             :   {
     236          70 :     long ipn=i*pn+1,ipnpn=ipn+phipn;
     237         546 :     for (s=0; s<phipn; s++)
     238             :     {
     239         476 :       *++px = (y*gs[s+ipn])%d;  /* gs[s+ipn] may > d */
     240         476 :       *++px = (y*gs[(s%pn_1)+ipnpn])%d;
     241             :     }
     242             :   }
     243          35 :   B = vecsmall_uniq(x);
     244          35 :   A = set_A(B, chi);
     245         273 :   for (s=0; s<phipn; s++)
     246             :   {
     247         238 :     long a=0, ipn=1, spn1=s%pn_1;
     248         714 :     for (i=0; i<p1; i++)
     249             :     {
     250         476 :       if ((j=zv_search(B, (y*gs[s+ipn])%d))<=0)
     251           0 :         pari_err_BUG("zv_search failed\n");
     252         476 :       a+=A[j];
     253         476 :       if ((j=zv_search(B, (y*gs[spn1+ipn+phipn])%d))<=0)
     254           0 :         pari_err_BUG("zv_search failed\n");
     255         476 :       a-=A[j];
     256         476 :       ipn+=pn;
     257             :     }
     258         238 :     cs[2+s] = a;
     259             :   }
     260          35 :   cs = zx_renormalize(cs, lg(cs));
     261          35 :   y = (lg(cs)==3) ? phipn*z_lval(cs[2], p) : zx_p_val(cs, p, n);
     262          35 :   return gc_long(av, y);
     263             : }
     264             : 
     265             : static GEN
     266           0 : quadteichstk(GEN Chi, int *chi, GEN Gam, long p, long m, long n)
     267             : {
     268           0 :   GEN Gam1 = Gam+1, xi;
     269             :   long i, j, j0, d, f0, pn, pn1, deg, pn1d;
     270             : 
     271           0 :   d = ((m&3)==1)?m:m<<2;
     272           0 :   f0 = ulcm(p, d)/p;
     273           0 :   pn = upowuu(p, n); pn1 = p*pn; pn1d = pn1%d;
     274           0 :   xi = cgetg(pn+2, t_POL); xi[1] = evalsigne(1) | evalvarn(0);
     275           0 :   for (i=0; i<pn; i++) gel(xi, 2+i) = const_vecsmall(p, 0);
     276           0 :   for (j=1; j<pn1; j++)
     277             :   {
     278             :     long jp, ipn1d, *xij0;
     279           0 :     if ((j0 = Gam1[j])<0) continue;
     280           0 :     jp = j%p; ipn1d = j%d; xij0 = gel(xi, 2+j0)+2;
     281           0 :     for (i=1; i<f0; i++)
     282             :     {
     283             :       int sgn;
     284           0 :       if ((ipn1d += pn1d) >= d) ipn1d -= d;
     285           0 :       if ((sgn = chi[ipn1d])==0) continue;
     286           0 :       deg = Chi[jp];  /* jp!=0 because j0>=0 */
     287           0 :       if (sgn>0) xij0[deg] += i;
     288           0 :       else xij0[deg] -= i;
     289             :     }
     290             :   }
     291           0 :   for (i=0; i<pn; i++) gel(xi, 2+i) = zx_renormalize(gel(xi, 2+i), p+1);
     292           0 :   return FlxX_renormalize(xi, pn+2);  /* zxX_renormalize does not exist */
     293             : }
     294             : 
     295             : #ifdef DEBUG_QUADSTK
     296             : /* return f0*xi_n */
     297             : static GEN
     298             : quadstkp_by_def(int *chi, GEN gam, long n, long p, long f, long f0)
     299             : {
     300             :   long i, a, a1, pn, pn1, qn;
     301             :   GEN x, x2, gam1 = gam+1;
     302             :   pn = upowuu(p, n); pn1 = p*pn; qn = f0*pn1;
     303             :   x = const_vecsmall(pn+1, 0); x2 = x+2;
     304             :   for (a=1; a<qn; a++)
     305             :   {
     306             :     int sgn;
     307             :     if ((a1=gam1[a%pn1])<0 || (sgn=chi[a%f])==0) continue;
     308             :     if (sgn>0) x2[a1]+=a;
     309             :     else x2[a1]-=a;
     310             :   }
     311             :   for (i=0; i<pn; i++)
     312             :   {
     313             :     if (x2[i]%pn1) pari_err_BUG("stickel. ele. is not integral.\n");
     314             :     else x2[i]/=pn1;
     315             :   }
     316             :   return zx_renormalize(x, pn+2);
     317             : }
     318             : #endif
     319             : 
     320             : /* f!=p
     321             :  * xi_n = f0^(-1)*
     322             :  *   sum_{0<=j<pn1,(j,p)=1}(Q_n/Q,j)^(-1)*(sum_{0<=i<f0}i*chi^(-1)(pn1*i+j)) */
     323             : static GEN
     324          14 : quadstkp1(int *chi, GEN gam, long n, long p, long f, long f0)
     325             : {
     326             :   long i, j, j0, pn, pn1, pn1f, den;
     327             :   GEN x, x2;
     328          14 :   pn = upowuu(p, n); pn1 = p*pn; pn1f = pn1%f;
     329          14 :   x = const_vecsmall(pn+1, 0); x2 = x+2;
     330          14 :   if (f==3) den = (chi[p%f]>0)?f0<<1:2;
     331          14 :   else if (f==4) den = (chi[p%f]>0)?f0<<1:f0;
     332          14 :   else den = f0<<1;
     333     1653372 :   for (j=1; j<pn1; j++)
     334             :   {
     335             :     long ipn1;
     336     1653358 :     if (j%p==0) continue;
     337     1102248 :     j0 = gam[1+j]; ipn1 = j%f;
     338   263437272 :     for (i=1; i<f0; i++)
     339             :     {
     340             :       int sgn;
     341   262335024 :       if ((ipn1+=pn1f)>=f) ipn1-=f;
     342   262335024 :       if ((sgn = chi[ipn1])>0) x2[j0]+=i;
     343   131716319 :       else if (sgn<0) x2[j0]-=i;
     344             :     }
     345             :   }
     346      551138 :   for (i=0; i<pn; i++)
     347             :   {
     348      551124 :     if (x2[i]%den) pari_err_BUG("stickel. ele. is not integral.\n");
     349      551124 :     else x2[i]/=den;
     350             :   }
     351          14 :   return zx_renormalize(x, pn+2);
     352             : }
     353             : 
     354             : /* f==p */
     355             : static GEN
     356           0 : quadstkp2(int *chi, GEN gam, long n, long p)
     357             : {
     358             :   long a, a1, i, pn, pn1, amodp;
     359           0 :   GEN x, x2, gam1 = gam+1;
     360           0 :   pn = upowuu(p, n); pn1 = p*pn;
     361           0 :   x = const_vecsmall(pn+1, 0); x2 = x+2;
     362           0 :   for (a=1,amodp=0; a<pn1; a++)
     363             :   {
     364             :     int sgn;
     365           0 :     if (++amodp==p) {amodp = 0; continue; }
     366           0 :     if ((sgn = chi[amodp])==0) continue;
     367           0 :     a1=gam1[a];
     368           0 :     if (sgn>0) x2[a1]+=a;
     369           0 :     else x2[a1]-=a;
     370             :   }
     371           0 :   for (i=0; i<pn; i++)
     372             :   {
     373           0 :     if (x2[i]%pn1) pari_err_BUG("stickel. ele. is not integral.\n");
     374           0 :     else x2[i]/=pn1;
     375             :   }
     376           0 :   return zx_renormalize(x, pn+2);
     377             : }
     378             : 
     379             : /*  p>=3
     380             :  *  f = conductor of Q(sqrt(m))
     381             :  *  q0 = lcm(f,p) = f0*p
     382             :  *  qn = q0*p^n = f0*p^(n+1)
     383             :  *  xi_n = qn^(-1)*sum_{1<=a<=qn,(a,qn)=1} a*chi(a)^(-1)*(Q_n/Q,a)^(-1) */
     384             : static GEN
     385          14 : quadstkp(long p, long m, long n, int *chi)
     386             : {
     387             :   long f, f0, pn, pn1, q0;
     388             :   GEN gam;
     389          14 :   f = ((m-1)%4==0)?labs(m):4*labs(m);
     390          14 :   pn = upowuu(p, n); pn1 = p*pn;
     391          14 :   if (f % p) { q0 = f * p; f0 = f; } else { q0 = f; f0 = f / p; }
     392          14 :   gam = set_gam((1+q0)%pn1, p, n);
     393             : #ifdef DEBUG_QUADSTK
     394             :   return quadstkp_by_def(chi, gam, n, p, f, f0);
     395             : #else
     396          14 :   return (f0!=1)?quadstkp1(chi, gam, n, p, f, f0):quadstkp2(chi, gam, n, p);
     397             : #endif
     398             : }
     399             : 
     400             : /* p=2 */
     401             : static GEN
     402          14 : quadstk2(long m, long n, int *chi)
     403             : {
     404             :   long i, j, j0, f, f0, pn, pn1, pn2, pn2f, q0;
     405             :   GEN x, x2, gam;
     406          14 :   f = ((m-1)%4==0)?labs(m):4*labs(m);
     407          14 :   pn = 1L<<n; pn1 = pn<<1; pn2 = pn1<<1; pn2f = pn2%f;
     408          14 :   q0 = (f&1)?f*4:f;
     409          14 :   f0 = (f&1)?f:f/4;
     410          14 :   x = const_vecsmall(pn+1, 0); x2 = x+2;
     411          14 :   gam = set_gam2((1+q0)%pn2, n);
     412       57344 :   for (j=1; j<pn2; j++)
     413             :   {
     414             :     long ipn2;
     415       57330 :     if (!(j&1)) continue;
     416       28672 :     j0 = gam[1+j];
     417       28672 :     ipn2 = j%f;
     418             :     /* for (i=1; i<f0; i++) x2[j0]+=i*chi[(i*pn2+j)%f]; */
     419     1691648 :     for (i=1; i<f0; i++)
     420             :     {
     421             :       int sgn;
     422     1662976 :       if ((ipn2+=pn2f)>=f) ipn2-=f;
     423     1662976 :       if ((sgn=chi[ipn2])>0) x2[j0]+=i;
     424      845390 :       else if (sgn<0) x2[j0]-=i;
     425             :     }
     426             :   }
     427       14350 :   for (f0<<=1, i=0; i<pn; i++)
     428             :   {
     429       14336 :     if (x2[i]%f0) pari_err_BUG("stickel. ele. is not integral.\n");
     430       14336 :     else x2[i]/=f0;
     431             :   }
     432          14 :   return zx_renormalize(x, pn+2);
     433             : }
     434             : 
     435             : /* Chin is a generator of the group of the characters of G(Q_n/Q).
     436             :  * chin[1+a]=k, k<0  ==> Chin(a)=0
     437             :  *              k>=0 ==> Chin(a)=zeta_{p^n}^k */
     438             : static GEN
     439          28 : set_chin(long p, long n)
     440             : {
     441          28 :   long i, j, x = 1, g, gpn, pn, pn1;
     442             :   GEN chin, chin1;
     443          28 :   pn = upowuu(p, n); pn1 = p*pn;
     444          28 :   chin = const_vecsmall(pn1, -1); chin1 = chin+1;
     445          28 :   g = pgener_Zl(p); gpn = Fl_powu(g, pn, pn1);
     446         322 :   for (i=0; i<pn; i++)
     447             :   {
     448         294 :     long y = x;
     449         882 :     for (j=1; j<p; j++)
     450             :     {
     451         588 :       chin1[y] = i;
     452         588 :       y = Fl_mul(y, gpn, pn1);
     453             :     }
     454         294 :     x = Fl_mul(x, g, pn1);
     455             :   }
     456          28 :   return chin;
     457             : }
     458             : 
     459             : /* k=Q(sqrt(m)), A_n=p-class gr. of k_n, |A_n|=p^(e_n), p|m
     460             :  * return e_n-e_(n-1).
     461             :  * There is an another method using the Stickelberger element based on
     462             :  * Coates-Lichtenbaum, Ann. Math. vol.98 No.3 (1973), 498-550, Lemma 2.15.
     463             :  * If kro(m,p)!=1, then orders of two groups coincide.
     464             :  * ediff_ber is faster than the Stickelberger element. */
     465             : static long
     466          28 : ediff_ber(ulong p, long m, ulong n, int *chi)
     467             : {
     468          28 :   pari_sp av = avma;
     469             :   long a, d, e, x, y, pn, pn1, qn1;
     470          28 :   GEN B, B2, chin = set_chin(p, n)+1;
     471             : 
     472          28 :   d = ((m-1)%4==0)?labs(m):4*labs(m);
     473          28 :   pn = upowuu(p, n); pn1 = p*pn; qn1 = (d*pn)>>1;
     474          28 :   B = const_vecsmall(pn+1, 0); B2 = B+2;
     475   522106886 :   for (a=x=y=1; a <= qn1; a++) /* x=a%d, y=a%pn1 */
     476             :   {
     477   522106858 :     int sgn = chi[x];
     478   522106858 :     if (sgn)
     479             :     {
     480   172938444 :       long k = chin[y];
     481   172938444 :       if (k >= 0) { if (sgn > 0) B2[k]++; else B2[k]--; }
     482             :     }
     483   522106858 :     if (++x == d) x = 0;
     484   522106858 :     if (++y == pn1) y = 0;
     485             :   }
     486          28 :   B = zx_renormalize(B, pn+2);
     487          14 :   e = (n==1)? zx_p_val(B, p, n)
     488          28 :             : ZX_p_val(ZX_rem(zx_to_ZX(B), polcyclo(pn, 0)), p, n);
     489          28 :   if (p==3 && chi[2] < 0) e--;  /* 2 is a primitive root of 3^n (n>=1) */
     490          28 :   return gc_long(av, e);
     491             : }
     492             : 
     493             : #ifdef DEBUG
     494             : static int*
     495             : set_quad_chi_slow(long m)
     496             : {
     497             :   long a, d = (m-1) % 4? 4*m: m, f = labs(d);
     498             :   int *chi = (int*)stack_calloc(sizeof(int)*f);
     499             :   for (a=1; a<f; a++) chi[a]=kross(d, a);
     500             :   return chi;
     501             : }
     502             : #endif
     503             : 
     504             : /* chi[a] = kross(d, a), 0 <= a < f, d = disc Q(sqrt(m)), f = abs(d)
     505             :  * Algorithm: m=-p1*p2*...*pr ==> kross(d,gi)=-1 (1<=i<=r), gi=proot(pi) */
     506             : static int*
     507          56 : set_quad_chi(long m)
     508             : {
     509          56 :   long d = (m-1) % 4? 4*m: m, f = labs(d);
     510          56 :   GEN fa = factoru(f), P = gel(fa, 1), E = gel(fa,2), u, v;
     511          56 :   long i, j, np, nm, l = lg(P);
     512          56 :   int *chi = (int*)stack_calloc(sizeof(int)*f);
     513          56 :   pari_sp av = avma;
     514          56 :   int *plus = (int*)stack_calloc(sizeof(int)*f), *p0 = plus;
     515          56 :   int *minus = (int*)stack_calloc(sizeof(int)*f), *p1 = minus;
     516             : 
     517          56 :   u = cgetg(32, t_VECSMALL);
     518          56 :   v = cgetg(32, t_VECSMALL);
     519         168 :   for (i = 1; i < l; i++)
     520             :   {
     521         112 :     ulong q = upowuu(P[i], E[i]), fq = f / q;
     522         112 :     u[i] = q * Fl_inv(q % fq, fq); /* 1 mod f/q, 0 mod q */
     523         112 :     v[i] = Fl_sub(1, u[i], f); /* => gv + u is 1 mod f/q and g mod q */
     524             :   }
     525          56 :   if (E[1]==2) /* f=4*(-m) */
     526             :   {
     527          14 :     *p0++ = Fl_add(v[1], u[1], f);
     528          14 :     *p1++ = Fl_add(Fl_mul(3, v[1], f), u[1], f);
     529          14 :     i = 2;
     530             :   }
     531          42 :   else if (E[1]==3) /* f=8*(-m) */
     532             :   {
     533             :     ulong a;
     534           7 :     *p0++ = Fl_add(v[1], u[1], f);
     535           7 :     a = Fl_add(Fl_mul(3, v[1], f), u[1], f);
     536           7 :     if (kross(d, a) > 0) *p0++ = a; else *p1++ = a;
     537           7 :     a = Fl_add(Fl_mul(5, v[1], f), u[1], f);
     538           7 :     if (kross(d, a) > 0) *p0++ = a; else *p1++ = a;
     539           7 :     a = Fl_add(Fl_mul(7, v[1], f), u[1], f);
     540           7 :     if (kross(d, a) > 0) *p0++ = a; else *p1++ = a;
     541           7 :     i = 2;
     542             :   }
     543             :   else /* f=-m */
     544          35 :   {*p0++ = 1; i = 1; }
     545         147 :   for (; i < l; i++)
     546             :   {
     547          91 :     ulong gn, g = pgener_Fl(P[i]);
     548          91 :     gn = g = Fl_add(Fl_mul(g, v[i], f), u[i], f);
     549          91 :     np = p0-plus;
     550          91 :     nm = p1-minus;
     551             :     for (;;)
     552             :     {
     553     4802616 :       for (j = 0; j < np; j++) *p1++ = Fl_mul(plus[j], gn, f);
     554     4799858 :       for (j = 0; j < nm; j++) *p0++ = Fl_mul(minus[j], gn, f);
     555        7784 :       gn = Fl_mul(gn, g, f); if (gn == 1) break;
     556     4763927 :       for (j= 0; j < np; j++) *p0++ = Fl_mul(plus[j], gn, f);
     557     4761204 :       for (j = 0; j < nm; j++) *p1++ = Fl_mul(minus[j], gn, f);
     558        7693 :       gn = Fl_mul(gn, g, f); if (gn == 1) break;
     559             :     }
     560             :   }
     561          56 :   np = p0-plus;
     562          56 :   nm = p1-minus;
     563     9548427 :   for (i = 0; i < np; i++) chi[plus[i]] = 1;
     564     9548427 :   for (i = 0; i < nm; i++) chi[minus[i]] = -1;
     565          56 :   set_avma(av); return chi;
     566             : }
     567             : 
     568             : static long
     569        8995 : srh_x(GEN T, long n, long x)
     570             : {
     571       30086 :   for (; x<n; x++) if (!T[x]) return x;
     572         623 :   return -1;
     573             : }
     574             : 
     575             : /* G is a cyclic group of order d. hat(G)=<chi>.
     576             :  * chi, chi^p, ... , chi^(p^(d_chi-1)) are conjugate.
     577             :  * {chi^j | j in C} are repre. of Q_p-congacy classes of inj. chars.
     578             :  *
     579             :  * C is a set of representatives of H/<p>, where H=(Z/dZ)^* */
     580             : static GEN
     581        1134 : set_C(long p, long d, long d_chi, long n_conj)
     582             : {
     583        1134 :   long i, j, x, y, pmodd = p%d;
     584        1134 :   GEN T = const_vecsmall(d, 0)+1;
     585        1134 :   GEN C = cgetg(1+n_conj, t_VECSMALL);
     586        1134 :   if (n_conj==1) { C[1] = 1; return C; }
     587        9618 :   for (i=0, x=1; x >= 0; x = srh_x(T, d, x))
     588             :   {
     589        8995 :     if (cgcd(x, d)==1) C[++i] = x;
     590       40929 :     for (j=0, y=x; j<d_chi; j++) T[y = Fl_mul(y, pmodd, d)] = 1;
     591             :   }
     592         623 :   return C;
     593             : }
     594             : 
     595             : static GEN
     596         343 : FpX_one_cyclo(long n, GEN p)
     597             : {
     598         343 :   if (lgefint(p)==3)
     599         301 :     return Flx_to_ZX(Flx_factcyclo(n, p[2], 1));
     600             :   else
     601          42 :     return FpX_factcyclo(n, p, 1);
     602             : }
     603             : 
     604             : static void
     605       17094 : Flx_red_inplace(GEN x, ulong p)
     606             : {
     607       17094 :   long i, l = lg(x);
     608      274540 :   for (i=2; i<l; i++) x[i] = uel(x, i)%p;
     609       17094 :   Flx_renormalize(x, l);
     610       17094 : }
     611             : 
     612             : /* x[i], T[i] < pn */
     613             : static GEN
     614       39046 : Flxq_xi_conj(GEN x, GEN T, long j, long d, long pn)
     615             : {
     616       39046 :   long i, deg = degpol(x);
     617       39046 :   GEN z = const_vecsmall(d+1, 0);
     618     1032304 :   for (i=0; i<=deg; i++) z[2+Fl_mul(i, j, d)] = x[2+i];
     619       39046 :   return Flx_rem(Flx_renormalize(z, d+2), T, pn);
     620             : }
     621             : 
     622             : static GEN
     623         966 : FlxqX_xi_conj(GEN x, GEN T, long j, long d, long pn)
     624             : {
     625         966 :   long i, l = lg(x);
     626             :   GEN z;
     627         966 :   z = cgetg(l, t_POL); z[1] = evalsigne(1) | evalvarn(0);
     628       40012 :   for (i=2; i<l; i++) gel(z, i) = Flxq_xi_conj(gel(x, i), T, j, d, pn);
     629         966 :   return z;
     630             : }
     631             : 
     632             : static GEN
     633           0 : FlxqX_xi_norm(GEN x, GEN T, long p, long d, long pn)
     634             : {
     635           0 :   long i, d_chi = degpol(T);
     636           0 :   GEN z = x, z1 = x;
     637           0 :   for (i=1; i<d_chi; i++)
     638             :   {
     639           0 :     z1 = FlxqX_xi_conj(z1, T, p, d, pn);
     640           0 :     z = FlxqX_mul(z, z1, T, pn);
     641             :   }
     642           0 :   return z;
     643             : }
     644             : 
     645             : /* assume 0 <= x[i], y[j] <= m-1 */
     646             : static GEN
     647          15 : FpV_shift_add(GEN x, GEN y, GEN m, long start, long end)
     648             : {
     649             :   long i, j;
     650      222300 :   for (i=start, j=1; i<=end; i++, j++)
     651             :   {
     652      222285 :     pari_sp av = avma;
     653      222285 :     GEN z = addii(gel(x, i), gel(y, j));
     654      222285 :     gel(x, i) = (cmpii(z, m) >= 0)? gerepileuptoint(av, subii(z, m)): z;
     655             :   }
     656          15 :   return x;
     657             : }
     658             : 
     659             : /* assume 0 <= x[i], y[j] <= m-1 */
     660             : static GEN
     661          10 : FpV_shift_sub(GEN x, GEN y, GEN m, long start, long end)
     662             : {
     663             :   long i, j;
     664      112430 :   for (i=start, j=1; i<=end; i++, j++)
     665             :   {
     666      112420 :     pari_sp av = avma;
     667      112420 :     GEN z = subii(gel(x, i), gel(y, j));
     668      112420 :     gel(x, i) = (signe(z) < 0)? gerepileuptoint(av, addii(z, m)): z;
     669             :   }
     670          10 :   return x;
     671             : }
     672             : 
     673             : /* assume 0 <= x[i], y[j] <= m-1 */
     674             : static GEN
     675         173 : Flv_shift_add(GEN x, GEN y, ulong m, long start, long end)
     676             : {
     677             :   long i, j;
     678     2320113 :   for (i=start, j=1; i<=end; i++, j++)
     679             :   {
     680     2319940 :     ulong xi = x[i], yj = y[j];
     681     2319940 :     x[i] = Fl_add(xi, yj, m);
     682             :   }
     683         173 :   return x;
     684             : }
     685             : 
     686             : /* assume 0 <= x[i], y[j] <= m-1 */
     687             : static GEN
     688         102 : Flv_shift_sub(GEN x, GEN y, ulong m, long start, long end)
     689             : {
     690             :   long i, j;
     691     1165182 :   for (i=start, j=1; i<=end; i++, j++)
     692             :   {
     693     1165080 :     ulong xi = x[i], yj = y[j];
     694     1165080 :     x[i] = Fl_sub(xi, yj, m);
     695             :   }
     696         102 :   return x;
     697             : }
     698             : 
     699             : /* return 0 if p|x. else return 1 */
     700             : INLINE long
     701         896 : Flx_divcheck(GEN x, ulong p)
     702             : {
     703         896 :   long i, l = lg(x);
     704         910 :   for (i=2; i<l; i++) if (uel(x, i)%p) return 1;
     705         448 :   return 0;
     706             : }
     707             : 
     708             : static long
     709         448 : FlxX_weier_deg(GEN pol, long p)
     710             : {
     711         448 :   long i, l = lg(pol);
     712         896 :   for (i=2; i<l && Flx_divcheck(gel(pol, i), p)==0; i++);
     713         448 :   return (i<l)?i-2:-1;
     714             : }
     715             : 
     716             : static long
     717        1582 : Flx_weier_deg(GEN pol, long p)
     718             : {
     719        1582 :   long i, l = lg(pol);
     720        3997 :   for (i=2; i<l && pol[i]%p==0; i++);
     721        1582 :   return (i<l)?i-2:-1;
     722             : }
     723             : 
     724             : static GEN
     725         308 : Flxn_shift_mul(GEN g, long n, GEN p, long d, long m)
     726             : {
     727         308 :   return Flx_shift(Flxn_mul(g, p, d, m), n);
     728             : }
     729             : 
     730             : INLINE long
     731        1057 : deg_trunc(long lam, long p, long n, long pn)
     732             : {
     733             :   long r, x, d;
     734        1260 :   for (r=1,x=p; x<lam; r++) x *= p;  /* r is min int s.t. lam<=p^r */
     735        1057 :   if ((d = (n-r+2)*lam+1)>=pn) d = pn;
     736        1057 :   return d;
     737             : }
     738             : 
     739             : /*  Flx_translate1_basecase(g, pn) becomes slow when degpol(g)>1000.
     740             :  *  So I wrote Flxn_translate1().
     741             :  *  I need lambda to truncate pol.
     742             :  *  But I need to translate T --> 1+T to know lambda.
     743             :  *  Though the code has a little overhead, it is still fast. */
     744             : static GEN
     745         756 : Flxn_translate1(GEN g, long p, long n)
     746             : {
     747             :   long i, j, d, lam, pn, start;
     748         756 :   if (n==1) start = 3;
     749          70 :   else if (n==2) start = 9;
     750          70 :   else start = 10;
     751         756 :   pn = upowuu(p, n);
     752         756 :   for (lam=start; lam; lam<<=1)  /* least upper bound is 3 */
     753             :   {
     754             :     GEN z;
     755         756 :     d = deg_trunc(lam, p, n, pn);
     756         756 :     z = const_vecsmall(d+1, 0);  /* z[2],...,z[d+1] <--> a_0,...,a_{d-1} */
     757       44954 :     for (i=degpol(g); i>=0; i--)
     758             :     {
     759     1683486 :       for (j=d+1; j>2; j--) z[j] = Fl_add(z[j], z[j-1], pn);  /* z = z*(1+T) */
     760       44198 :       z[2] = Fl_add(z[2], g[2+i], pn);
     761             :     }
     762         756 :     z = Flx_renormalize(z, d+2);
     763         756 :     if (Flx_weier_deg(z, p) <= lam) return z;
     764             :   }
     765             :   return NULL; /*LCOV_EXCL_LINE*/
     766             : }
     767             : 
     768             : static GEN
     769         224 : FlxXn_translate1(GEN g, long p, long n)
     770             : {
     771             :   long i, j, d, lam, pn, start;
     772             :   GEN z;
     773         224 :   if (n==1) start = 3;
     774           0 :   else if (n==2) start = 9;
     775           0 :   else start = 10;
     776         224 :   pn = upowuu(p, n);
     777         224 :   for (lam=start; lam; lam<<=1)  /* least upper bound is 3 */
     778             :   {
     779         224 :     d = deg_trunc(lam, p, n, pn);
     780         224 :     z = const_vec(d+1, pol0_Flx(0));  /* z[2],...,z[d+1] <--> a_0,...,a_{d-1} */
     781         224 :     settyp(z, t_POL); z[1] = evalsigne(1) | evalvarn(0);
     782        9408 :     for (i=degpol(g); i>=0; i--)
     783             :     {
     784       64288 :       for (j=d+1; j>2; j--) gel(z, j) = Flx_add(gel(z, j), gel(z, j-1), pn);
     785        9184 :       gel(z, 2) = Flx_add(gel(z, 2), gel(g, 2+i), pn);
     786             :     }
     787         224 :     z = FlxX_renormalize(z, d+2);
     788         224 :     if (FlxX_weier_deg(z, p) <= lam) return z;
     789             :   }
     790             :   return NULL; /*LCOV_EXCL_LINE*/
     791             : }
     792             : 
     793             : /* lam < 0 => error (lambda can't be determined)
     794             :  * lam = 0 => return 1
     795             :  * lam > 0 => return dist. poly. of degree lam. */
     796             : static GEN
     797          84 : Flxn_Weierstrass_prep(GEN g, long p, long n, long d_chi)
     798             : {
     799          84 :   long i, r0, d, dg = degpol(g), lam, pn, t;
     800             :   ulong lam0;
     801             :   GEN U, UINV, P, PU, g0, g1, gp, gU;
     802          84 :   if ((lam = Flx_weier_deg(g, p))==0) return(pol1_Flx(0));
     803          77 :   else if (lam<0)
     804           0 :     pari_err(e_MISC,"Flxn_Weierstrass_prep: precision too low. Increase n!");
     805          77 :   lam0 = lam/d_chi;
     806          77 :   pn = upowuu(p, n);
     807          77 :   d = deg_trunc(lam, p, n, pn);
     808          77 :   if (d>dg) d = dg;
     809          77 :   if (d<=lam) d=1+lam;
     810         140 :   for (r0=1; upowuu(p, r0)<lam0; r0++);
     811          77 :   g = Flxn_red(g, d);
     812          77 :   t = Fl_inv(g[2+lam], pn);
     813          77 :   g = Flx_Fl_mul(g, t, pn);  /* normalized so as g[2+lam]=1 */
     814          77 :   U = Flx_shift(g, -lam);
     815          77 :   UINV = Flxn_inv(U, d, pn);
     816          77 :   P = zx_z_divexact(Flxn_red(g, lam), p);  /* assume g[i] <= LONG_MAX */
     817          77 :   PU = Flxn_mul(P, UINV, d, pn);
     818          77 :   gU = Flxn_mul(g, UINV, d, pn);
     819          77 :   g0 = pol1_Flx(0);
     820          77 :   g1 = pol1_Flx(0);
     821         385 :   for (i=1; i<n; i++)
     822             :   {
     823         308 :     g1 = Flxn_shift_mul(g1, -lam, PU, d, pn);
     824         308 :     gp = Flx_Fl_mul(g1, upowuu(p, i), pn);
     825         308 :     g0 = (i&1)?Flx_sub(g0, gp, pn):Flx_add(g0, gp, pn);
     826             :   }
     827          77 :   g0 = Flxn_mul(g0, gU, lam+1, pn);
     828          77 :   g0 = Flx_red(g0, upowuu(p, (p==2)?n-r0:n+1-r0));
     829          77 :   return g0;
     830             : }
     831             : 
     832             : /* xi_n and Iwasawa pol. for Q(sqrt(m)) and p
     833             :  *
     834             :  * (flag&1)!=0 ==> output xi_n
     835             :  * (flag&2)!=0 ==> output power series
     836             :  * (flag&4)!=0 ==> output Iwasawa polynomial */
     837             : static GEN
     838          14 : imagquadstkpol(long p, long m, long n)
     839             : {
     840          14 :   long pn = upowuu(p, n);
     841             :   GEN pol, stk, stk2;
     842             :   int *chi;
     843          14 :   if (p==2 && (m==-1 || m==-2 || m==-3 || m==-6)) return nullvec();
     844          14 :   if (p==3 && m==-3) return nullvec();
     845          14 :   if (p==2 && m%2==0) m /= 2;
     846          14 :   chi = set_quad_chi(m);
     847          14 :   stk = (p==2)? quadstk2(m, n, chi): quadstkp(p, m, n, chi);
     848          14 :   stk2 = zx_to_Flx(stk, pn);
     849          14 :   pol = Flxn_Weierstrass_prep(zlx_translate1(stk2, p, n), p, n, 1);
     850          14 :   return degpol(pol)? mkvec(Flx_to_ZX(pol)): nullvec();
     851             : }
     852             : 
     853             : /* a mod p == g^i mod p ==> omega(a)=zeta_(p-1)^(-i)
     854             :  *  Chi[g^i mod p]=i (0 <= i <= p-2) */
     855             : static GEN
     856           0 : get_teich(long p, long g)
     857             : {
     858           0 :   long i, gi = 1, p1 = p-1;
     859           0 :   GEN Chi = cgetg(p, t_VECSMALL);
     860           0 :   for (i=0; i<p1; i++) { Chi[gi] = i; gi = Fl_mul(gi, g, p); }
     861           0 :   return Chi;
     862             : }
     863             : 
     864             : /* Ichimura-Sumida criterion for Greenberg conjecture for real quadratic field.
     865             :  * chi: character of Q(sqrt(m)), omega: Teichmuller character mod p or 4.
     866             :  * Get Stickelberger element from chi^* = omega*chi^(-1) and convert it to
     867             :  * power series by the correspondence (Q_n/Q,1+q0)^(-1) <-> (1+T)(1+q0)^(-1) */
     868             : static GEN
     869          14 : realquadstkpol(long p, long m, long n)
     870             : {
     871             :   int *chi;
     872          14 :   long pnm1 = upowuu(p, n-1),pn = p*pnm1, pn1 = p*pn, d, q0;
     873             :   GEN stk, ser, pol;
     874          14 :   if (m==1) pari_err_DOMAIN("quadstkpol", "m", "=", gen_1, gen_1);
     875          14 :   if (p==2 && (m&1)==0) m>>=1;
     876          14 :   d = ((m&3)==1)?m:m<<2;
     877          14 :   q0 = ulcm((p==2)?4:p, d);
     878          14 :   if (p==2)
     879             :   {
     880          14 :     chi = set_quad_chi(-m);
     881          14 :     stk = quadstk2(-m, n, chi);
     882          14 :     stk = zx_to_Flx_inplace(stk, pn);
     883             :   }
     884           0 :   else if (p==3 && m%3==0 && kross(-m/3,3)==1)
     885           0 :   {
     886           0 :     long m3 = m/3;
     887           0 :     chi = set_quad_chi(-m3);
     888           0 :     stk = quadstkp(3, -m3, n, chi);
     889           0 :     stk = zx_to_Flx_inplace(stk, pn);
     890             :   }
     891             :   else
     892             :   {
     893           0 :     long g = pgener_Zl(p);
     894           0 :     long x = Fl_powu(Fl_inv(g, p), pnm1, pn);
     895           0 :     GEN Chi = get_teich(p, g);
     896           0 :     GEN Gam = set_gam((1+q0)%pn1, p, n);
     897           0 :     chi = set_quad_chi(m);
     898           0 :     stk = quadteichstk(Chi, chi, Gam, p, m, n);  /* exact */
     899           0 :     stk = zxX_to_FlxX(stk, pn);  /* approx. */
     900           0 :     stk = FlxY_evalx(stk, x, pn);
     901             :   }
     902          14 :   stk = Flx_rescale_inplace(Flx_recip1_inplace(stk, pn), (1+q0)%pn, pn);
     903          14 :   ser = Flxn_translate1(stk, p, n);
     904          14 :   pol = Flxn_Weierstrass_prep(ser, p, n, 1);
     905          14 :   return degpol(pol)? mkvec(Flx_to_ZX(pol)): nullvec();
     906             : }
     907             : 
     908             : /* m > 0 square-free. lambda_2(Q(sqrt(-m)))
     909             :  * Kida, Tohoku Math. J. vol.31 (1979), 91-96, Theorem 1. */
     910             : static GEN
     911           0 : quadlambda2(long m)
     912             : {
     913             :   long i, l, L;
     914             :   GEN P;
     915           0 :   if ((m&1)==0) m >>= 1;  /* lam_2(Q(sqrt(-m)))=lam_2(Q(sqrt(-2*m))) */
     916           0 :   if (m <= 3) return mkvecs(0);
     917           0 :   P = gel(factoru(m), 1); l = lg(P);
     918           0 :   for (L = -1,i = 1; i < l; i++) L += 1L << (-3 + vals(P[i]-1) + vals(P[i]+1));
     919           0 :   return mkvecs(L);
     920             : }
     921             : 
     922             : /* Iwasawa lambda invariant of Q(sqrt(m)) (m<0) for p
     923             :  * |A_n|=p^(e[n])
     924             :  * kross(m,p)!=1 : e[n]-e[n-1]<eulerphi(p^n)  ==> lambda=e[n]-e[n-1]
     925             :  * kross(m,p)==1 : e[n]-e[n-1]<=eulerphi(p^n) ==> lambda=e[n]-e[n-1]
     926             :  * Gold, Acta Arith. XXVI (1974), p.25, Cor. 3
     927             :  * Gold, Acta Arith. XXVI (1975), p.237, Cor. */
     928             : static GEN
     929          28 : quadlambda(long p, long m)
     930             : {
     931             :   long flag, n, phipn;
     932          28 :   GEN e = cgetg(31, t_VECSMALL);
     933             :   int *chi;
     934          28 :   if (m>0) pari_err_IMPL("plus part of lambda invariant in quadlambda()");
     935          28 :   if (p==2) return quadlambda2(-m);
     936          28 :   if (p==3 && m==-3) return mkvec3(gen_0, gen_0, nullvec());
     937          28 :   flag = kross(m, p);
     938          28 :   e[1] = Z_lval(quadclassno(quaddisc(stoi(m))), p);
     939          28 :   if (flag!=1 && e[1]==0) return mkvec3(gen_0, gen_0, nullvec());
     940          28 :   chi = set_quad_chi(m);
     941          28 :   phipn = p-1;  /* phipn=phi(p^n) */
     942          63 :   for (n=1; n; n++, phipn *= p)
     943             :   {
     944          63 :     long L = flag? ediff(p, m, n, chi): ediff_ber(p, m, n, chi);
     945          63 :     e[n+1] = e[n] + L;
     946          63 :     if ((flag!=1 && (L < phipn))|| (flag==1 && (L <= phipn))) break;
     947             :   }
     948          28 :   return vecsmall2vec2(e, n);
     949             : }
     950             : 
     951             : /* factor n-th cyclotomic polynomial mod p^r and return a minimal
     952             :  *  polynomial of zeta_n over Q_p.
     953             :  *  phi(n)=deg*n_conj, n_conj == 1 <=> polcyclo(n) is irred mod p. */
     954             : static GEN
     955         945 : set_minpol(ulong n, GEN p, ulong r, long n_conj)
     956             : {
     957             :   GEN z, v, pol, pr;
     958             :   pari_timer ti;
     959         945 :   if (umodiu(p, n)==1) /* zeta_n in Z_p, faster than polcyclo() */
     960             :   {
     961         420 :     GEN prm1 = powiu(p, r-1), pr = mulii(prm1, p); /* pr=p^r */
     962         420 :     GEN prn = diviuexact(subii(pr, prm1), n);      /* prn=phi(p^r)/n */
     963         420 :     z = Fp_pow(pgener_Fp(p), prn, pr);
     964         420 :     return deg1pol_shallow(gen_1, Fp_neg(z, pr), 0);
     965             :   }
     966         525 :   pr = powiu(p, r);
     967         525 :   pol = polcyclo(n, 0);
     968         525 :   if (n_conj==1) return FpX_red(pol, pr);
     969         343 :   if (DEBUGLEVEL>3) timer_start(&ti);
     970         343 :   z = FpX_one_cyclo(n, p);
     971         343 :   if (DEBUGLEVEL>3) timer_printf(&ti, "FpX_one_cyclo:n=%ld  ", n);
     972         343 :   v = ZpX_liftfact(pol, mkvec2(z, FpX_div(pol, z, p)), pr, p, r);
     973         343 :   return gel(v, 1);
     974             : }
     975             : 
     976             : static GEN
     977          91 : set_minpol_teich(ulong g_K, GEN p, ulong r)
     978             : {
     979          91 :   GEN prm1 = powiu(p, r-1), pr = mulii(prm1, p), z;
     980          91 :   z = Fp_pow(Fp_inv(utoi(g_K), p), prm1, pr);
     981          91 :   return deg1pol_shallow(gen_1, Fp_neg(z, pr), 0);
     982             : }
     983             : 
     984             : static long
     985       18963 : srh_1(GEN H)
     986             : {
     987       18963 :   GEN bits = gel(H, 3);
     988       18963 :   ulong f = bits[1];
     989       18963 :   return F2v_coeff(bits, f-1);
     990             : }
     991             : 
     992             : /* (1/f)sum_{1<=a<=f}a*chi^{-1}(a) = -(1/(2-chi(a)))sum_{1<=a<=f/2} chi^{-1}(a)
     993             :  *  does not overflow */
     994             : static GEN
     995       13146 : zx_ber_num(GEN Chi, long f, long d)
     996             : {
     997       13146 :   long i, f2 = f>>1;
     998       13146 :   GEN x = const_vecsmall(d+1, 0), x2 = x+2;
     999    51965081 :   for (i = 1; i <= f2; i++)
    1000    51951935 :     if (Chi[i] >= 0) x2[Chi[i]] ++;
    1001       13146 :   return zx_renormalize(x, d+2);
    1002             : }
    1003             : 
    1004             : /* x a zx
    1005             :  * zx_ber_num is O(f). ZX[FpX,Flx]_ber_conj is O(d). Sometimes d<<f. */
    1006             : static GEN
    1007       26257 : ZX_ber_conj(GEN x, long j, long d)
    1008             : {
    1009       26257 :   long i, deg = degpol(x);
    1010       26257 :   GEN y = pol_zero(d), x2 = x+2, y2 = y+2;
    1011      818202 :   for (i=0; i<=deg; i++) gel(y2, Fl_mul(i, j, d)) = stoi(x2[i]);
    1012       26257 :   return ZX_renormalize(y, d+2);
    1013             : }
    1014             : 
    1015             : /* x a zx */
    1016             : static GEN
    1017         252 : FpX_ber_conj(GEN x, long j, long d, GEN p)
    1018             : {
    1019         252 :   long i, deg = degpol(x);
    1020         252 :   GEN y = pol_zero(d), x2 = x+2, y2 = y+2;
    1021         756 :   for (i=0; i<=deg; i++) gel(y2, Fl_mul(i, j, d)) = modsi(x2[i], p);
    1022         252 :   return FpX_renormalize(y, d+2);
    1023             : }
    1024             : 
    1025             : /* x a zx */
    1026             : static GEN
    1027       21756 : Flx_ber_conj(GEN x, long j, long d, ulong p)
    1028             : {
    1029       21756 :   long i, deg = degpol(x);
    1030       21756 :   GEN y = const_vecsmall(d+1, 0), x2 = x+2, y2 = y+2;
    1031     1076565 :   for (i=0; i<=deg; i++) y2[Fl_mul(i, j, d)] = umodsu(x2[i], p);
    1032       21756 :   return Flx_renormalize(y, d+2);
    1033             : }
    1034             : 
    1035             : static GEN
    1036       26257 : ZX_ber_den(GEN Chi, long j, long d)
    1037             : {
    1038       26257 :   GEN x = pol_zero(d), x2 = x+2;
    1039       26257 :   if (Chi[2]>=0) gel(x2, Fl_neg(Fl_mul(Chi[2], j, d), d)) = gen_1;
    1040       26257 :   gel(x2, 0) = subiu(gel(x2, 0), 2);
    1041       26257 :   return ZX_renormalize(x, d+2);
    1042             : }
    1043             : 
    1044             : static GEN
    1045       14490 : Flx_ber_den(GEN Chi, long j, long d, ulong p)
    1046             : {
    1047       14490 :   GEN x = const_vecsmall(d+1, 0), x2 = x+2;
    1048       14490 :   if (Chi[2]>=0) x2[Fl_neg(Fl_mul(Chi[2], j, d), d)] = 1;
    1049       14490 :   x2[0] = Fl_sub(x2[0], 2, p);
    1050       14490 :   return Flx_renormalize(x, d+2);
    1051             : }
    1052             : 
    1053             : /* x is ZX of deg <= d-1 */
    1054             : static GEN
    1055         196 : ber_conj(GEN x, long k, long d)
    1056             : {
    1057         196 :   long i, deg = degpol(x);
    1058         196 :   GEN z = pol_zero(d);
    1059         196 :   if (k==1)
    1060           0 :     for (i=0; i<=deg; i++) gel(z, 2+i) = gel(x, 2+i);
    1061             :   else
    1062      122206 :     for (i=0; i<=deg; i++) gel(z, 2+Fl_mul(i, k, d)) = gel(x, 2+i);
    1063         196 :   return ZX_renormalize(z, d+2);
    1064             : }
    1065             : 
    1066             : /* The computation is fast when p^n and el=1+k*f*p^n are less than 2^64
    1067             :  *  for m <= n <= M
    1068             :  *  We believe M>=3 is enough when f%p=0 and M>=2 is enough for other case
    1069             :  *  because we expect that p^2 does not divide |A_{K,psi}| for a large p.
    1070             :  *  FIXME: M should be set according to p and f. */
    1071             : static void
    1072         196 : set_p_f(GEN pp, ulong f, long *pm, long *pM)
    1073             : {
    1074         196 :   ulong p = itou_or_0(pp);
    1075         196 :   if (!p || p >= 2000000) { *pm=2; *pM = dvdui(f, pp)? 3: 2; }
    1076         189 :   else if (p == 3)      { *pm=5; *pM=20; }
    1077         119 :   else if (p == 5)      { *pm=5; *pM=13; }
    1078          56 :   else if (p == 7)      { *pm=5; *pM=11; }
    1079          42 :   else if (p == 11)     { *pm=5; *pM=9; }
    1080          35 :   else if (p == 13)     { *pm=5; *pM=8; }
    1081          28 :   else if (p < 400)     { *pm=5; *pM=7; }
    1082           0 :   else if (p < 5000)    { *pm=3; *pM=5; }
    1083           0 :   else if (p < 50000)   { *pm=2; *pM=4; }
    1084           0 :   else                  { *pm=2; *pM=3; }
    1085         196 : }
    1086             : 
    1087             : static GEN
    1088       18795 : subgp2ary(GEN H, long n)
    1089             : {
    1090       18795 :   GEN v = gel(H, 3), w = cgetg(n+1, t_VECSMALL);
    1091       18795 :   long i, j, f = v[1];
    1092   399982464 :   for (i = 1, j = 0; i <= f; i++)
    1093   399963669 :     if (F2v_coeff(v,i)) w[++j] = i;
    1094       18795 :   return w;
    1095             : }
    1096             : 
    1097             : static GEN
    1098       19124 : Flv_FlvV_factorback(GEN g, GEN x, ulong q)
    1099       90216 : { pari_APPLY_ulong(Flv_factorback(g, gel(x,i), q)) }
    1100             : 
    1101             : /* lift chi character on G/H to character on G */
    1102             : static GEN
    1103       18795 : zncharlift(GEN chi, GEN ncycGH, GEN U, GEN cycG)
    1104             : {
    1105       18795 :   GEN nchi = char_normalize(chi, ncycGH);
    1106       18795 :   GEN c = ZV_ZM_mul(gel(nchi, 2), U), d = gel(nchi, 1);
    1107       18795 :   return char_denormalize(cycG, d, c);
    1108             : }
    1109             : 
    1110             : /* 0 <= c[i] < d, i=1..r; (c[1],...,c[r], d) = 1; find e[i] such that
    1111             :  * sum e[i]*c[i] = 1 mod d */
    1112             : static GEN
    1113       18795 : Flv_extgcd(GEN c, ulong d)
    1114             : {
    1115       18795 :   long i, j, u, f, l = lg(c);
    1116       18795 :   GEN e = zero_zv(l-1);
    1117       18795 :   if (l == 1) return e;
    1118       46053 :   for (f = d, i = 1; f != 1 && i < l; i++)
    1119             :   {
    1120       27258 :     f = cbezout(f, itou(gel(c,i)), &u, &e[i]);
    1121       27258 :     if (!e[i]) continue;
    1122       25004 :     e[i] = umodsu(e[i], d);
    1123       25004 :     u = umodsu(u, d);
    1124       32998 :     if (u != 1) for (j = 1; j < i; j++) e[j] = Fl_mul(e[j], u, d);
    1125             :   }
    1126       18795 :   return e;
    1127             : }
    1128             : 
    1129             : /* f!=p; return exact xi. */
    1130             : static GEN
    1131         462 : get_xi_1(GEN Chi, GEN Gam, long p, long f, long n, long d, ulong pm)
    1132             : {
    1133         462 :   GEN Gam1 = Gam+1, xi;
    1134             :   long i, j, j0, f0, pn, pn1, deg, pn1f;
    1135             : 
    1136         462 :   f0 = (f%p)?f:f/p;
    1137         462 :   pn = upowuu(p, n); pn1 = p*pn; pn1f = pn1%f;
    1138         462 :   xi = cgetg(pn+2, t_POL); xi[1] = evalsigne(1) | evalvarn(0);
    1139       17556 :   for (i=0; i<pn; i++) gel(xi, 2+i) = const_vecsmall(d+1, 0);
    1140      432754 :   for (j=1; j<pn1; j++)
    1141             :   {
    1142             :     long ipn1,*xij0;
    1143      432292 :     if ((j0 = Gam1[j])<0) continue;
    1144      415660 :     ipn1 = j%f; xij0 = gel(xi, 2+j0)+2;
    1145  4027401588 :     for (i=1; i<f0; i++)
    1146             :     {
    1147  4026985928 :       if ((ipn1 += pn1f) >= f) ipn1 -= f;
    1148  4026985928 :       if (ipn1==0 || (deg = Chi[ipn1])<0) continue;
    1149  2203029612 :       xij0[deg] += i;
    1150             :     }
    1151             :   }
    1152       17556 :   for (i=0; i<pn; i++) Flx_red_inplace(gel(xi, 2+i), pm);
    1153         462 :   return FlxX_renormalize(xi, pn+2);
    1154             : }
    1155             : 
    1156             : /* f=p; return p^(n+1)*xi mod pm. */
    1157             : static GEN
    1158           0 : get_xi_2(GEN Chi, GEN Gam, long p, long f, long n, long d, ulong pm)
    1159             : {
    1160             :   long a, amodf, i, j0, pn, pn1, deg;
    1161           0 :   GEN Gam1 = Gam+1, xi;
    1162             : 
    1163           0 :   pn = upowuu(p, n); pn1 = p*pn;
    1164           0 :   xi = cgetg(pn+2, t_POL); xi[1] = evalsigne(1) | evalvarn(0);
    1165           0 :   for (i=0; i<pn; i++) gel(xi, 2+i) = const_vecsmall(d+1, 0);
    1166           0 :   for (a=1,amodf=0; a<pn1; a++)  /* xi is exact */
    1167             :   {
    1168           0 :     if (++amodf==f) amodf = 0;
    1169           0 :     if ((j0=Gam1[a])<0 || amodf==0 || (deg=Chi[amodf])<0) continue;
    1170           0 :     mael(xi, 2+j0, 2+deg) += a;
    1171             :   }
    1172           0 :   for (i=0; i<pn; i++) Flx_red_inplace(gel(xi, 2+i), pm);
    1173           0 :   return FlxX_renormalize(xi, pn+2);
    1174             : }
    1175             : 
    1176             : static GEN
    1177          56 : pol_chi_xi(GEN K, long p, long j, long n)
    1178             : {
    1179          56 :   pari_sp av = avma;
    1180          56 :   GEN MinPol2 = gel(K, 7), xi = gel(K, 8);
    1181          56 :   long d = K_get_d(K), f = K_get_f(K), d_chi = K_get_dchi(K);
    1182          56 :   long wd, minpolpow = (f==p)?2*n+1:n, pm = upowuu(p, minpolpow);
    1183             :   GEN ser, pol, xi_conj;
    1184             :   pari_timer ti;
    1185             : 
    1186             :   /* xi is FlxX mod p^m, MinPol2 is Flx mod p^m, xi_conj is FlxqX. */
    1187          56 :   xi_conj = FlxqX_xi_conj(xi, MinPol2, j, d, pm);
    1188          56 :   if (d_chi==1)  /* d_chi==1 if f==p */
    1189             :   {
    1190          56 :     xi_conj = FlxX_to_Flx(xi_conj);
    1191          56 :     if (f==p) xi_conj = zx_z_divexact(xi_conj, upowuu(p, n+1));
    1192             :   }
    1193             :   /* Now xi_conj is mod p^n */
    1194          56 :   if (DEBUGLEVEL>1) timer_start(&ti);
    1195          56 :   ser = (d_chi==1) ? Flxn_translate1(xi_conj, p, n)
    1196          56 :     : FlxXn_translate1(xi_conj, p, n);
    1197          56 :   if (DEBUGLEVEL>1) timer_printf(&ti, "Flx%sn_translate1",(d_chi==1)?"":"X");
    1198          56 :   wd = (d_chi==1)?Flx_weier_deg(ser, p):FlxX_weier_deg(ser, p);
    1199          56 :   if (wd<0) pari_err(e_MISC,"pol_chi_xi: precision too low. Increase n!\n");
    1200          56 :   else if (wd==0) return pol_1(0);
    1201             :   /* wd>0, convert to dist. poly. */
    1202          56 :   if (d_chi>1)  /* f!=p. minpolpow==n */
    1203             :   {
    1204           0 :     ser = FlxqX_xi_norm(ser, MinPol2, p, d, upowuu(p, n));
    1205           0 :     ser = FlxX_to_Flx(ser);
    1206             :   }
    1207          56 :   pol = Flx_to_ZX(Flxn_Weierstrass_prep(ser, p, n, d_chi));
    1208          56 :   setvarn(pol, fetch_user_var("T"));
    1209             : #ifdef DEBUG
    1210             :   if (wd>0 && d_chi>1)
    1211             :     err_printf("(wd,d_chi,p,f,d,j,H)=(%ld,%ld,%ld,%ld,%ld,%ld,%Ps)\n",
    1212             :         wd,d_chi,p,f,d,j,gmael3(K, 1, 1, 1));
    1213             : #endif
    1214          56 :   return gerepilecopy(av, pol);
    1215             : }
    1216             : 
    1217             : /* return 0 if lam_psi (psi=chi^j) is determined to be zero.
    1218             :  * else return -1.
    1219             :  * If psi(p)!=1, then N_{Q(zeta_d)/Q}(1-psi(p))!=0 (mod p) */
    1220             : static long
    1221       14504 : lam_chi_ber(GEN K, long p, long j)
    1222             : {
    1223       14504 :   pari_sp av = avma;
    1224       14504 :   GEN B1, B2, Chi = gel(K, 2), MinPol2 = gel(K, 7), B_num = gel(K, 8);
    1225       14504 :   long x, p2 = p*p, d = K_get_d(K), f = K_get_f(K);
    1226             : 
    1227       14504 :   if (f == d+1 && p == f && j == 1) return 0;  /* Teichmuller */
    1228             : 
    1229       14490 :   B1 = Flx_rem(Flx_ber_conj(B_num, j, d, p2), MinPol2, p2);
    1230       14490 :   B2 = Flx_rem(Flx_ber_den(Chi, j, d, p2), MinPol2, p2);
    1231       14490 :   if (degpol(B1)<0 || degpol(B2)<0)
    1232          49 :     return gc_long(av, -1); /* 0 mod p^2 */
    1233       14441 :   x = zx_lval(B1, p) - zx_lval(B2, p);
    1234       14441 :   if (x<0) pari_err_BUG("subcycloiwasawa [Bernoulli number]");
    1235       14441 :   return gc_long(av, x==0 ? 0: -1);
    1236             : }
    1237             : 
    1238             : static long
    1239         910 : lam_chi_xi(GEN K, long p, long j, long n)
    1240             : {
    1241         910 :   pari_sp av = avma;
    1242         910 :   GEN xi_conj, z, MinPol2 = gel(K, 7), xi = gel(K, 8);
    1243         910 :   long d = K_get_d(K), f = K_get_f(K), d_chi = K_get_dchi(K);
    1244         910 :   long wd, minpolpow = (f==p)?n+2:1, pm = upowuu(p, minpolpow);
    1245             : 
    1246             :   /* xi is FlxX mod p^m, MinPol2 is Flx mod p^m, xi_conj is FlxqX. */
    1247         910 :   xi_conj = FlxqX_xi_conj(xi, MinPol2, j, d, pm);
    1248         910 :   if (d_chi==1)  /* d_chi==1 if f==p */
    1249             :   {
    1250         686 :     xi_conj = FlxX_to_Flx(xi_conj);
    1251         686 :     if (f==p) xi_conj = zx_z_divexact(xi_conj, upowuu(p, n+1));
    1252             :   }
    1253             :   /* Now xi_conj is mod p^n */
    1254         686 :   z = (d_chi==1) ? Flxn_translate1(xi_conj, p, n)
    1255         910 :     : FlxXn_translate1(xi_conj, p, n);
    1256         910 :   wd = (d_chi==1)?Flx_weier_deg(z, p):FlxX_weier_deg(z, p);
    1257             : #ifdef DEBUG
    1258             :   if (wd>0 && d_chi>1)
    1259             :     err_printf("(wd,d_chi,p,f,d,j,H)=(%ld,%ld,%ld,%ld,%ld,%ld,%Ps)\n",
    1260             :         wd,d_chi,p,f,d,j,gmael3(K, 1, 1, 1));
    1261             : #endif
    1262         910 :   return gc_long(av, wd<0 ? -1 : wd*d_chi);
    1263             : }
    1264             : 
    1265             : /* K = [H1, Chi, Minpol, C, [d_chi, n_conj]] */
    1266             : static GEN
    1267          56 : imag_cyc_pol(GEN K, long p, long n)
    1268             : {
    1269          56 :   pari_sp av = avma;
    1270          56 :   GEN Chi = gel(K, 2), MinPol = gel(K, 3), C = gel(K, 4), MinPol2;
    1271          56 :   long d_K = K_get_d(K), f = K_get_f(K), n_conj = K_get_nconj(K);
    1272          56 :   long i, q0, pn1, pM, pmodf = p%f, n_done = 0;
    1273          56 :   GEN z = nullvec(), Gam, xi, Lam, K2;
    1274             : 
    1275          56 :   Lam = const_vecsmall(n_conj, -1);
    1276          56 :   if (pmodf==0 || Chi[pmodf]) /* mark trivial chi-part using Bernoulli number */
    1277             :   {
    1278          14 :     MinPol2 = ZX_to_Flx(MinPol, p*p);  /* p^2 for B_{1,chi} */
    1279          14 :     K2 = shallowconcat(K, mkvec2(MinPol2, zx_ber_num(Chi, f, d_K)));
    1280          42 :     for (i=1; i<=n_conj; i++)
    1281          28 :       if ((Lam[i] = lam_chi_ber(K2, p, C[i])) == 0) n_done++;
    1282          14 :     if (n_conj==n_done) return gerepilecopy(av, z); /* all chi-parts trivial */
    1283             :   }
    1284          49 :   q0 = (f%p)? f*p: f;
    1285          49 :   pn1 = upowuu(p, n+1);
    1286          49 :   Gam = set_gam((1+q0)%pn1, p, n);
    1287          49 :   pM = upowuu(p, (f==p)? 2*n+1: n);
    1288          49 :   MinPol2 = ZX_to_Flx(MinPol, pM);
    1289           0 :   xi = (f==p)? get_xi_2(Chi, Gam, p, f, n, d_K, pM)
    1290          49 :              : get_xi_1(Chi, Gam, p, f, n, d_K, pM);
    1291          49 :   K2 = shallowconcat(K, mkvec2(MinPol2, xi));
    1292         105 :   for (i=1; i<=n_conj; i++)
    1293             :   {
    1294             :     GEN z1;
    1295          56 :     if (Lam[i]>=0) continue;
    1296          56 :     z1 = pol_chi_xi(K2, p, C[i], n);
    1297          56 :     if (degpol(z1)) z = vec_append(z, z1);  /* degpol(z1) may be zero */
    1298             :   }
    1299          49 :   return gerepilecopy(av, z);
    1300             : }
    1301             : 
    1302             : /* K is an imaginary cyclic extension of degree d contained in Q(zeta_f)
    1303             :  * H is the subgr of G=(Z/fZ)^* corresponding to K
    1304             :  * h=|H|, d*h=phi(f)
    1305             :  * G/H=<g> i.e. g^d \in H
    1306             :  * d_chi=[Q_p(zeta_d):Q_p], i.e. p^d_chi=1 (mod d)
    1307             :  * An inj. char. of G(K/Q) is automatically imaginary.
    1308             :  *
    1309             :  * G(K/Q)=G/H=<g>, chi:G(K/Q) -> overline{Q_p} s.t. chi(g)=zeta_d^(-1)
    1310             :  * Chi[a]=k, k<0  => chi(a)=0
    1311             :  *           k>=0 => chi(a)=zeta_d^(-k)
    1312             :  * psi=chi^j, j in C : repre. of inj. odd char.
    1313             :  * psi(p)==1 <=> chi(p)^j==0 <=> j*Chi[p]=0 (mod d) <=> Chi[p]==0
    1314             :  *
    1315             :  * K = [H1, Chi, Minpol, C, [d_chi, n_conj]] */
    1316             : static long
    1317        3262 : imag_cyc_lam(GEN K, long p)
    1318             : {
    1319        3262 :   pari_sp av = avma;
    1320        3262 :   GEN Chi = gel(K, 2), MinPol = gel(K, 3), C = gel(K, 4), MinPol2;
    1321        3262 :   long d_K = K_get_d(K), f = K_get_f(K), n_conj = K_get_nconj(K);
    1322        3262 :   long i, q0, n, pmodf = p%f, n_done = 0;
    1323             :   ulong pn1, pM;
    1324        3262 :   GEN p0 = utoi(p), Gam, Lam, xi, K2;
    1325             : 
    1326        3262 :   q0 = (f%p)? f*p: f;
    1327        3262 :   Lam = const_vecsmall(n_conj, -1);
    1328        3262 :   if (pmodf==0 || Chi[pmodf])  /* 1st trial is Bernoulli number */
    1329             :   {
    1330        3052 :     MinPol2 = ZX_to_Flx(MinPol, p*p);  /* p^2 for B_{1,chi} */
    1331        3052 :     K2 = shallowconcat(K, mkvec2(MinPol2, zx_ber_num(Chi, f, d_K)));
    1332       17528 :     for (i=1; i<=n_conj; i++)
    1333       14476 :       if ((Lam[i] = lam_chi_ber(K2, p, C[i])) == 0) n_done++;
    1334        3052 :     if (n_conj==n_done) return gc_long(av, 0);  /* all chi-parts trivial */
    1335             :   }
    1336         413 :   pM = pn1 = p;
    1337         413 :   for (n=1; n>=0; n++)  /* 2nd trial is Stickelberger element */
    1338             :   {
    1339         413 :     pn1 *= p; /* p^(n+1) */
    1340         413 :     if (f == p)
    1341             :     { /* do not use set_minpol: it returns a new pol for each call */
    1342           0 :       GEN fac, cofac, v, pol = polcyclo(d_K, 0);
    1343           0 :       pM = pn1 * p; /* p^(n+2) */
    1344           0 :       fac = FpX_red(MinPol, p0); cofac = FpX_div(pol, fac, p0);
    1345           0 :       v = ZpX_liftfact(pol, mkvec2(fac, cofac), utoipos(pM), p0, n+2);
    1346           0 :       MinPol2 = gel(v, 1);
    1347             :     }
    1348         413 :     Gam = set_gam((1+q0)%pn1, p, n);
    1349         413 :     MinPol2 = ZX_to_Flx(MinPol, pM);
    1350           0 :     xi = (f==p)? get_xi_2(Chi, Gam, p, f, n, d_K, pM)
    1351         413 :                : get_xi_1(Chi, Gam, p, f, n, d_K, pM);
    1352         413 :     K2 = shallowconcat(K, mkvec2(MinPol2, xi));
    1353        2205 :     for (i=1; i<=n_conj; i++)
    1354        1792 :       if (Lam[i]<0 && (Lam[i] = lam_chi_xi(K2, p, C[i], n)) >= 0) n_done++;
    1355         413 :     if (n_conj==n_done) break;
    1356             :   }
    1357         413 :   return gc_long(av, zv_sum(Lam));
    1358             : }
    1359             : static GEN
    1360         329 : GHinit(long f, GEN HH, GEN *pcycGH)
    1361             : {
    1362         329 :   GEN G = znstar0(utoipos(f), 1);
    1363         329 :   GEN U, Ui, cycG, cycGH, ncycGH, gG, gGH, vChar, vH1, P, gH = gel(HH, 1);
    1364         329 :   long i, expG, n_f, lgH = lg(gH); /* gens. of H */
    1365         329 :   P = cgetg(lgH, t_MAT);
    1366         805 :   for (i = 1; i < lgH; i++) gel(P,i) = Zideallog(G, utoi(gH[i]));
    1367             : 
    1368             :   /* group structure of G/H */
    1369         329 :   cycG = znstar_get_cyc(G);
    1370         329 :   expG = itou(gel(cycG, 1));
    1371             :   /* gG generators of G, gGH generators of G/H: gGH = g.Ui, g = gGH.U */
    1372         329 :   cycGH = ZM_snf_group(hnfmodid(P, cycG), &U, &Ui);
    1373         329 :   ncycGH = cyc_normalize(cycGH);
    1374         329 :   gG = ZV_to_Flv(znstar_get_gen(G), f); /* gens. of G */
    1375             :   /* generators of G/H */
    1376         329 :   gGH = Flv_FlvV_factorback(gG, ZM_to_Flm(Ui, expG), f);
    1377         329 :   vChar = chargalois(cycGH, NULL); n_f = lg(vChar)-2;
    1378         329 :   vH1 = cgetg(n_f+1, t_VEC);
    1379       19124 :   for (i = 1; i <= n_f; i++)
    1380             :   { /* skip trivial character */
    1381       18795 :     GEN chi = gel(vChar,i+1), nchi = char_normalize(chi, ncycGH);
    1382       18795 :     GEN chiG, E, H1, C = gel(nchi, 2);
    1383       18795 :     long e, he, gen, d = itou(gel(nchi, 1));
    1384             :     /* chi(prod g[i]^e[i]) = e(sum e[i]*C[i] / d), chi has order d = #(G/H1)*/
    1385       18795 :     E = Flv_extgcd(C, d); /* \sum C[i]*E[i] = 1 in Z/dZ */
    1386             : 
    1387       18795 :     chiG = zncharlift(chi, ncycGH, U, cycG);
    1388       18795 :     H1 =  charker(cycG, chiG); /* H1 < G with G/H1 cyclic */
    1389       18795 :     e = itou( zncharconductor(G, chiG) ); /* cond H1 = cond chi */
    1390       18795 :     H1 = Flv_FlvV_factorback(zv_to_Flv(gG, e), ZM_to_Flm(H1, expG), e);
    1391       18795 :     gen = Flv_factorback(zv_to_Flv(gGH, e), E, e);
    1392       18795 :     H1 = znstar_generate(e, H1); /* G/H1 = <gen>, chi(gen) = e(1/d) */
    1393       18795 :     he = eulerphiu(e) / d;
    1394             :     /* G/H1 = <gen> cyclic of index d, e = cond(H1) */
    1395       18795 :     gel(vH1,i) = mkvec3(H1, mkvecsmall5(d,e,he,srh_1(H1), gen),
    1396             :                         subgp2ary(H1, he));
    1397             :   }
    1398         329 :   if (pcycGH) *pcycGH = cycGH;
    1399         329 :   return vH1;
    1400             : }
    1401             : 
    1402             : /* aH=g^iH ==> chi(a)=zeta_n^(-i); Chi[g^iH]=i; Chi[0] is never accessed */
    1403             : static GEN
    1404       13419 : get_chi(GEN H1)
    1405             : {
    1406       13419 :   GEN H = _get_H(H1);
    1407       13419 :   long i, j, gi, d = _get_d(H1), f = _get_f(H1), h = _get_h(H1), g = _get_g(H1);
    1408       13419 :   GEN Chi = const_vecsmall(f-1, -1);
    1409             : 
    1410     5584159 :   for (j=1; j<=h; j++) Chi[H[j]] = 0; /* i = 0 */
    1411      396592 :   for (i = 1, gi = g; i < d; i++)
    1412             :   {
    1413    40492081 :     for (j=1; j<=h; j++) Chi[Fl_mul(gi, H[j], f)] = i;
    1414      383173 :     gi = Fl_mul(gi, g, f);
    1415             :   }
    1416       13419 :   return Chi;
    1417             : }
    1418             : 
    1419             : static void
    1420          14 : errpdiv(const char *f, GEN p, long d)
    1421             : {
    1422          14 :   pari_err_DOMAIN(f, "p", "divides",
    1423          14 :                   strtoGENstr(stack_sprintf("[F:Q] = %ld", d)), p);
    1424           0 : }
    1425             : /* p odd doesn't divide degF; return lambda invariant if n==0 and
    1426             :  * iwasawa polynomials if n>=1 */
    1427             : static GEN
    1428          49 : abeliwasawa(long p, long f, GEN HH, long degF, long n)
    1429             : {
    1430          49 :   long lam = 0, i, n_f;
    1431          49 :   GEN vH1, vData, z = nullvec(), p0 = utoi(p) ;
    1432             : 
    1433          49 :   vH1 = GHinit(f, HH, NULL); n_f = lg(vH1)-1;
    1434          49 :   vData = const_vec(degF, NULL);
    1435        6608 :   for (i=1; i<=n_f; i++) /* prescan. set Teichmuller */
    1436             :   {
    1437        6573 :     GEN H1 = gel(vH1,i);
    1438        6573 :     long d_K = _get_d(H1), f_K = _get_f(H1), g_K = _get_g(H1);
    1439             : 
    1440        6573 :     if (f_K == d_K+1 && p == f_K)  /* found K=Q(zeta_p) */
    1441             :     {
    1442          14 :       long d_chi = 1, n_conj = eulerphiu(d_K);
    1443          14 :       GEN C = set_C(p, d_K, d_chi, n_conj);
    1444          14 :       long minpow = n? 2*n+1: 2;
    1445          14 :       GEN MinPol = set_minpol_teich(g_K, p0, minpow);
    1446          14 :       gel(vData, d_K) = mkvec4(MinPol, C, NULL, mkvecsmall2(d_chi, n_conj));
    1447          14 :       break;
    1448             :     }
    1449             :   }
    1450             : 
    1451        6664 :   for (i=1; i<=n_f; i++)
    1452             :   {
    1453        6615 :     GEN H1 = gel(vH1,i), z1, Chi, K;
    1454        6615 :     long d_K = _get_d(H1), s = _get_s(H1);
    1455             : 
    1456        6615 :     if (s) continue;  /* F is real */
    1457             : #ifdef DEBUG
    1458             :     err_printf("  handling %s cyclic subfield K, deg(K)=%ld, cond(K)=%ld\n",
    1459             :         s? "a real": "an imaginary", d_K, _get_f(H1));
    1460             : #endif
    1461        3318 :     if (!gel(vData, d_K))
    1462             :     {
    1463         126 :       long d_chi = order_f_x(d_K, p), n_conj = eulerphiu(d_K)/d_chi;
    1464         126 :       GEN C = set_C(p, d_K, d_chi, n_conj);
    1465         126 :       long minpow = n? n+1: 2;
    1466         126 :       GEN MinPol = set_minpol(d_K, p0, minpow, n_conj);
    1467         126 :       gel(vData, d_K) = mkvec4(MinPol, C, NULL, mkvecsmall2(d_chi, n_conj));
    1468             :     }
    1469        3318 :     Chi = get_chi(H1);
    1470        3318 :     K = shallowconcat(mkvec2(H1, Chi), gel(vData, d_K));
    1471        3318 :     if (n==0) lam += imag_cyc_lam(K, p);
    1472          56 :     else if (lg(z1 = imag_cyc_pol(K, p, n)) > 1) z = shallowconcat(z, z1);
    1473             :   }
    1474          49 :   return n? z: mkvecs(lam);
    1475             : }
    1476             : 
    1477             : static GEN
    1478          77 : ary2mat(GEN x, long n)
    1479             : {
    1480             :   long i, j;
    1481          77 :   GEN z = cgetg(n+1,t_MAT);
    1482         182 :   for (i=1; i<=n; i++)
    1483             :   {
    1484         105 :     gel(z,i) = cgetg(n+1,t_COL);
    1485         280 :     for (j=1; j<=n; j++) gmael(z, i, j) = utoi(x[(i-1)*n+j-1]);
    1486             :   }
    1487          77 :   return z;
    1488             : }
    1489             : 
    1490             : static long
    1491           0 : is_cyclic(GEN x)
    1492             : {
    1493           0 :   GEN y = gel(x, 2);
    1494           0 :   long i, l = lg(y), n = 0;
    1495           0 :   for (i = 1; i < l; i++) if (signe(gel(y,i))) n++;
    1496           0 :   return n <= 1;
    1497             : }
    1498             : 
    1499             : static GEN
    1500          77 : make_p_part(GEN y, ulong p, long d_pow)
    1501             : {
    1502          77 :   long i, l = lg(y);
    1503          77 :   GEN z = cgetg(l, t_VECSMALL);
    1504         182 :   for (i = 1; i < l; i++) z[i] = signe(gel(y,i))? Z_lval(gel(y,i), p): d_pow;
    1505          77 :   return z;
    1506             : }
    1507             : 
    1508             : static GEN
    1509          77 : structure_MLL(GEN y, long d_pow)
    1510             : {
    1511          77 :   long y0, i, l = lg(y);
    1512          77 :   GEN x = gen_0, E = cgetg(l, t_VEC);
    1513         182 :   for (i = 1; i < l; i++)
    1514             :   {
    1515         105 :     if ((y0 = d_pow-y[i]) < 0) y0 = 0;
    1516         105 :     x = addiu(x, y0);
    1517         105 :     gel(E, l-i) = utoi(y0);
    1518             :   }
    1519          77 :   return mkvec2(x, E);
    1520             : }
    1521             : 
    1522             : static long
    1523          14 : find_del_el(GEN *oldgr, GEN newgr, long n, long n_el, long d_chi)
    1524             : {
    1525          14 :   if (n_el==1) return 1;
    1526          14 :   if (equalis(gmael(newgr, 2, 1), n_el)) return n;
    1527          14 :   if (cmpii(gel(*oldgr, 1), gel(newgr, 1)) >= 0) return n;
    1528          14 :   if (n > 1 && is_cyclic(newgr)) { *oldgr = newgr; return n-1; }
    1529          14 :   if (n == n_el) return n;
    1530          14 :   if (cmpis(gel(newgr, 1), n*d_chi) < 0) return n;
    1531          14 :   return 0;
    1532             : }
    1533             : 
    1534             : static GEN
    1535           7 : subgr2vecsmall(GEN H, long h, long f)
    1536             : {
    1537             :   long i;
    1538           7 :   GEN z = const_vecsmall(f-1, 0); /* f=lg(z) */
    1539        2023 :   for (i=1; i<=h; i++) z[H[i]] = 1; /* H[i]!=0 */
    1540           7 :   return z;
    1541             : }
    1542             : 
    1543             : /* K is the subfield of Q(zeta_f) with degree d corresponding to the subgroup
    1544             :  * H in (Z/fZ)^*; for a divisor e of f, zeta_e \in K <=> H \subset He. */
    1545             : static long
    1546         119 : root_of_1(long f, GEN H)
    1547             : {
    1548         119 :   GEN g = gel(H, 1); /* generators */
    1549         119 :   long e = f, i, l = lg(g);
    1550         119 :   for (i = 1; i < l; i++)
    1551             :   {
    1552          98 :     e = cgcd(e, g[i] - 1);
    1553          98 :     if (e <= 2) return 2;
    1554             :   }
    1555          21 :   return odd(e)? (e<<1): e;
    1556             : }
    1557             : 
    1558             : static long
    1559         259 : find_ele(GEN H)
    1560             : {
    1561         259 :   long i, f=lg(H);
    1562      369852 :   for (i=1; i<f; i++) if (H[i]) return i;
    1563           7 :   return 0;
    1564             : }
    1565             : 
    1566             : static void
    1567         252 : delete_ele(GEN H, long j, long el)
    1568             : {
    1569         252 :   long f = lg(H), x = 1;
    1570        2016 :   do H[Fl_mul(j,x,f)] = 0;
    1571        2016 :   while ((x=Fl_mul(x,el,f))!=1);
    1572         252 : }
    1573             : 
    1574             : static GEN
    1575           7 : get_coset(GEN H, long h, long f, long el)
    1576             : {
    1577           7 :   long i, j, k = h/order_f_x(f, el);
    1578           7 :   GEN H2, coset = const_vecsmall(k, 0);
    1579           7 :   H2 = subgr2vecsmall(H, h, f);
    1580         259 :   for (i=0; (j=find_ele(H2))>0; i++)
    1581             :   {
    1582         252 :     coset[1+i] = j;
    1583         252 :     delete_ele(H2, j, el);
    1584             :   }
    1585           7 :   if (i != k) pari_err_BUG("failed to find coset\n");
    1586           7 :   return coset;
    1587             : }
    1588             : 
    1589             : static long
    1590        3024 : srh_pol(GEN xpows, GEN vn, GEN pols, long el, long k, long f)
    1591             : {
    1592        3024 :   pari_sp av = avma;
    1593        3024 :   long i, j, l = lg(pols), d = degpol(gel(pols, 1));
    1594        3024 :   GEN pol = gel(pols, 1);
    1595             : 
    1596      654696 :   for (i=1; i<l; i++)
    1597             :   {
    1598             :     GEN x, y, z;
    1599      654696 :     if (vn[i]==0) continue;
    1600      331170 :     y = gel(pols, vn[i]);
    1601      331170 :     z = pol0_Flx(0);
    1602     3311700 :     for (j=0; j<=d; j++)
    1603     2980530 :       z = Flx_add(z, Flx_Fl_mul(gel(xpows, 1+Fl_mul(j, k, f)), y[2+j], el), el);
    1604      331170 :     x = Flx_rem(z, pol, el);
    1605      331170 :     if (lg(x)==2)
    1606        3024 :     {vn[i] = 0; return gc_long(av, i); }  /* pols[i] is min pol of zeta_f^k */
    1607             :   }
    1608           0 :   pari_err_BUG("subcyclopclgp [minimal polinomial]");
    1609           0 :   return 0;  /* to suppress warning */
    1610             : }
    1611             : 
    1612             : /* e_chi[i mod dK] = chi(i*j), i = 0..dK-1; beware: e_chi is translated ! */
    1613             : static GEN
    1614       35857 : get_e_chi(GEN K, ulong j, ulong d, ulong *pdK)
    1615             : {
    1616       35857 :   ulong i, dK = K_get_d(K);
    1617       35857 :   GEN TR = gel(K,4) + 2, e_chi = cgetg(dK+1, t_VECSMALL) + 1;
    1618       35857 :   if (j == 1)
    1619      286902 :     for (i = 0; i < dK; i++) e_chi[i] = umodiu(gel(TR, i), d);
    1620             :   else
    1621      983345 :     for (i = 0; i < dK; i++) e_chi[i] = umodiu(gel(TR, Fl_mul(i, j, dK)), d);
    1622       35857 :   *pdK = dK; return e_chi;
    1623             : }
    1624             : static GEN
    1625        5145 : get_e_chi_ll(GEN K, ulong j, GEN d)
    1626             : {
    1627        5145 :   ulong i, dK = umael3(K, 1, 2, 1);
    1628        5145 :   GEN TR = gel(K,4) + 2, e_chi = cgetg(dK+1, t_VEC) + 1;
    1629      246785 :   for (i = 0; i < dK; i++) gel(e_chi,i) = modii(gel(TR, Fl_mul(i, j, dK)), d);
    1630        5145 :   return e_chi;
    1631             : }
    1632             : 
    1633             : /* el=1 (mod f) */
    1634             : static long
    1635           0 : chk_el_real_f(GEN K, ulong p, ulong d_pow, ulong el, ulong j0)
    1636             : {
    1637           0 :   pari_sp av = avma;
    1638           0 :   GEN H = K_get_H(K);
    1639           0 :   ulong d_K, f = K_get_f(K), h = K_get_h(K), g_K = K_get_g(K);
    1640           0 :   ulong  i, j, gi, d = upowuu(p, d_pow), dp = d*p;
    1641           0 :   ulong g_el, z_f, flag = 0, el1f = (el-1)/f, el1dp = (el-1)/dp;
    1642           0 :   GEN e_chi = get_e_chi(K, j0, dp, &d_K);
    1643           0 :   GEN vz_f, xi_el = cgetg(d_K+1, t_VECSMALL)+1;
    1644             : 
    1645           0 :   g_el = pgener_Fl(el);
    1646           0 :   z_f = Fl_powu(g_el, el1f, el);
    1647           0 :   vz_f = Fl_powers(z_f, f-1, el)+1;
    1648             : 
    1649           0 :   for (gi = 1, i = 0; i < d_K; i++)
    1650             :   {
    1651           0 :     ulong x = 1;
    1652           0 :     for (j = 1; j <= h; j++)
    1653             :     {
    1654           0 :       ulong y = Fl_mul(H[j], gi, f);
    1655           0 :       x = Fl_mul(x, vz_f[y]-1, el); /* vz_f[y] = z_f^y  */
    1656             :     }
    1657           0 :     gi = Fl_mul(gi, g_K, f);
    1658           0 :     xi_el[i] = x;  /* xi_el[i] = xi^{g_K^i} mod el */
    1659             :   }
    1660           0 :   for (i=0; i<d_K; i++)
    1661             :   {
    1662           0 :     ulong x = 1;
    1663           0 :     for (j=0; j<d_K; j++)
    1664           0 :       x = Fl_mul(x, Fl_powu(xi_el[j], e_chi[(i+j)%d_K], el), el);
    1665           0 :     if ((x = Fl_powu(x, el1dp, el))!=1) flag = 1;
    1666           0 :     if (Fl_powu(x, p, el)!=1) return gc_long(av,0);
    1667             :   }
    1668           0 :   return gc_long(av, flag?1:0);
    1669             : }
    1670             : 
    1671             : /* For a cyclic field K contained in Q(zeta_f),
    1672             :  * computes minimal polynomial T of theta=Tr_{Q(zeta_f)/K}(zeta_f) over Q
    1673             :  * and conjugates of theta */
    1674             : static GEN
    1675          14 : minpol_theta(GEN K)
    1676             : {
    1677          14 :   GEN HH = gmael3(K,1,1,1);
    1678          14 :   return galoissubcyclo(utoi(K_get_f(K)), HH, 0, 0);
    1679             : }
    1680             : 
    1681             : /*  xi[1+i] = i-th conj of xi = Tr_{Q(zeta_f)/K}(1-zeta_f).
    1682             :  * |1-(cos(x)+i*sin(x))|^2 = 2(1-cos(x)) */
    1683             : static GEN
    1684           0 : xi_approx(GEN K, long prec)
    1685             : {
    1686           0 :   pari_sp av = avma;
    1687           0 :   GEN H = K_get_H(K);
    1688           0 :   long d_K = K_get_d(K), f = K_get_f(K), h = K_get_h(K), g_K = K_get_g(K);
    1689           0 :   GEN xi = cgetg(d_K+1, t_COL), vz_f = grootsof1(f, prec);
    1690           0 :   long i, j, g = 1, h2 = h>>1;
    1691           0 :   for (i=1; i<=d_K; i++)
    1692             :   {
    1693           0 :     GEN y = real_1(prec);
    1694           0 :     for (j=1; j<=h2; j++)
    1695             :     {
    1696           0 :       GEN z = gmael(vz_f, 1+Fl_mul(H[j], g, f), 1);
    1697           0 :       y = mulrr(y, shiftr(subsr(1, z), 1));
    1698             :     }
    1699           0 :     gel(xi, i) = y;
    1700           0 :     g = Fl_mul(g, g_K, f);
    1701             :   }
    1702           0 :   return gerepilecopy(av, xi);
    1703             : }
    1704             : 
    1705             : static GEN
    1706          47 : theta_xi_el(GEN K, ulong el)
    1707             : {
    1708          47 :   pari_sp av = avma;
    1709          47 :   GEN H = K_get_H(K);
    1710          47 :   ulong d_K = K_get_d(K), f = K_get_f(K), h = K_get_h(K), g_K = K_get_g(K);
    1711          47 :   GEN theta = cgetg(d_K+1, t_VECSMALL), xi = cgetg(d_K+1, t_VECSMALL), vz_f;
    1712          47 :   ulong i, j, g = 1, x, y, g_el, z_f;
    1713             : 
    1714          47 :   g_el = pgener_Fl(el);
    1715          47 :   z_f = Fl_powu(g_el, (el-1)/f, el);
    1716          47 :   vz_f = Fl_powers(z_f, f-1, el);
    1717        1109 :   for (i=1; i<=d_K; i++)
    1718             :   {
    1719        1062 :     x = 0; y = 1;
    1720       28254 :     for (j=1; j<=h; j++)
    1721             :     {
    1722       27192 :       ulong z = vz_f[1+Fl_mul(H[j], g, f)];
    1723       27192 :       x = Fl_add(x, z, el);
    1724       27192 :       y = Fl_mul(y, z-1, el);
    1725             :     }
    1726        1062 :     theta[i] = x;
    1727        1062 :     xi[i] = y;
    1728        1062 :     g = Fl_mul(g, g_K, f);
    1729             :   }
    1730          47 :   return gerepilecopy(av, mkvec2(theta, xi));
    1731             : }
    1732             : 
    1733             : static GEN
    1734          47 : make_Xi(GEN xi, long d)
    1735             : {
    1736             :   long i, j;
    1737          47 :   GEN Xi = cgetg(d+1, t_MAT);
    1738        1109 :   for (j=0; j<d; j++)
    1739             :   {
    1740        1062 :     GEN x = cgetg(d+1, t_VECSMALL);
    1741        1062 :     gel(Xi, 1+j) = x;
    1742       27714 :     for (i=0; i<d; i++) x[1+i] = xi[1+(i+j)%d];
    1743             :   }
    1744          47 :   return Xi;
    1745             : }
    1746             : 
    1747             : static GEN
    1748          47 : make_Theta(GEN theta, ulong d, ulong el)
    1749             : {
    1750             :   ulong i;
    1751          47 :   GEN Theta = cgetg(d+1, t_MAT);
    1752        1109 :   for (i=1; i<=d; i++) gel(Theta, i) = Fl_powers(theta[i], d-1, el);
    1753          47 :   return Flm_inv(Theta, el);
    1754             : }
    1755             : 
    1756             : static GEN
    1757          47 : Xi_el(GEN K, GEN tInvA, ulong el)
    1758             : {
    1759          47 :   pari_sp av = avma;
    1760          47 :   ulong d_K = K_get_d(K);
    1761          47 :   GEN tx = theta_xi_el(K, el), Theta, Xi, X;
    1762             : 
    1763          47 :   if ((Theta = make_Theta(gel(tx, 1), d_K, el))==NULL) return NULL;
    1764          47 :   Xi = make_Xi(gel(tx, 2), d_K);
    1765          47 :   X = Flm_mul(Flm_mul(Xi, Theta, el), ZM_to_Flm(tInvA, el), el);
    1766          47 :   return gerepilecopy(av, X);
    1767             : }
    1768             : 
    1769             : static GEN
    1770           0 : pol_xi_el(GEN K, ulong el)
    1771             : {
    1772           0 :   pari_sp av = avma;
    1773           0 :   ulong d_K = K_get_d(K), f = K_get_f(K), h = K_get_h(K), g_K = K_get_g(K);
    1774           0 :   GEN H = K_get_H(K), xi = cgetg(d_K+1, t_VECSMALL), vz_f;
    1775           0 :   ulong i, j, g = 1, y, g_el, z_f;
    1776             : 
    1777           0 :   g_el = pgener_Fl(el);
    1778           0 :   z_f = Fl_powu(g_el, (el-1)/f, el);
    1779           0 :   vz_f = Fl_powers(z_f, f-1, el);
    1780           0 :   for (i=1; i<=d_K; i++)
    1781             :   {
    1782           0 :     y = 1;
    1783           0 :     for (j=1; j<=h; j++)
    1784             :     {
    1785           0 :       ulong z = vz_f[1+Fl_mul(H[j], g, f)];
    1786           0 :       y = Fl_mul(y, z-1, el);
    1787             :     }
    1788           0 :     xi[i] = y;
    1789           0 :     g = Fl_mul(g, g_K, f);
    1790             :   }
    1791           0 :   return gerepilecopy(av, Flv_roots_to_pol(xi, el, 0));
    1792             : }
    1793             : 
    1794             : /* theta[1+i] = i-th conj of theta; xi[1+i] = i-th conj of xi. */
    1795             : static GEN
    1796          14 : theta_xi_approx(GEN K, long prec)
    1797             : {
    1798          14 :   pari_sp av = avma;
    1799          14 :   GEN H = K_get_H(K);
    1800          14 :   long d_K = K_get_d(K), f = K_get_f(K), h = K_get_h(K), g_K = K_get_g(K);
    1801          14 :   GEN theta = cgetg(d_K+1, t_VEC), xi = cgetg(d_K+1, t_VEC);
    1802          14 :   GEN vz_f = grootsof1(f, prec);
    1803          14 :   long i, j, g = 1, h2 = h>>1;
    1804             : 
    1805         238 :   for (i=1; i<=d_K; i++)
    1806             :   {
    1807         224 :     GEN x = real_0(prec), y = real_1(prec);
    1808        5068 :     for (j=1; j<=h2; j++)
    1809             :     {
    1810        4844 :       GEN z = gmael(vz_f, 1+Fl_mul(H[j], g, f), 1);
    1811        4844 :       x = addrr(x, z);
    1812        4844 :       y = mulrr(y, shiftr(subsr(1, z), 1));
    1813             :     }
    1814         224 :     gel(theta, i) = shiftr(x, 1);
    1815         224 :     gel(xi, i) = y;
    1816         224 :     g = Fl_mul(g, g_K, f);
    1817             :   }
    1818          14 :   return gerepilecopy(av, mkvec2(theta, xi));
    1819             : }
    1820             : 
    1821             : static GEN
    1822          14 : bound_coeff_xi(GEN K, GEN tInvA)
    1823             : {
    1824          14 :   pari_sp av = avma;
    1825          14 :   long d_K = K_get_d(K), prec = MEDDEFAULTPREC, i;
    1826          14 :   GEN tInvV, R = cgetg(d_K+1, t_MAT), theta_xi = theta_xi_approx(K, prec);
    1827          14 :   GEN theta = gel(theta_xi, 1), xi = gel(theta_xi, 2), x1, x2, bound;
    1828             : 
    1829         238 :   for (i=1; i<=d_K; i++)
    1830             :   {
    1831         224 :     GEN z = gpowers(gel(theta, i), d_K-1);
    1832         224 :     settyp(z, t_COL);
    1833         224 :     gel(R, i) = z;
    1834             :   }
    1835          14 :   tInvV = RgM_mul(RgM_inv(R), tInvA);
    1836          14 :   x1 = gsupnorm(tInvV, prec); x2 = gsupnorm(xi, prec);
    1837          14 :   bound = mulrs(mulrr(x1, x2), 3*d_K);
    1838          14 :   return gerepilecopy(av, bound);
    1839             : }
    1840             : 
    1841             : static GEN
    1842          14 : get_Xi(GEN K, GEN tInvA)
    1843             : {
    1844          14 :   pari_sp av = avma;
    1845             :   GEN M0, XI, EL, Xi;
    1846          14 :   ulong f = K_get_f(K), el, e, n, i;
    1847             :   forprime_t T0;
    1848             : 
    1849          14 :   M0 = bound_coeff_xi(K, tInvA);
    1850          14 :   e = expo(M0)+1; n = e/(BITS_IN_LONG-1); n++;
    1851          14 :   EL = cgetg(1+n, t_VECSMALL);
    1852          14 :   XI = cgetg(1+n, t_VEC);
    1853          14 :   u_forprime_arith_init(&T0, LONG_MAX, ULONG_MAX, 1, f);
    1854             : 
    1855          61 :   for (i=1; i<=n; i++)
    1856             :   {
    1857          47 :     el = u_forprime_next(&T0);
    1858          47 :     while ((Xi=Xi_el(K, tInvA, el))==NULL) el = u_forprime_next(&T0);
    1859          47 :     gel(XI, i) = Xi;
    1860          47 :     EL[i] = el;
    1861             :   }
    1862          14 :   return gerepileupto(av, nmV_chinese_center(XI, EL, NULL));
    1863             : }
    1864             : 
    1865             : /* K is a cyclic field of conductor f with degree d=d_K
    1866             :  * xi = Norm_{Q(zeta_f)/K}(1-zeta_f)
    1867             :  * 1: T, min poly of a=Tr_{Q(zeta_f)/K}(zeta_f) over Q
    1868             :  * 2: B, power basis of K with respect to a
    1869             :  * 3: A, rational matrix s.t. t(v_1,...v_d) = A*t(1,a,...,a^{d-1})
    1870             :  * 4: Xi, integer matrix s.t. t(xi^(1),...,xi^(d)) = Xi*t(v_1,...,v_d) */
    1871             : static GEN
    1872          14 : xi_data_basis(GEN K)
    1873             : {
    1874          14 :   pari_sp av = avma;
    1875          14 :   GEN T = minpol_theta(K);
    1876             :   GEN InvA, A, M, Xi, A_den;
    1877          14 :   GEN D, B = nfbasis(T, &D);
    1878             :   pari_timer ti;
    1879          14 :   if (DEBUGLEVEL>1) timer_start(&ti);
    1880          14 :   A = RgXV_to_RgM(B, lg(B)-1);
    1881          14 :   M = gmael(A, 1, 1);
    1882          14 :   if (!equali1(M)) A = RgM_Rg_div(A, M);
    1883          14 :   InvA = QM_inv(A);
    1884          14 :   A = Q_remove_denom(A, &A_den);
    1885          14 :   if (A_den==NULL) A_den = gen_1;
    1886          14 :   Xi = get_Xi(K, shallowtrans(InvA));
    1887          14 :   if (DEBUGLEVEL>1) timer_printf(&ti, "xi_data_basis");
    1888          14 :   return gerepilecopy(av, mkvec5(T, B, shallowtrans(A), Xi, A_den));
    1889             : }
    1890             : 
    1891             : /* When factorization of polcyclo mod el is difficult, one can try to
    1892             :  * check the condition of el using an integral basis of K.
    1893             :  * This is useful when d_K is small. */
    1894             : static long
    1895          14 : chk_el_real_basis(GEN K, long p, long d_pow, long el, long j0)
    1896             : {
    1897          14 :   pari_sp av = avma;
    1898          14 :   GEN xi = gel(K, 7), T = gel(xi, 1), A = gel(xi, 3), Xi = gel(xi, 4);
    1899          14 :   GEN A_den = gel(xi, 5);
    1900          14 :   ulong i, j, x, found = 0;
    1901             :   GEN v_el, xi_el;
    1902             :   GEN e_chi, xi_e_chi;
    1903             :   ulong d_K, d, dp, el1dp;
    1904             : 
    1905          14 :   if (dvdiu(A_den, el)) return 0;
    1906             : 
    1907          14 :   d = upowuu(p, d_pow); dp = d*p; el1dp = (el-1)/dp;
    1908          14 :   e_chi = get_e_chi(K, j0, dp, &d_K);
    1909          14 :   xi_e_chi = cgetg(d_K+1, t_VECSMALL)+1;
    1910             : 
    1911          14 :   if (DEBUGLEVEL>1) err_printf("chk_el_real_basis: d_K=%ld el=%ld\n",d_K,el);
    1912          14 :   A = ZM_to_Flm(A, el);
    1913          14 :   A = Flm_Fl_mul(A, Fl_inv(umodiu(A_den, el), el), el);
    1914          14 :   x = Flx_oneroot_split(ZX_to_Flx(T, el), el);
    1915          14 :   v_el = Flm_Flc_mul(A, Fl_powers(x, d_K-1, el), el);
    1916          14 :   xi_el = Flm_Flc_mul(ZM_to_Flm(Xi, el), v_el, el);
    1917          14 :   if (DEBUGLEVEL>2) err_printf("el=%ld xi_el=%Ps\n", el, xi_el);
    1918         238 :   for (i=0; i<d_K; i++)
    1919             :   {
    1920         224 :     long z = 1;
    1921        5208 :     for (j=0; j<d_K; j++)
    1922        4984 :       z = Fl_mul(z, Fl_powu(xi_el[1+j], e_chi[(i+j)%d_K], el), el);
    1923         224 :     xi_e_chi[i] = z;
    1924             :   }
    1925          14 :   if (DEBUGLEVEL>2) err_printf("xi_e_chi=%Ps\n", xi_e_chi-1);
    1926         238 :   for (i=0; i<d_K; i++)
    1927             :   {
    1928         224 :     long x = Fl_powu(xi_e_chi[i], el1dp, el);
    1929         224 :     if (x!=1) found = 1;
    1930         224 :     if (Fl_powu(x, p, el)!=1) return gc_long(av, 0);
    1931             :   }
    1932          14 :   return gc_long(av, found);
    1933             : }
    1934             : 
    1935             : static GEN
    1936           0 : bound_pol_xi(GEN K)
    1937             : {
    1938           0 :   pari_sp av = avma;
    1939           0 :   GEN xi = xi_approx(K, MEDDEFAULTPREC);
    1940           0 :   GEN M = real_1(MEDDEFAULTPREC), one = rtor(dbltor(1.0001), MEDDEFAULTPREC);
    1941           0 :   long i, n = lg(xi);
    1942             : 
    1943           0 :   for (i=1; i<n; i++) M = mulrr(M, addrr(one, gel(xi, i)));
    1944           0 :   M = mulrs(M, 3);
    1945           0 :   return gerepilecopy(av, M);
    1946             : }
    1947             : 
    1948             : static GEN
    1949           0 : minpol_xi(GEN K)
    1950             : {
    1951           0 :   pari_sp av = avma;
    1952             :   GEN M0, POL, EL;
    1953           0 :   ulong f = K_get_f(K), el, e, n, i;
    1954             :   forprime_t T0;
    1955             : 
    1956           0 :   M0 = bound_pol_xi(K);
    1957           0 :   e = expo(M0)+1; n = e/(BITS_IN_LONG-1); n++;
    1958           0 :   EL = cgetg(1+n, t_VECSMALL);
    1959           0 :   POL = cgetg(1+n, t_VEC);
    1960           0 :   u_forprime_arith_init(&T0, LONG_MAX, ULONG_MAX, 1, f);
    1961           0 :   for (i=1; i<=n; i++)
    1962             :   {
    1963           0 :     el = u_forprime_next(&T0);
    1964           0 :     gel(POL, i) = pol_xi_el(K, el);
    1965           0 :     EL[i] = el;
    1966             :   }
    1967           0 :   return gerepileupto(av, nxV_chinese_center(POL, EL, NULL));
    1968             : }
    1969             : 
    1970             : static long
    1971           0 : find_conj_el(GEN K, GEN pol, GEN Den)
    1972             : {
    1973           0 :   pari_sp av = avma;
    1974           0 :   GEN H = K_get_H(K);
    1975           0 :   ulong d_K = K_get_d(K), f = K_get_f(K), h = K_get_h(K), g = K_get_g(K);
    1976           0 :   ulong j, k, el, g_el, z_f, xi = 1, xi_g = 1;
    1977           0 :   GEN T = NULL, vz_f;
    1978             : 
    1979           0 :   for (el=f+1; el; el+=f)
    1980           0 :     if (uisprime(el) && dvdiu(Den, el)==0)
    1981             :     {
    1982           0 :       T = ZX_to_Flx(pol, el);
    1983           0 :       T = Flx_Fl_mul(T, Fl_inv(umodiu(Den, el), el), el);
    1984           0 :       break;
    1985             :     }
    1986           0 :   g_el = pgener_Fl(el);
    1987           0 :   z_f = Fl_powu(g_el, (el-1)/f, el);
    1988           0 :   vz_f = Fl_powers(z_f, f-1, el);
    1989           0 :   for (j=1; j<=h; j++)
    1990           0 :     xi = Fl_mul(xi, vz_f[1+H[j]]-1, el);
    1991           0 :   for (j=1; j<=h; j++)
    1992           0 :     xi_g = Fl_mul(xi_g, vz_f[1+Fl_mul(H[j], g, f)]-1, el);
    1993           0 :   for (k=1; k<=d_K; k++)
    1994             :   {
    1995           0 :     xi = Flx_eval(T, xi, el);
    1996           0 :     if (xi == xi_g) break;
    1997             :   }
    1998           0 :   if (xi != xi_g) pari_err_BUG("find_conj_el");
    1999           0 :   return gc_long(av, k);
    2000             : }
    2001             : 
    2002             : /* G = H_1*H_2*...*H_m is cyclic of order n, H_i=<perms[i]>
    2003             :  * G is not necessarily a direct product.
    2004             :  * If p^e || n, then p^e || |H_i| for some i.
    2005             :  * return a generator of G. */
    2006             : static GEN
    2007           0 : find_gen(GEN perms, long n)
    2008             : {
    2009           0 :   GEN fa = factoru(n), P = gel(fa, 1), E = gel(fa, 2);
    2010           0 :   long i, j, l = lg(perms), r = lg(P);
    2011           0 :   GEN gen = cgetg(r, t_VEC), orders = cgetg(l, t_VECSMALL), perm;
    2012             : 
    2013           0 :   for (i=1; i<l; i++) orders[i] = perm_orderu(gel(perms, i));
    2014           0 :   for (i=1; i<r; i++)
    2015             :   {
    2016           0 :     long pe = upowuu(P[i], E[i]);
    2017           0 :     for (j=1; j<l; j++) if (orders[j]%pe==0) break;
    2018           0 :     gel(gen, i) = perm_powu(gel(perms, j), orders[j]/pe);
    2019             :   }
    2020           0 :   perm = gel(gen, 1);
    2021           0 :   for (i=2; i<l; i++) perm = perm_mul(perm, gel(gen, i));
    2022           0 :   return perm;
    2023             : }
    2024             : 
    2025             : /* R is the roots of T. R[1+i] = R[1]^(g_K^i), 0 <= i <= d_K-1
    2026             :  * 1: min poly T of xi over Q
    2027             :  * 2: F(x)\in Q[x] s.t. xi^(g_K)=F(xi) */
    2028             : static GEN
    2029           0 : xi_data_galois(GEN K)
    2030             : {
    2031           0 :   pari_sp av = avma;
    2032             :   GEN T, G, perms, perm, pol, pol2, Den;
    2033           0 :   ulong k, d_K = K_get_d(K);
    2034             :   pari_timer ti;
    2035             : 
    2036           0 :   if (DEBUGLEVEL>1) timer_start(&ti);
    2037           0 :   T = minpol_xi(K);
    2038           0 :   if (DEBUGLEVEL>1) timer_printf(&ti, "minpol_xi");
    2039           0 :   G = galoisinit(T, NULL);
    2040           0 :   if (DEBUGLEVEL>1) timer_printf(&ti, "galoisinit");
    2041           0 :   perms = gal_get_gen(G);
    2042           0 :   perm = (lg(perms)==2)?gel(perms, 1):find_gen(perms, d_K);
    2043           0 :   if (DEBUGLEVEL>1) timer_start(&ti);
    2044           0 :   pol = galoispermtopol(G, perm);
    2045           0 :   if (DEBUGLEVEL>1) timer_printf(&ti, "galoispermtopol");
    2046           0 :   pol = Q_remove_denom(pol, &Den);
    2047           0 :   if (Den==NULL) Den = gen_1;
    2048           0 :   k = find_conj_el(K, pol, Den);
    2049           0 :   if (DEBUGLEVEL>1) timer_printf(&ti,"find_conj");
    2050           0 :   pol2 = galoispermtopol(G, perm_powu(perm, k));
    2051           0 :   pol2 = Q_remove_denom(pol2, &Den);
    2052           0 :   if (Den==NULL) Den = gen_1;
    2053           0 :   return gerepilecopy(av, mkvec3(T, pol2, Den));
    2054             : }
    2055             : 
    2056             : /* If g(X)\in Q[X] s.t. g(xi)=xi^{g_K} was found,
    2057             :  * then we fix an integer x_0 s.t. xi=x_0 (mod el) and construct x_i
    2058             :  * s.t. xi^{g_K^i}=x_i (mod el) using g(X). */
    2059             : static long
    2060           0 : chk_el_real_galois(GEN K, long p, long d_pow, long el, long j0)
    2061             : {
    2062           0 :   pari_sp av = avma;
    2063           0 :   GEN xi = gel(K, 7), T = gel(xi, 1), F = gel(xi, 2), Den = gel(xi, 3);
    2064             :   GEN Fel, xi_el, xi_e_chi, e_chi;
    2065           0 :   ulong i, j, x, found = 0;
    2066             :   ulong d_K, d, dp, el1dp;
    2067             : 
    2068           0 :   if (dvdiu(Den, el)) return 0;
    2069             : 
    2070           0 :   d = upowuu(p, d_pow); dp = d*p; el1dp = (el-1)/dp;
    2071           0 :   e_chi = get_e_chi(K, j0, dp, &d_K);
    2072           0 :   xi_el = cgetg(d_K+1, t_VECSMALL)+1;
    2073           0 :   xi_e_chi = cgetg(d_K+1, t_VECSMALL)+1;
    2074             : 
    2075           0 :   if (DEBUGLEVEL>1) err_printf("chk_el_real_galois: d_K=%ld el=%ld\n",d_K,el);
    2076           0 :   Fel = ZX_to_Flx(F, el);
    2077           0 :   Fel = Flx_Fl_mul(Fel, Fl_inv(umodiu(Den, el), el), el);
    2078           0 :   x = Flx_oneroot_split(ZX_to_Flx(T, el), el);
    2079           0 :   for (i=0; i<d_K; i++) { xi_el[i] = x; x = Flx_eval(Fel, x, el); }
    2080           0 :   if (DEBUGLEVEL>2) err_printf("el=%ld xi_el=%Ps\n", el, xi_el-1);
    2081           0 :   for (i=0; i<d_K; i++)
    2082             :   {
    2083           0 :     long z = 1;
    2084           0 :     for (j=0; j<d_K; j++)
    2085           0 :       z = Fl_mul(z, Fl_powu(xi_el[j], e_chi[(i+j)%d_K], el), el);
    2086           0 :     xi_e_chi[i] = z;
    2087             :   }
    2088           0 :   if (DEBUGLEVEL>2) err_printf("xi_e_chi=%Ps\n", xi_e_chi-1);
    2089           0 :   for (i=0; i<d_K; i++)
    2090             :   {
    2091           0 :     long x = Fl_powu(xi_e_chi[i], el1dp, el);
    2092           0 :     if (x!=1) found = 1;
    2093           0 :     if (Fl_powu(x, p, el)!=1) return gc_long(av, 0);
    2094             :   }
    2095           0 :   return gc_long(av, found);
    2096             : }
    2097             : 
    2098             : /* checks the condition of el using the irreducible polynomial G_K(X) of zeta_f
    2099             :  * over K. G_K(X) mod el is enough for our purpose and it is obtained by
    2100             :  * factoring polcyclo(f) mod el */
    2101             : static long
    2102           7 : chk_el_real_factor(GEN K, long p, long d_pow, long el, long j0)
    2103             : {
    2104           7 :   pari_sp av = avma;
    2105           7 :   GEN H = K_get_H(K);
    2106           7 :   ulong f = K_get_f(K), h = K_get_h(K), g_K = K_get_g(K);
    2107           7 :   ulong i, j, k, d_K, d = upowuu(p, d_pow), dp = d*p, found = 0;
    2108             :   GEN pols, coset, vn_g, polnum, xpows, G_K;
    2109           7 :   ulong el1dp = (el-1)/dp, n_coset, n_g, gi;
    2110           7 :   GEN e_chi = get_e_chi(K, j0, dp, &d_K);
    2111             :   pari_timer ti;
    2112             : 
    2113           7 :   if (DEBUGLEVEL>1) err_printf("chk_el_real_factor: f=%ld el=%ld\n",f,el);
    2114           7 :   coset = get_coset(H, h, f, el);
    2115           7 :   if (DEBUGLEVEL>1)
    2116             :   {
    2117           0 :     timer_start(&ti);
    2118           0 :     err_printf("factoring polyclo(d) (mod %ld)\n",f, el);
    2119             :   }
    2120           7 :   pols = Flx_factcyclo(f, el, 0);
    2121           7 :   if (DEBUGLEVEL>1) timer_printf(&ti,"Flx_factcyclo(%lu,%lu)",f,el);
    2122           7 :   n_coset = lg(coset)-1;
    2123           7 :   n_g = lg(pols)-1;
    2124           7 :   vn_g = identity_perm(n_g);
    2125             : 
    2126           7 :   polnum = const_vec(d_K, NULL);
    2127          91 :   for (i=1; i<=d_K; i++) gel(polnum, i) = const_vecsmall(n_coset, 0);
    2128           7 :   xpows = Flxq_powers(polx_Flx(0), f-1, gel(pols, 1), el);
    2129          91 :   for (gi=1,i=1; i<=d_K; i++)
    2130             :   {
    2131        3108 :     for (j=1; j<=n_coset; j++)
    2132             :     {
    2133        3024 :       long x, conj = Fl_mul(gi, coset[j], f);
    2134        3024 :       x = srh_pol(xpows, vn_g, pols, el, conj, f);
    2135        3024 :       gel(polnum, i)[j] = x;
    2136             :     }
    2137          84 :     gi = Fl_mul(gi, g_K, f);
    2138             :   }
    2139           7 :   G_K = const_vec(d_K, NULL);
    2140          91 :   for (i=1; i<=d_K; i++)
    2141             :   {
    2142          84 :     GEN z = pol1_Flx(0);
    2143        3108 :     for (j=1; j<=n_coset; j++) z = Flx_mul(z, gel(pols, gel(polnum,i)[j]), el);
    2144          84 :     gel(G_K, i) = z;
    2145             :   }
    2146           7 :   if (DEBUGLEVEL>2) err_printf("G_K(x)=%Ps\n",Flx_to_ZX(gel(G_K, 1)));
    2147          91 :   for (k=0; k<d_K; k++)
    2148             :   {
    2149          84 :     long x = 1;
    2150        1092 :     for (i = 0; i < d_K; i++)
    2151             :     {
    2152             :       long x0, t;
    2153        1008 :       x0 = Flx_eval(gel(G_K, 1+i), 1, el);
    2154        1008 :       t = Fl_powu(x0, e_chi[(i+k)%d_K], el);
    2155        1008 :       x = Fl_mul(x, t, el);
    2156             :     }
    2157          84 :     x = Fl_powu(x, el1dp, el);
    2158          84 :     if (x!=1) found = 1;
    2159          84 :     if (Fl_powu(x, p, el)!=1) return gc_long(av, 0);
    2160             :   }
    2161           7 :   return gc_long(av, found);
    2162             : }
    2163             : 
    2164             : static long
    2165          21 : chk_el_real_chi(GEN K, ulong p, ulong d_pow, ulong el, ulong j0, long flag)
    2166             : {
    2167          21 :   ulong f = K_get_f(K);
    2168             : 
    2169          21 :   if (el%f == 1)
    2170           0 :     return chk_el_real_f(K, p, d_pow, el, j0); /* must be faster */
    2171          21 :   if (flag&USE_BASIS)
    2172          14 :     return chk_el_real_basis(K, p, d_pow, el, j0);
    2173           7 :   if (flag&USE_GALOIS_POL)
    2174           0 :     return chk_el_real_galois(K, p, d_pow, el, j0);
    2175           7 :   return chk_el_real_factor(K, p, d_pow, el, j0);
    2176             : }
    2177             : 
    2178             : static long
    2179         616 : chk_ell_real(GEN K, long d2, GEN ell, long j0)
    2180             : {
    2181         616 :   pari_sp av = avma;
    2182         616 :   GEN H = K_get_H(K);
    2183         616 :   ulong f = K_get_f(K), h = K_get_h(K), g_K = K_get_g(K);
    2184             :   ulong d_K, i, j, gi;
    2185         616 :   GEN e_chi = get_e_chi(K, j0, d2, &d_K);
    2186         616 :   GEN g_ell, z_f, vz_f, xi_el = cgetg(d_K+1, t_VEC)+1;
    2187         616 :   GEN ell_1 = subiu(ell,1), ell1d2 = diviuexact(ell_1, d2);
    2188             : 
    2189         616 :   g_ell = pgener_Fp(ell);
    2190         616 :   z_f = Fp_pow(g_ell, diviuexact(subiu(ell, 1), f), ell);
    2191         616 :   vz_f = Fp_powers(z_f, f-1, ell)+1;
    2192       12964 :   for (gi=1, i=0; i<d_K; i++)
    2193             :   {
    2194       12348 :     GEN x = gen_1;
    2195     1074948 :     for (j = 1; j <= h; j++)
    2196             :     {
    2197     1062600 :       ulong y = Fl_mul(H[j], gi, f);
    2198     1062600 :       x = Fp_mul(x, subiu(gel(vz_f, y), 1), ell); /* vz_f[y] = z_f^y  */
    2199             :     }
    2200       12348 :     gi = Fl_mul(gi, g_K, f);
    2201       12348 :     gel(xi_el, i) = x;  /* xi_el[i]=xi^{g_K^i} */
    2202             :   }
    2203        1421 :   for (i=0; i<d_K; i++)
    2204             :   {
    2205        1386 :     GEN x = gen_1;
    2206       31976 :     for (j=0; j<d_K; j++)
    2207       30590 :       x = Fp_mul(x, Fp_powu(gel(xi_el, j), e_chi[(i+j)%d_K], ell), ell);
    2208        1386 :     x = Fp_pow(x, ell1d2, ell);
    2209        1386 :     if (!equali1(x)) return gc_long(av, 0);
    2210             :   }
    2211          35 :   return gc_long(av, 1);
    2212             : }
    2213             : 
    2214             : static GEN
    2215          21 : next_el_real(GEN K, long p, long d_pow, GEN elg, long j0, long flag)
    2216             : {
    2217          21 :   GEN Chi = gel(K, 2);
    2218          21 :   ulong f = K_get_f(K), h = K_get_h(K), d = upowuu(p, d_pow), d2 = d*d;
    2219          21 :   ulong D = (flag & USE_F)? d2*f: d2<<1, el = elg[1] + D;
    2220             : 
    2221             :   /* O(el*h) -> O(el*log(el)) by FFT */
    2222          21 :   if (1000 < h && el < h) { el = (h/D)*D+1; if (el < h) el += D; }
    2223          21 :   if (flag&USE_F)  /* el=1 (mod f) */
    2224             :   {
    2225           0 :     for (;; el += D)
    2226           0 :       if (uisprime(el) && chk_el_real_f(K, p, d_pow, el, j0)) break;
    2227             :   }
    2228             :   else
    2229             :   {
    2230         413 :     for (;; el += D)
    2231         455 :       if (Chi[el%f]==0 && uisprime(el) &&
    2232          42 :           chk_el_real_chi(K, p, d_pow, el, j0, flag)) break;
    2233             :   }
    2234          21 :   return mkvecsmall2(el, pgener_Fl(el));
    2235             : }
    2236             : 
    2237             : static GEN
    2238          35 : next_ell_real(GEN K, GEN ellg, long d2, GEN df0l, long j0)
    2239             : {
    2240          35 :   GEN ell = gel(ellg, 1);
    2241        7602 :   for (ell = addii(ell, df0l);; ell = addii(ell, df0l))
    2242        7602 :     if (BPSW_psp(ell) && chk_ell_real(K, d2, ell, j0))
    2243          35 :       return mkvec2(ell, pgener_Fp(ell));
    2244             : }
    2245             : 
    2246             : /* #velg >= n */
    2247             : static long
    2248           0 : delete_el(GEN velg, long n)
    2249             : {
    2250             :   long i, l;
    2251           0 :   if (DEBUGLEVEL>1) err_printf("deleting %ld ...\n", gmael(velg, n, 1));
    2252           0 :   for (l = lg(velg)-1; l >= 1; l--) if (gel(velg, l)) break;
    2253           0 :   for (i = n; i < l; i++) gel(velg, i) = gel(velg, i+1);
    2254           0 :   return l;
    2255             : }
    2256             : 
    2257             : /* velg has n components */
    2258             : static GEN
    2259          21 : set_ell_real(GEN K, GEN velg, long n, long d_chi, long d2, long f0, long j0)
    2260             : {
    2261          21 :   long i, n_ell = n*d_chi;
    2262          21 :   GEN z = cgetg(n_ell + 1, t_VEC);
    2263          21 :   GEN df0l = muluu(d2, f0), ellg = mkvec2(gen_1, gen_1);
    2264          42 :   for (i=1; i<=n; i++) df0l = muliu(df0l, gel(velg, i)[1]);
    2265          56 :   for (i=1; i<=n_ell; i++) ellg = gel(z, i)= next_ell_real(K, ellg, d2, df0l, j0);
    2266          21 :   return z;
    2267             : }
    2268             : 
    2269             : static GEN
    2270         182 : G_K_vell(GEN K, GEN vellg, ulong gk)
    2271             : {
    2272         182 :   pari_sp av = avma;
    2273         182 :   GEN H = K_get_H(K);
    2274         182 :   ulong f = K_get_f(K), h = K_get_h(K);
    2275         182 :   GEN z_f, vz_f, A, P, M, z =  cgetg(h+1, t_VEC);
    2276         182 :   ulong i, lv = lg(vellg);
    2277             : 
    2278         182 :   A=cgetg(lv, t_VEC);
    2279         182 :   P=cgetg(lv, t_VEC);
    2280         728 :   for (i=1; i<lv; i++)
    2281             :   {
    2282         546 :     GEN ell = gmael(vellg, i, 1), g_ell = gmael(vellg, i, 2);
    2283         546 :     gel(A, i) = Fp_pow(g_ell, diviuexact(subiu(ell, 1), f), ell);
    2284         546 :     gel(P, i) = ell;
    2285             :   }
    2286         182 :   z_f = ZV_chinese(A, P, &M);
    2287         182 :   vz_f = Fp_powers(z_f, f-1, M)+1;
    2288        3822 :   for (i=1; i<=h; i++) gel(z, i) = gel(vz_f, Fl_mul(H[i], gk, f));
    2289         182 :   return gerepilecopy(av, FpV_roots_to_pol(z, M, 0));
    2290             : }
    2291             : 
    2292             : /* f=cond(K), M=product of ell in vell, G(K/Q)=<g_K>
    2293             :  * G_K[1+i]=minimal polynomial of zeta_f^{g_k^i} over K mod M, 0 <= i < d_K */
    2294             : static GEN
    2295           7 : make_G_K(GEN K, GEN vellg)
    2296             : {
    2297           7 :   ulong d_K = K_get_d(K), f = K_get_f(K), g_K = K_get_g(K);
    2298           7 :   GEN G_K = cgetg(d_K+1, t_VEC);
    2299           7 :   ulong i, g = 1;
    2300             : 
    2301         189 :   for (i=0; i<d_K; i++)
    2302             :   {
    2303         182 :     gel(G_K, 1+i) = G_K_vell(K, vellg, g);
    2304         182 :     g = Fl_mul(g, g_K, f);
    2305             :   }
    2306           7 :   return G_K;
    2307             : }
    2308             : 
    2309             : static GEN
    2310          12 : G_K_p(GEN K, GEN ellg, ulong gk)
    2311             : {
    2312          12 :   pari_sp av = avma;
    2313          12 :   ulong i, f = K_get_f(K), h = K_get_h(K);
    2314          12 :   GEN ell = gel(ellg, 1), g_ell = gel(ellg, 2);
    2315          12 :   GEN H = K_get_H(K), z_f, vz_f, z = cgetg(h+1, t_VEC);
    2316             : 
    2317          12 :   z_f = Fp_pow(g_ell, diviuexact(subiu(ell, 1), f), ell);
    2318          12 :   vz_f = Fp_powers(z_f, f-1, ell)+1;
    2319        3468 :   for (i=1; i<=h; i++) gel(z, i) = gel(vz_f, Fl_mul(H[i], gk, f));
    2320          12 :   return gerepilecopy(av, FpV_roots_to_pol(z, ell, 0));
    2321             : }
    2322             : 
    2323             : static GEN
    2324         114 : G_K_l(GEN K, GEN ellg, ulong gk)
    2325             : {
    2326         114 :   pari_sp av = avma;
    2327         114 :   ulong ell = itou(gel(ellg, 1)), g_ell = itou(gel(ellg, 2));
    2328         114 :   ulong f = K_get_f(K), h = K_get_h(K), i, z_f;
    2329         114 :   GEN H = K_get_H(K), vz_f, z = cgetg(h+1, t_VEC);
    2330             : 
    2331         114 :   z_f = Fl_powu(g_ell, (ell-1) / f, ell);
    2332         114 :   vz_f = Fl_powers(z_f, f-1, ell)+1;
    2333       26898 :   for (i=1; i<=h; i++) z[i] = vz_f[Fl_mul(H[i], gk, f)];
    2334         114 :   return gerepilecopy(av, Flv_roots_to_pol(z, ell, 0));
    2335             : }
    2336             : 
    2337             : static GEN
    2338           6 : vz_vell(long d, GEN vellg, GEN *pM)
    2339             : {
    2340           6 :   long i, l = lg(vellg);
    2341           6 :   GEN A = cgetg(l, t_VEC), P = cgetg(l, t_VEC), z;
    2342             : 
    2343          18 :   for (i = 1; i < l; i++)
    2344             :   {
    2345          12 :     GEN ell = gmael(vellg, i, 1), g_ell = gmael(vellg, i, 2);
    2346          12 :     gel(A, i) = Fp_pow(g_ell, diviuexact(subiu(ell, 1), d), ell);
    2347          12 :     gel(P, i) = ell;
    2348             :   }
    2349           6 :   z = ZV_chinese(A, P, pM); return Fp_powers(z, d-1, *pM);
    2350             : }
    2351             : 
    2352             : static GEN
    2353           0 : D_xi_el_vell_FFT(GEN K, GEN elg, GEN vellg, ulong d, ulong j0, GEN vG_K)
    2354             : {
    2355           0 :   pari_sp av = avma;
    2356           0 :   ulong d_K, h = K_get_h(K), d_chi = K_get_dchi(K);
    2357           0 :   ulong el = elg[1], g_el = elg[2], el_1 = el-1;
    2358             :   ulong i, j, i2, k, dwel;
    2359           0 :   GEN u = cgetg(el+2, t_POL) , v = cgetg(h+3, t_POL);
    2360           0 :   GEN w = cgetg(el+1, t_VEC), ww;
    2361           0 :   GEN M, vz_el, G_K, z = const_vec(d_chi, gen_1);
    2362           0 :   GEN e_chi = get_e_chi(K, j0, d, &d_K);
    2363             : 
    2364           0 :   vz_el = vz_vell(el, vellg, &M);
    2365           0 :   u[1] = evalsigne(1) | evalvarn(0);
    2366           0 :   v[1] = evalsigne(1) | evalvarn(0);
    2367             : 
    2368           0 :   for (i=i2=0; i<el; i++)
    2369             :   {
    2370           0 :     ulong j2 = i2?el-i2:i2; /* i2=(i*i)%el */
    2371           0 :     gel(u, 2+i) = gel(vz_el, 1+j2);
    2372           0 :     if ((i2+=i+i+1)>=el) i2%=el;
    2373             :   }
    2374           0 :   for (k=0; k<d_K; k++)
    2375             :   {
    2376           0 :     pari_sp av = avma;
    2377             :     pari_timer ti;
    2378             :     long gd, gi;
    2379           0 :     GEN x1 = gen_1;
    2380           0 :     G_K = gel(vG_K, 1+k);
    2381           0 :     for (i=i2=0; i<=h; i++)
    2382             :     {
    2383           0 :       gel(v, 2+i) = Fp_mul(gel(G_K, 2+i), gel(vz_el, 1+i2), M);
    2384           0 :       if ((i2+=i+i+1)>=el) i2%=el;
    2385             :     }
    2386           0 :     if (DEBUGLEVEL>2) timer_start(&ti);
    2387           0 :     ww = ZX_mul(u, v);
    2388           0 :     if (DEBUGLEVEL>2)
    2389           0 :       timer_printf(&ti, "ZX_mul:%ld/%ld h*el=%ld*%ld", k, d_K, h, el);
    2390           0 :     dwel = degpol(ww)-el;
    2391           0 :     for (i=0; i<=dwel; i++) gel(w, 1+i) = addii(gel(ww, 2+i), gel(ww, 2+i+el));
    2392           0 :     for (; i<el; i++) gel(w, 1+i) = gel(ww, 2+i);
    2393           0 :     for (i=i2=1; i<el; i++)  /* w[i]=G_K(z_el^(2*i)) */
    2394             :     {
    2395           0 :       gel(w, i) = Fp_mul(gel(w, 1+i), gel(vz_el, 1+i2), M);
    2396           0 :       if ((i2+=i+i+1)>=el) i2%=el;
    2397             :     }
    2398           0 :     gd = Fl_powu(g_el, d, el);  /* a bit faster */
    2399           0 :     gi = g_el;
    2400           0 :     for (i=1; i<d; i++)
    2401             :     {
    2402           0 :       GEN xi = gen_1;
    2403           0 :       long gdi = gi;
    2404           0 :       for (j=0; i+j<el_1; j+=d)
    2405             :       {
    2406           0 :         xi = Fp_mul(xi, gel(w, (gdi+gdi)%el), M);
    2407           0 :         gdi = Fl_mul(gdi, gd, el);
    2408             :       }
    2409           0 :       gi = Fl_mul(gi, g_el, el);
    2410           0 :       xi = Fp_powu(xi, i, M);
    2411           0 :       x1 = Fp_mul(x1, xi, M);
    2412             :     }
    2413           0 :     for (i=1; i<=d_chi; i++)
    2414             :     {
    2415           0 :       GEN x2 = Fp_powu(x1, e_chi[(k+i-1)%d_K], M);
    2416           0 :       gel(z, i) = Fp_mul(gel(z, i), x2, M);
    2417             :     }
    2418           0 :     z = gerepilecopy(av, z);
    2419             :   }
    2420           0 :   return gerepilecopy(av, z);
    2421             : }
    2422             : 
    2423             : static GEN
    2424           0 : D_xi_el_vell(GEN K, GEN elg, GEN vellg, ulong d, ulong j0)
    2425             : {
    2426           0 :   pari_sp av = avma;
    2427           0 :   GEN H = K_get_H(K);
    2428           0 :   ulong f = K_get_f(K), h = K_get_h(K), g_K = K_get_g(K);
    2429             :   GEN z_f, z_el, vz_f, vz_el;
    2430           0 :   ulong el = elg[1], g_el = elg[2], el_1 = el-1;
    2431           0 :   ulong i, j, k, d_K, lv = lg(vellg), d_chi = K_get_dchi(K);
    2432           0 :   GEN A, B, P, M, z = const_vec(d_chi, gen_1);
    2433           0 :   GEN e_chi = get_e_chi(K, j0, d, &d_K);
    2434             : 
    2435           0 :   A=cgetg(lv, t_VEC);
    2436           0 :   B=cgetg(lv, t_VEC);
    2437           0 :   P=cgetg(lv, t_VEC);
    2438           0 :   for (i = 1; i < lv; i++)
    2439             :   {
    2440           0 :     GEN ell = gmael(vellg, i, 1), g_ell = gmael(vellg, i, 2);
    2441           0 :     GEN ell_1 = subiu(ell, 1);
    2442           0 :     gel(A, i) = Fp_pow(g_ell, diviuexact(ell_1, f), ell);
    2443           0 :     gel(B, i) = Fp_pow(g_ell, diviuexact(ell_1, el), ell);
    2444           0 :     gel(P, i) = ell;
    2445             :   }
    2446           0 :   z_f = ZV_chinese(A, P, &M);
    2447           0 :   z_el = ZV_chinese(B, P, NULL);
    2448           0 :   vz_f = Fp_powers(z_f, f-1, M);
    2449           0 :   vz_el = Fp_powers(z_el, el-1, M);
    2450           0 :   for (k = 0; k < d_K; k++)
    2451             :   {
    2452           0 :     pari_sp av = avma;
    2453           0 :     GEN x0 = gen_1;
    2454           0 :     long gk = Fl_powu(g_K, k, f);
    2455           0 :     for (i=1; i<el_1; i++)
    2456             :     {
    2457           0 :       long gi = Fl_powu(g_el, i, el);
    2458           0 :       GEN x1 = gen_1;
    2459           0 :       GEN x2 = gel(vz_el, 1+gi);
    2460           0 :       for (j=1; j<=h; j++)
    2461             :       {
    2462           0 :         long y = Fl_mul(H[j], gk, f);
    2463           0 :         x1 = Fp_mul(x1, Fp_sub(x2, gel(vz_f, 1+y), M), M);
    2464             :       }
    2465           0 :       x1 = Fp_powu(x1, i%d, M);
    2466           0 :       x0 = Fp_mul(x0, x1, M);
    2467             :     }
    2468           0 :     for (i=1; i<=d_chi; i++)
    2469             :     {
    2470           0 :       GEN x2 = Fp_powu(x0, e_chi[(k+i-1)%d_K], M);
    2471           0 :       gel(z, i) = Fp_mul(gel(z, i), x2, M);
    2472             :     }
    2473           0 :     z = gerepilecopy(av, z);
    2474             :   }
    2475           0 :   return gerepilecopy(av, z);
    2476             : }
    2477             : 
    2478             : static GEN
    2479          34 : D_xi_el_Flx_mul(GEN K, GEN elg, GEN ellg, GEN vG_K, ulong d, ulong j0)
    2480             : {
    2481          34 :   pari_sp av = avma;
    2482          34 :   ulong d_K, f = K_get_f(K), h = K_get_h(K), g_K = K_get_g(K);
    2483          34 :   ulong el = elg[1], g_el = elg[2], el_1 = el-1, d_chi = K_get_dchi(K);
    2484          34 :   ulong ell = itou(gel(ellg, 1)), g_ell = itou(gel(ellg, 2)), z_el;
    2485          34 :   GEN u = cgetg(el+2, t_VECSMALL), v = cgetg(h+3, t_VECSMALL);
    2486          34 :   GEN w = cgetg(el+1, t_VECSMALL), ww;
    2487          34 :   GEN vz_el, G_K, z = const_vecsmall(d_chi, 1);
    2488          34 :   GEN e_chi = get_e_chi(K, j0, d, &d_K);
    2489             :   ulong i, j, i2, k, dwel;
    2490             : 
    2491          34 :   u[1] = evalvarn(0);
    2492          34 :   v[1] = evalvarn(0);
    2493          34 :   z_el = Fl_powu(g_ell, (ell - 1) / el, ell);
    2494          34 :   vz_el = Fl_powers(z_el, el_1, ell)+1;
    2495             : 
    2496      700376 :   for (i=i2=0; i<el; i++)
    2497             :   {
    2498      700342 :     ulong j2 = i2?el-i2:i2;
    2499      700342 :     u[2+i] = vz_el[j2];
    2500      700342 :     if ((i2+=i+i+1)>=el) i2%=el;  /* i2=(i*i)%el */
    2501             :   }
    2502         694 :   for (k=0; k<d_K; k++)
    2503             :   {
    2504         660 :     pari_sp av = avma;
    2505             :     pari_timer ti;
    2506         660 :     ulong gk = Fl_powu(g_K, k, f);
    2507         660 :     long gd, gi, x1 = 1;
    2508         660 :     if (DEBUGLEVEL>2) timer_start(&ti);
    2509         660 :     G_K = (vG_K==NULL)?G_K_l(K, ellg, gk):ZX_to_Flx(gel(vG_K, 1+k), ell);
    2510         660 :     if (DEBUGLEVEL>2) timer_printf(&ti, "G_K_l");
    2511       39024 :     for (i=i2=0; i<=h; i++)
    2512             :     {
    2513       38364 :       v[2+i] = Fl_mul(G_K[2+i], vz_el[i2], ell);
    2514       38364 :       if ((i2+=i+i+1)>=el) i2%=el;  /* i2=(i*i)%el */
    2515             :     }
    2516         660 :     if (DEBUGLEVEL>2) timer_start(&ti);
    2517         660 :     ww = Flx_mul(u, v, ell);
    2518         660 :     if (DEBUGLEVEL>2)
    2519           0 :       timer_printf(&ti, "Flx_mul:%ld/%ld h*el=%ld*%ld", k, d_K, h, el);
    2520         660 :     dwel=degpol(ww)-el; /* dwel=h-1 */
    2521       38364 :     for (i=0; i<=dwel; i++) w[1+i] = Fl_add(ww[2+i], ww[2+i+el], ell);
    2522     8388480 :     for (; i<el; i++) w[1+i] = ww[2+i];
    2523     8425524 :     for (i=i2=1; i<el; i++)  /* w[i]=G_K(z_el^(2*i)) */
    2524             :     {
    2525     8424864 :       w[i] = Fl_mul(w[1+i], vz_el[i2], ell);
    2526     8424864 :       if ((i2+=i+i+1)>=el) i2%=el;  /* i2=(i*i)%el */
    2527             :     }
    2528         660 :     gd = Fl_powu(g_el, d, el);  /* a bit faster */
    2529         660 :     gi = g_el;
    2530        4596 :     for (i=1; i<d; i++)
    2531             :     {
    2532        3936 :       long xi = 1, gdi = gi;
    2533     8163696 :       for (j=0; i+j<el_1; j+=d)
    2534             :       {
    2535     8159760 :         xi = Fl_mul(xi, w[(gdi+gdi)%el], ell);
    2536     8159760 :         gdi = Fl_mul(gdi, gd, el);
    2537             :       }
    2538        3936 :       gi = Fl_mul(gi, g_el, el);
    2539        3936 :       xi = Fl_powu(xi, i, ell);
    2540        3936 :       x1 = Fl_mul(x1, xi, ell);
    2541             :     }
    2542        2412 :     for (i=1; i<=d_chi; i++)
    2543             :     {
    2544        1752 :       long x2 = Fl_powu(x1, e_chi[(k+i-1)%d_K], ell);
    2545        1752 :       z[i] = Fl_mul(z[i], x2, ell);
    2546             :     }
    2547         660 :     set_avma(av);
    2548             :   }
    2549          34 :   return gerepilecopy(av, Flv_to_ZV(z));
    2550             : }
    2551             : 
    2552             : static GEN
    2553          35 : D_xi_el_ZX_mul(GEN K, GEN elg, GEN ellg, GEN vG_K, ulong d, ulong j0)
    2554             : {
    2555          35 :   pari_sp av = avma;
    2556          35 :   GEN ell = gel(ellg,1), g_ell, u, v, w, ww, z_el, vz_el, G_K, z, e_chi;
    2557             :   ulong d_K, f, h, g_K, el, g_el, el_1, d_chi, i, j, i2, k, dwel;
    2558             : 
    2559          35 :   if (lgefint(ell) == 3) return D_xi_el_Flx_mul(K, elg, ellg, vG_K, d, j0);
    2560           1 :   f = K_get_f(K); h = K_get_h(K); g_K = K_get_g(K);
    2561           1 :   el = elg[1]; g_el = elg[2]; el_1 = el-1; d_chi = K_get_dchi(K);
    2562           1 :   g_ell = gel(ellg, 2);
    2563           1 :   z = const_vec(d_chi, gen_1);
    2564           1 :   e_chi = get_e_chi(K, j0, d, &d_K);
    2565             : 
    2566           1 :   u = cgetg(el+2,t_POL); u[1] = evalsigne(1) | evalvarn(0);
    2567           1 :   v = cgetg(h+3, t_POL); v[1] = evalsigne(1) | evalvarn(0);
    2568           1 :   w = cgetg(el+1, t_VEC);
    2569           1 :   z_el = Fp_pow(g_ell, diviuexact(subiu(ell, 1), el), ell);
    2570           1 :   vz_el = Fp_powers(z_el, el_1, ell)+1;
    2571             : 
    2572      114998 :   for (i=i2=0; i<el; i++)
    2573             :   {
    2574      114997 :     ulong j2 = i2?el-i2:i2; /* i2=(i*i)%el */
    2575      114997 :     gel(u, 2+i) = gel(vz_el, j2);
    2576      114997 :     if ((i2+=i+i+1)>=el) i2%=el;
    2577             :   }
    2578          13 :   for (k=0; k<d_K; k++)
    2579             :   {
    2580          12 :     pari_sp av = avma;
    2581             :     pari_timer ti;
    2582          12 :     long gd, gi, gk = Fl_powu(g_K, k, f);
    2583          12 :     GEN x1 = gen_1;
    2584          12 :     if (DEBUGLEVEL>2) timer_start(&ti);
    2585          12 :     G_K = (vG_K==NULL) ? G_K_p(K, ellg, gk):RgX_to_FpX(gel(vG_K, 1+k), ell);
    2586          12 :     if (DEBUGLEVEL>2) timer_printf(&ti, "G_K_p");
    2587        3480 :     for (i=i2=0; i<=h; i++)
    2588             :     {
    2589        3468 :       gel(v, 2+i) = Fp_mul(gel(G_K, 2+i), gel(vz_el, i2), ell);
    2590        3468 :       if ((i2+=i+i+1)>=el) i2%=el;
    2591             :     }
    2592          12 :     if (DEBUGLEVEL>2) timer_start(&ti);
    2593          12 :     ww = ZX_mul(u, v);
    2594          12 :     if (DEBUGLEVEL>2)
    2595           0 :       timer_printf(&ti, "ZX_mul:%ld/%ld h*el=%ld*%ld", k, d_K, h, el);
    2596          12 :     dwel = degpol(ww)-el;
    2597        3468 :     for (i=0; i<=dwel; i++) gel(w, 1+i) = addii(gel(ww, 2+i), gel(ww, 2+i+el));
    2598     1376520 :     for (; i<el; i++) gel(w, 1+i) = gel(ww, 2+i);
    2599     1379964 :     for (i=i2=1; i<el; i++)  /* w[i]=G_K(z_el^(2*i)) */
    2600             :     {
    2601     1379952 :       gel(w, i) = Fp_mul(gel(w, 1+i), gel(vz_el, i2), ell);
    2602     1379952 :       if ((i2+=i+i+1)>=el) i2%=el;
    2603             :     }
    2604          12 :     gd = Fl_powu(g_el, d, el);  /* a bit faster */
    2605          12 :     gi = g_el;
    2606         444 :     for (i=1; i<d; i++)
    2607             :     {
    2608         432 :       GEN xi = gen_1;
    2609         432 :       long gdi = gi;
    2610     1343088 :       for (j=0; i+j<el_1; j+=d)
    2611             :       {
    2612     1342656 :         xi = Fp_mul(xi, gel(w, (gdi+gdi)%el), ell);
    2613     1342656 :         gdi = Fl_mul(gdi, gd, el);
    2614             :       }
    2615         432 :       gi = Fl_mul(gi, g_el, el);
    2616         432 :       xi = Fp_powu(xi, i, ell);
    2617         432 :       x1 = Fp_mul(x1, xi, ell);
    2618             :     }
    2619          24 :     for (i=1; i<=d_chi; i++)
    2620             :     {
    2621          12 :       GEN x2 = Fp_powu(x1, e_chi[(k+i-1)%d_K], ell);
    2622          12 :       gel(z, i) = Fp_mul(gel(z, i), x2, ell);
    2623             :     }
    2624          12 :     z = gerepilecopy(av, z);
    2625             :   }
    2626           1 :   return gerepilecopy(av, z);
    2627             : }
    2628             : 
    2629             : static GEN
    2630           0 : D_xi_el_ss(GEN K, GEN elg, GEN ellg, ulong d, ulong j0)
    2631             : {
    2632           0 :   pari_sp av = avma;
    2633           0 :   GEN H = K_get_H(K);
    2634           0 :   ulong d_K, f = K_get_f(K), h = K_get_h(K), g_K = K_get_g(K);
    2635           0 :   ulong el = elg[1], g_el = elg[2], el_1 = el-1;
    2636           0 :   ulong ell = itou(gel(ellg, 1)), g_ell = itou(gel(ellg, 2));
    2637           0 :   ulong i, j, k, gk, z_f, z_el, d_chi = K_get_dchi(K);
    2638           0 :   GEN vz_f, vz_el, z = const_vecsmall(d_chi, 1);
    2639           0 :   GEN e_chi = get_e_chi(K, j0, d, &d_K);
    2640             : 
    2641           0 :   z_f = Fl_powu(g_ell, (ell - 1) / f, ell);
    2642           0 :   z_el = Fl_powu(g_ell, (ell - 1) / el, ell);
    2643           0 :   vz_f = Fl_powers(z_f, f-1, ell)+1;
    2644           0 :   vz_el = Fl_powers(z_el, el-1, ell)+1;
    2645           0 :   gk = 1; /* g_K^k */
    2646           0 :   for (k = 0; k < d_K; k++)
    2647             :   {
    2648           0 :     ulong x0 = 1, gi = g_el; /* g_el^i */
    2649           0 :     for (i = 1; i < el_1; i++)
    2650             :     {
    2651           0 :       ulong x1 = 1, x2 = vz_el[gi];
    2652           0 :       for (j=1; j<=h; j++)
    2653             :       {
    2654           0 :         ulong y = Fl_mul(H[j], gk, f);
    2655           0 :         x1 = Fl_mul(x1, Fl_sub(x2, vz_f[y], ell), ell);
    2656             :       }
    2657           0 :       x1 = Fl_powu(x1, i%d, ell);
    2658           0 :       x0 = Fl_mul(x0, x1, ell);
    2659           0 :       gi = Fl_mul(gi, g_el, el);
    2660             :     }
    2661           0 :     for (i = 1; i <= d_chi; i++)
    2662             :     {
    2663           0 :       ulong x2 = Fl_powu(x0, e_chi[(k+i-1)%d_K], ell);
    2664           0 :       z[i] = Fl_mul(z[i], x2, ell);
    2665             :     }
    2666           0 :     gk = Fl_mul(gk, g_K, f);
    2667             :   }
    2668           0 :   return gerepileupto(av, Flv_to_ZV(z));
    2669             : }
    2670             : 
    2671             : static GEN
    2672           0 : D_xi_el_sl(GEN K, GEN elg, GEN ellg, ulong d, ulong j0)
    2673             : {
    2674           0 :   pari_sp av = avma;
    2675           0 :   GEN ell = gel(ellg, 1), H;
    2676             :   GEN g_ell, ell_1, z_f, z_el, vz_f, vz_el, z, e_chi;
    2677             :   ulong d_K, f, h, g_K, el, g_el, el_1, d_chi, i, j, k, gk;
    2678             : 
    2679           0 :   if (lgefint(ell) == 3) return D_xi_el_ss(K, elg, ellg, d, j0);
    2680           0 :   H = K_get_H(K);
    2681           0 :   f = K_get_f(K); h = K_get_h(K); g_K = K_get_g(K);
    2682           0 :   el = elg[1]; g_el = elg[2]; el_1 = el-1; d_chi = K_get_dchi(K);
    2683           0 :   g_ell = gel(ellg, 2); ell_1 = subiu(ell, 1);
    2684           0 :   z = const_vec(d_chi, gen_1);
    2685           0 :   e_chi = get_e_chi(K, j0, d, &d_K);
    2686             : 
    2687           0 :   z_f = Fp_pow(g_ell, diviuexact(ell_1, f), ell);
    2688           0 :   z_el = Fp_pow(g_ell, diviuexact(ell_1, el), ell);
    2689           0 :   vz_f = Fp_powers(z_f, f-1, ell) + 1;
    2690           0 :   vz_el = Fp_powers(z_el, el-1, ell) + 1;
    2691           0 :   gk = 1; /* g_K^k */
    2692           0 :   for (k = 0; k < d_K; k++)
    2693             :   {
    2694           0 :     pari_sp av2 = avma;
    2695           0 :     GEN x0 = gen_1;
    2696           0 :     ulong gi = g_el; /* g_el^i */
    2697           0 :     for (i = 1; i < el_1; i++)
    2698             :     {
    2699           0 :       pari_sp av3 = avma;
    2700           0 :       GEN x1 = gen_1, x2 = gel(vz_el, gi);
    2701           0 :       for (j = 1; j <= h; j++)
    2702             :       {
    2703           0 :         ulong y = Fl_neg(Fl_mul(H[j], gk, f), f);
    2704           0 :         x1 = Fp_mul(x1, Fp_sub(x2, gel(vz_f, y), ell), ell);
    2705             :       }
    2706           0 :       x1 = Fp_powu(x1, i%d, ell);
    2707           0 :       x0 = gerepileuptoint(av3, Fp_mul(x0, x1, ell));
    2708           0 :       gi = Fl_mul(gi, g_el, el);
    2709             :     }
    2710           0 :     for (i=1; i<=d_chi; i++)
    2711             :     {
    2712           0 :       GEN x2 = Fp_powu(x0, e_chi[(k+i-1)%d_K], ell);
    2713           0 :       gel(z, i) = Fp_mul(gel(z, i), x2, ell);
    2714             :     }
    2715           0 :     if (k == d_K-1) break;
    2716           0 :     z = gerepilecopy(av2, z);
    2717           0 :     gk = Fl_mul(gk, g_K, f);
    2718             :   }
    2719           0 :   return gerepilecopy(av, z);
    2720             : }
    2721             : 
    2722             : static long
    2723         175 : get_y(GEN z, GEN ellg, long d)
    2724             : {
    2725         175 :   GEN ell = gel(ellg, 1), g_ell = gel(ellg, 2);
    2726         175 :   GEN elld = diviuexact(subiu(ell, 1), d);
    2727         175 :   GEN g_elld = Fp_pow(g_ell, elld, ell);
    2728         175 :   GEN x = Fp_pow(modii(z, ell), elld, ell);
    2729             :   long k;
    2730      211778 :   for (k=0; k<d; k++)
    2731             :   {
    2732      211778 :     if (equali1(x)) break;
    2733      211603 :     x = Fp_mul(x, g_elld, ell);
    2734             :   }
    2735         175 :   if (k==0) k=d;
    2736         161 :   else if (d<=k) pari_err_BUG("subcyclopclgp [MLL]");
    2737         175 :   return k;
    2738             : }
    2739             : 
    2740             : static void
    2741           0 : real_MLLn(long *y, GEN K, ulong p, ulong d_pow, ulong n,
    2742             :     GEN velg, GEN vellg, GEN vG_K, ulong j0)
    2743             : {
    2744           0 :   pari_sp av = avma;
    2745           0 :   ulong i, j, k, d = upowuu(p, d_pow), h = gmael(K, 1, 2)[3];
    2746           0 :   ulong row = lg(vellg)-1;
    2747           0 :   for (i=1; i<=n; i++)
    2748             :   {
    2749           0 :     GEN elg = gel(velg, i), z;
    2750           0 :     ulong el = elg[1], nz;
    2751             :     pari_timer ti;
    2752           0 :     if (DEBUGLEVEL>1) timer_start(&ti);
    2753           0 :     z = (h<el) ? D_xi_el_vell_FFT(K, elg, vellg, d, j0, vG_K)
    2754           0 :                : D_xi_el_vell(K, elg, vellg, d, j0);
    2755           0 :     if (DEBUGLEVEL>1) timer_printf(&ti, "subcyclopclgp:[D_xi_el]");
    2756           0 :     if (DEBUGLEVEL>2) err_printf("z=%Ps\n", z);
    2757           0 :     nz = lg(z)-1;
    2758           0 :     for (k = 1; k <= nz; k++)
    2759           0 :       for (j=1; j<=row; j++)
    2760           0 :         y[(j-1)*row+(i-1)*nz+k-1] = get_y(gel(z, k), gel(vellg, j), d);
    2761           0 :     set_avma(av);
    2762             :   }
    2763           0 : }
    2764             : 
    2765             : static void
    2766          14 : real_MLL1(long *y, GEN K, ulong p, ulong d_pow, GEN velg, GEN vellg, ulong j0)
    2767             : {
    2768          14 :   ulong h = gmael(K, 1, 2)[3], d = upowuu(p, d_pow);
    2769          14 :   GEN elg = gel(velg, 1), ellg = gel(vellg, 1), z;
    2770          14 :   ulong el = elg[1];
    2771             :   pari_timer ti;
    2772             : 
    2773          14 :   if (DEBUGLEVEL>2) timer_start(&ti);
    2774          14 :   z = h < el? D_xi_el_ZX_mul(K, elg, ellg, NULL, d, j0)
    2775          14 :             : D_xi_el_sl(K, elg, ellg, d, j0);
    2776          14 :   if (DEBUGLEVEL>2) timer_printf(&ti, "subcyclopclgp:[D_xi_el]");
    2777          14 :   if (DEBUGLEVEL>2) err_printf("z=%Ps\n", z);
    2778          14 :   y[0] = get_y(gel(z, 1), ellg, d);
    2779          14 : }
    2780             : 
    2781             : static void
    2782           7 : real_MLL(long *y, GEN K, ulong p, ulong d_pow, ulong n,
    2783             :     GEN velg, GEN vellg, GEN vG_K, ulong j0)
    2784             : {
    2785           7 :   ulong i, j, k, d = upowuu(p, d_pow), h = gmael(K, 1, 2)[3];
    2786           7 :   ulong row = lg(vellg)-1;
    2787          28 :   for (j=1; j<=row; j++)
    2788             :   {
    2789          21 :     GEN ellg = gel(vellg, j);
    2790          42 :     for (i=1; i<=n; i++)
    2791             :     {
    2792          21 :       pari_sp av2 = avma;
    2793          21 :       GEN elg = gel(velg, i), z;
    2794          21 :       ulong el = elg[1], nz;
    2795             :       pari_timer ti;
    2796          21 :       if (DEBUGLEVEL>2) timer_start(&ti);
    2797          21 :       z = h < el? D_xi_el_ZX_mul(K, elg, ellg, vG_K, d, j0)
    2798          21 :                 : D_xi_el_sl(K, elg, ellg, d, j0);
    2799          21 :       if (DEBUGLEVEL>2) timer_printf(&ti, "subcyclopclgp:[D_xi_el]");
    2800          21 :       if (DEBUGLEVEL>3) err_printf("z=%Ps\n", z);
    2801          21 :       nz = lg(z)-1;
    2802          84 :       for (k = 1; k <= nz; k++)
    2803          63 :         y[(j-1)*row+(i-1)*nz+k-1] = get_y(gel(z, k), ellg, d);
    2804          21 :       set_avma(av2);
    2805             :     }
    2806             :   }
    2807           7 : }
    2808             : 
    2809             : static long
    2810          21 : use_basis(long d_K, long f) { return (d_K<=10 || (d_K<=30 && f<=5000)); }
    2811             : 
    2812             : static long
    2813           7 : use_factor(ulong f)
    2814           7 : { GEN fa = factoru(f), P = gel(fa, 1); return (P[lg(P)-1]<500); }
    2815             : 
    2816             : /* group structure, destroy gr */
    2817             : static GEN
    2818          63 : get_str(GEN gr)
    2819             : {
    2820          63 :   GEN z = gel(gr,2);
    2821          63 :   long i, j, l = lg(z);
    2822         154 :   for (i = j = 1; i < l; i++)
    2823          91 :     if (lgefint(gel(z, i)) > 2) gel(z,j++) = gel(z,i);
    2824          63 :   setlg(z, j); return z;
    2825             : }
    2826             : 
    2827             : static GEN
    2828          21 : cyc_real_MLL(GEN K, ulong p, long d_pow, long j0, long flag)
    2829             : {
    2830          21 :   ulong d_K = K_get_d(K), f = K_get_f(K), d_chi = K_get_dchi(K);
    2831          21 :   ulong n, n0 = 1, f0, n_el = d_pow, d = upowuu(p, d_pow), rank = n_el*d_chi;
    2832          21 :   GEN velg = const_vec(n_el, NULL), vellg = NULL;
    2833          21 :   GEN oldgr = mkvec2(gen_0, NULL), newgr = mkvec2(gen_0, NULL);
    2834          21 :   long *y0 = (long*)stack_calloc(sizeof(long)*rank*rank);
    2835             : 
    2836          21 :   if (DEBUGLEVEL>1)
    2837           0 :     err_printf("cyc_real_MLL:p=%ld d_pow=%ld deg(K)=%ld cond(K)=%ld g_K=%ld\n",
    2838             :         p, d_pow, d_K, f, K_get_g(K));
    2839          21 :   gel(K, 2) = get_chi(gel(K,1));
    2840          21 :   if (f-1 <= (d_K<<1)) flag |= USE_F;
    2841          21 :   else if (use_basis(d_K, f)) flag |= USE_BASIS;
    2842           7 :   else if (use_factor(f)) flag |= USE_FACTOR;
    2843           0 :   else flag |= USE_GALOIS_POL;
    2844          21 :   if (flag&USE_BASIS) K = vec_append(K, xi_data_basis(K));
    2845           7 :   else if (flag&USE_GALOIS_POL) K = vec_append(K, xi_data_galois(K));
    2846          21 :   f0 = f%p?f:f/p;
    2847          21 :   gel(velg, 1) = next_el_real(K, p, d_pow, mkvecsmall2(1, 1), j0, flag);
    2848          21 :   if (flag&USE_FULL_EL)
    2849             :   {
    2850           0 :     for (n=2; n<=n_el; n++)
    2851           0 :       gel(velg, n) = next_el_real(K, p, d_pow, gel(velg, n+1), j0, flag);
    2852           0 :     n0 = n_el;
    2853             :   }
    2854             : 
    2855          21 :   for (n=n0; n<=n_el; n++) /* loop while structure is unknown */
    2856             :   {
    2857          21 :     pari_sp av2 = avma;
    2858             :     long n_ell, m, M;
    2859             :     GEN y;
    2860             :     pari_timer ti;
    2861          21 :     if (DEBUGLEVEL>2) timer_start(&ti);
    2862          21 :     vellg = set_ell_real(K, velg, n, d_chi, d*d, f0, j0);
    2863          21 :     n_ell = lg(vellg) -1; /* equal to n*d_chi */
    2864          21 :     if (DEBUGLEVEL>2) timer_printf(&ti, "set_ell_real");
    2865          21 :     if (DEBUGLEVEL>3) err_printf("vel=%Ps\nvell=%Ps\n", velg, vellg);
    2866          21 :     if (n_ell==1)
    2867          14 :       real_MLL1(y0, K, p, d_pow, velg, vellg, j0);
    2868             :     else
    2869             :     {
    2870             :       GEN vG_K;
    2871           7 :       if (DEBUGLEVEL>2) timer_start(&ti);
    2872           7 :       vG_K = make_G_K(K, vellg);
    2873           7 :       if (DEBUGLEVEL>2) timer_printf(&ti, "make_G_K");
    2874           7 :       if (lgefint(gmael(vellg, n_ell, 1))<=3 || (flag&SAVE_MEMORY))
    2875           7 :         real_MLL(y0, K, p, d_pow, n, velg, vellg, vG_K, j0);
    2876             :       else
    2877           0 :         real_MLLn(y0, K, p, d_pow, n, velg, vellg, vG_K, j0);
    2878             :     }
    2879          21 :     set_avma(av2);
    2880          21 :     y = ary2mat(y0, n_ell);
    2881          21 :     if (DEBUGLEVEL>3) err_printf("y=%Ps\n", y);
    2882          21 :     y = ZM_snf(y);
    2883          21 :     if (DEBUGLEVEL>3) err_printf("y=%Ps\n", y);
    2884          21 :     y = make_p_part(y, p, d_pow);
    2885          21 :     if (DEBUGLEVEL>3) err_printf("y=%Ps\n", y);
    2886          21 :     newgr = structure_MLL(y, d_pow);
    2887          21 :     if (DEBUGLEVEL>3)
    2888           0 :       err_printf("d_pow=%ld d_chi=%ld old=%Ps new=%Ps\n",d_pow,d_chi,oldgr,newgr);
    2889          21 :     if (equalsi(d_pow*d_chi, gel(newgr, 1))) break;
    2890           0 :     if ((m = find_del_el(&oldgr, newgr, n, n_el, d_chi)))
    2891           0 :     { M = m = delete_el(velg, m); n--; }
    2892             :     else
    2893           0 :     { M = n+1; m = n; }
    2894           0 :     gel(velg, M) = next_el_real(K, p, d_pow, gel(velg, m), j0, flag);
    2895             :   }
    2896          21 :   return get_str(newgr);
    2897             : }
    2898             : 
    2899             : static GEN
    2900           0 : cyc_buch(long dK, GEN p, long d_pow)
    2901             : {
    2902           0 :   GEN z = Buchquad(stoi(dK), 0.0, 0.0, 0), cyc = gel(z,2);
    2903           0 :   long i, l = lg(cyc);
    2904           0 :   if (Z_pval(gel(z,1), p) != d_pow) pari_err_BUG("subcyclopclgp [Buchquad]");
    2905           0 :   for (i = 1; i < l; i++)
    2906             :   {
    2907           0 :     long x = Z_pval(gel(cyc, i), p); if (!x) break;
    2908           0 :     gel(cyc, i) = utoipos(x);
    2909             :   }
    2910           0 :   setlg(cyc, i); return cyc;
    2911             : }
    2912             : 
    2913             : static void
    2914           0 : verbose_output(GEN K, GEN p, long pow, long j)
    2915             : {
    2916           0 :   long d = K_get_d(K), f = K_get_f(K), s = K_get_s(K), d_chi = K_get_dchi(K);
    2917           0 :   err_printf("|A_K_psi|=%Ps^%ld, psi=chi^%ld, d_psi=%ld, %s,\n\
    2918             :     [K:Q]=%ld, [f,H]=[%ld, %Ps]\n",
    2919           0 :     p,pow*d_chi,j,d_chi,s?"real":"imaginary",d,f,zv_to_ZV(gmael3(K,1,1,1)));
    2920           0 : }
    2921             : 
    2922             : static int
    2923       35091 : cyc_real_pre(GEN K, GEN xi, ulong p, ulong j, long el)
    2924             : {
    2925       35091 :   pari_sp av = avma;
    2926       35091 :   ulong i, d_K, x = 1;
    2927       35091 :   GEN e_chi = get_e_chi(K, j, p, &d_K);
    2928             : 
    2929       35091 :   xi++;
    2930     1254827 :   for (i = 0; i < d_K; i++) x = Fl_mul(x, Fl_powu(xi[i], e_chi[i], el), el);
    2931       35091 :   return gc_ulong(av, Fl_powu(x, (el-1)/p, el));
    2932             : }
    2933             : 
    2934             : /* return vec[-1,[],0], vec[0,[],0], vec[1,[1],0], vec[2,[1,1],0] etc */
    2935             : static GEN
    2936       26579 : cyc_real_ss(GEN K, GEN xi, ulong p, long j, long pow, long el, ulong pn, long flag)
    2937             : {
    2938       26579 :   ulong d_chi = K_get_dchi(K);
    2939       26579 :   if (cyc_real_pre(K, xi, pn, j, el) == 1) return NULL; /* not determined */
    2940       21273 :   if (--pow==0) return mkvec3(gen_0, nullvec(), gen_0); /* trivial */
    2941         105 :   if (DEBUGLEVEL) verbose_output(K, utoi(p), pow, j);
    2942         105 :   if (flag&USE_MLL)
    2943             :   {
    2944          21 :     pari_sp av = avma;
    2945          21 :     GEN gr = (K_get_d(K) == 2)? cyc_buch(K_get_f(K), utoi(p), pow)
    2946          21 :                                : cyc_real_MLL(K, p, pow, j, flag);
    2947          21 :     return gerepilecopy(av, mkvec3(utoipos(d_chi * pow), gr, gen_0));
    2948             :   }
    2949          84 :   if (pow==1) return mkvec3(utoi(d_chi), onevec(d_chi), gen_0);
    2950          21 :   return mkvec3(utoi(pow*d_chi), nullvec(), gen_0);
    2951             : }
    2952             : 
    2953             : static GEN
    2954        5145 : cyc_real_ll(GEN K, GEN xi, GEN p, long j, long pow, GEN el, GEN pn, long flag)
    2955             : {
    2956        5145 :   pari_sp av = avma;
    2957        5145 :   ulong i, d_K = K_get_d(K), d_chi = K_get_dchi(K);
    2958        5145 :   GEN e_chi = get_e_chi_ll(K, j, pn), x = gen_1;
    2959             : 
    2960        5145 :   xi++;
    2961      246785 :   for (i = 0; i < d_K; i++)
    2962      241640 :     x = Fp_mul(x, Fp_pow(gel(xi, i), gel(e_chi, i), el), el);
    2963        5145 :   x = Fp_pow(x, diviiexact(subiu(el, 1), pn), el); /* x = x^(el-1)/pn mod el */
    2964        5145 :   set_avma(av); if (equali1(x)) return NULL; /* not determined */
    2965        5145 :   if (--pow==0) return mkvec3(gen_0, nullvec(), gen_0); /* trivial */
    2966           0 :   if (DEBUGLEVEL) verbose_output(K, p, pow, j);
    2967           0 :   if (flag&USE_MLL)
    2968           0 :     pari_err_IMPL(stack_sprintf("flag=%ld for large prime", USE_MLL));
    2969           0 :   if (pow==1) return mkvec3(utoi(d_chi), onevec(d_chi), gen_0);
    2970           0 :   return mkvec3(utoi(pow*d_chi), nullvec(), gen_0);
    2971             : }
    2972             : 
    2973             : /* xi[1+i] = xi^(g^i), 0 <= i <= d-1 */
    2974             : static GEN
    2975       13797 : xi_conj_s(GEN K, ulong el)
    2976             : {
    2977       13797 :   pari_sp av = avma;
    2978       13797 :   GEN  H = K_get_H(K);
    2979       13797 :   long d = K_get_d(K), f = K_get_f(K), h = K_get_h(K), g = K_get_g(K);
    2980       13797 :   long i, gi = 1, z = Fl_powu(pgener_Fl(el), (el-1)/f, el);
    2981       13797 :   GEN vz = Fl_powers(z, f, el)+1, xi = cgetg(d+1, t_VECSMALL);
    2982             : 
    2983      372953 :   for (i=1; i<=d; i++)
    2984             :   {
    2985      359156 :     long j, x = 1;
    2986   170777852 :     for (j=1; j<=h; j++)
    2987   170418696 :       x = Fl_mul(x, vz[Fl_mul(H[j], gi, f)]-1, el);
    2988      359156 :     xi[i] = x;
    2989      359156 :     gi = Fl_mul(gi, g, f);
    2990             :   }
    2991       13797 :   return gerepilecopy(av, xi);
    2992             : }
    2993             : 
    2994             : static GEN
    2995        1673 : xi_conj_l(GEN K, GEN el)
    2996             : {
    2997        1673 :   pari_sp av = avma;
    2998        1673 :   GEN  H = K_get_H(K);
    2999        1673 :   long d = K_get_d(K), f = K_get_f(K), h = K_get_h(K), g = K_get_g(K);
    3000        1673 :   long i, gi = 1;
    3001        1673 :   GEN z = Fp_pow(pgener_Fp(el), diviuexact(subiu(el, 1), f), el);
    3002        1673 :   GEN vz = Fp_powers(z, f, el)+1, xi = cgetg(d+1, t_VEC);
    3003             : 
    3004       70462 :   for (i=1; i<=d; i++)
    3005             :   {
    3006             :     long j;
    3007       68789 :     GEN x = gen_1;
    3008     6756043 :     for (j=1; j<=h; j++)
    3009     6687254 :       x = Fp_mul(x, subiu(gel(vz, Fl_mul(H[j], gi, f)), 1), el);
    3010       68789 :     gel(xi, i) = x;
    3011       68789 :     gi = Fl_mul(gi, g, f);
    3012             :   }
    3013        1673 :   return gerepilecopy(av, xi);
    3014             : }
    3015             : 
    3016             : static GEN
    3017       10164 : pclgp_cyc_real(GEN K, GEN p, long max_pow, long flag)
    3018             : {
    3019       10164 :   const long NUM_EL = 20;
    3020       10164 :   GEN C = gel(K, 5);
    3021       10164 :   long f_K = K_get_f(K), n_conj = K_get_nconj(K);
    3022       10164 :   long i, pow, n_el, n_done = 0;
    3023       10164 :   GEN gr = nullvec(), Done = const_vecsmall(n_conj, 0), xi;
    3024       10164 :   long first = 1;
    3025             : 
    3026       10290 :   for (pow=1; pow<=max_pow; pow++)
    3027             :   {
    3028       10290 :     GEN pn = powiu(p, pow), fpn = muliu(pn, f_K), el = addiu(fpn, 1);
    3029      109753 :     for (n_el = 0; n_el < NUM_EL; el = addii(el, fpn))
    3030             :     {
    3031             :       ulong uel;
    3032      109627 :       if (!BPSW_psp(el)) continue;
    3033       15470 :       n_el++; uel = itou_or_0(el);
    3034       15470 :       if (uel)
    3035             :       {
    3036       13797 :         xi = xi_conj_s(K, uel);
    3037       13797 :         if (first && n_conj > 10) /* mark trivial chi-part */
    3038             :         {
    3039        8715 :           for (i = 1; i <= n_conj; i++)
    3040             :           {
    3041        8512 :             if (cyc_real_pre(K, xi, p[2], C[i], uel) == 1) continue;
    3042        8260 :             Done[i] = 1;
    3043        8260 :             if (++n_done == n_conj) return gr;
    3044             :           }
    3045         203 :           first = 0; continue;
    3046             :         }
    3047             :       }
    3048             :       else
    3049        1673 :         xi = xi_conj_l(K, el);
    3050       55132 :       for (i = 1; i <= n_conj; i++)
    3051             :       {
    3052             :         GEN z;
    3053       50029 :         if (Done[i]) continue;
    3054       31724 :         if (uel)
    3055       26579 :           z = cyc_real_ss(K, xi, p[2], C[i], pow, uel, itou(pn), flag);
    3056             :         else
    3057        5145 :           z = cyc_real_ll(K, xi, p, C[i], pow, el, pn, flag);
    3058       31724 :         if (!z) continue;
    3059       26418 :         Done[i] = 1;
    3060       26418 :         if (!isintzero(gel(z, 1))) gr = vec_append(gr, z);
    3061       26418 :         if (++n_done == n_conj) return gr;
    3062             :       }
    3063             :     }
    3064             :   }
    3065           0 :   pari_err_BUG("pclgp_cyc_real: max_pow is not enough");
    3066             :   return NULL; /*LCOV_EXCL_LINE*/
    3067             : }
    3068             : 
    3069             : /* return (el, g_el) */
    3070             : static GEN
    3071          56 : next_el_imag(GEN elg, long f)
    3072             : {
    3073          56 :   long el = elg[1];
    3074          56 :   if (odd(f)) f<<=1;
    3075         140 :   while (!uisprime(el+=f));
    3076          56 :   return mkvecsmall2(el, pgener_Fl(el));
    3077             : }
    3078             : 
    3079             : /* return (ell, g_ell) */
    3080             : static GEN
    3081          70 : next_ell_imag(GEN ellg, GEN df0l)
    3082             : {
    3083          70 :   GEN ell = gel(ellg, 1);
    3084         770 :   while (!BPSW_psp(ell = addii(ell, df0l)));
    3085          70 :   return mkvec2(ell, pgener_Fp(ell));
    3086             : }
    3087             : 
    3088             : static GEN
    3089          56 : set_ell_imag(GEN velg, long n, long d_chi, GEN df0)
    3090             : {
    3091          56 :   long i, n_ell = n*d_chi;
    3092          56 :   GEN z = cgetg(n_ell + 1, t_VEC);
    3093          56 :   GEN df0l = shifti(df0, 1), ellg = mkvec2(gen_1, gen_1);
    3094         126 :   for (i=1; i<=n; i++) df0l = muliu(df0l, gel(velg, i)[1]);
    3095         126 :   for (i=1; i<=n_ell; i++) ellg = gel(z, i)= next_ell_imag(ellg, df0l);
    3096          56 :   return z;
    3097             : }
    3098             : 
    3099             : /* U(X)=u(x)+u(X)*X^f+...+f(X)*X^((m-1)f) or u(x)-u(X)*X^f+...
    3100             :  * U(X)V(X)=u(X)V(X)(1+X^f+...+X^((m-1)f))
    3101             :  *         =w_0+w_1*X+...+w_{f+el-3}*X^(f+el-3)
    3102             :  * w_i (1 <= i <= f+el-2) are needed.
    3103             :  * w_{f+el-2}=0 if el-1 == f.
    3104             :  * W_i = w_i + w_{i+el-1} (1 <= i <= f-1). */
    3105             : static GEN
    3106          85 : gauss_Flx_mul(ulong f, GEN elg, GEN ellg)
    3107             : {
    3108          85 :   pari_sp av = avma;
    3109          85 :   ulong el = elg[1], g_el= elg[2];
    3110          85 :   ulong el_1 = el-1, f2 = f<<1, lv = el_1, lu = f, m = el_1/f;
    3111          85 :   ulong ell = itou(gel(ellg, 1)), g_ell = itou(gel(ellg, 2));
    3112          85 :   ulong z_2f = Fl_powu(g_ell, (ell - 1) / f2, ell);
    3113          85 :   ulong z_el = Fl_powu(g_ell, (ell - 1) / el, ell);
    3114             :   ulong i, i2, gi;
    3115          85 :   GEN W = cgetg(f+1, t_VECSMALL), vz_2f, vz_el;
    3116          85 :   GEN u = cgetg(lu+2, t_VECSMALL), v = cgetg(lv+2, t_VECSMALL), w0;
    3117             : 
    3118          85 :   u[1] = evalsigne(1);
    3119          85 :   v[1] = evalsigne(1);
    3120          85 :   vz_2f = Fl_powers(z_2f, f2-1, ell);
    3121          85 :   vz_el = Fl_powers(z_el, el_1, ell);
    3122      528459 :   for (i=i2=0; i<lu; i++)
    3123             :   {
    3124             :     long j2; /* i2=(i*i)%f2, gi=g_el^i */
    3125      528374 :     j2 = i2?f2-i2:i2;
    3126      528374 :     u[2+i] = vz_2f[1+j2];
    3127      528374 :     if ((i2+=i+i+1)>=f2) i2-=f2; /* same as i2%=f2 */
    3128             :   }
    3129     1600537 :   for (gi=1,i=i2=0; i<lv; i++)
    3130             :   {
    3131     1600452 :     v[2+i] = Fl_mul(vz_2f[1+i2], vz_el[1+gi], ell);
    3132     1600452 :     gi = Fl_mul(gi, g_el, el);
    3133     1600452 :     if ((i2+=i+i+1)>=f2) i2%=f2; /* i2-=f2 does not work */
    3134             :   }
    3135          85 :   w0 = Flx_mul(u, v, ell) + 1;
    3136          85 :   if (m==1)
    3137             :   { /* el_1=f */
    3138           0 :     for (i=1; i<f; i++) W[i] = Fl_add(w0[i], w0[i+lv], ell);
    3139           0 :     W[f] = w0[f];
    3140             :   }
    3141             :   else
    3142             :   {
    3143          85 :     ulong start = 1+f, end = f+el-1;
    3144          85 :     GEN w = cgetg(end+1, t_VECSMALL);
    3145     2128826 :     for (i=1; i<end; i++) w[i] = w0[i];
    3146          85 :     w[end] = 0;
    3147         360 :     for (i=1; i<m; i++, start+=f)
    3148         550 :       w = both_odd(f,i)? Flv_shift_sub(w, w0, ell, start, end)
    3149         275 :                        : Flv_shift_add(w, w0, ell, start, end);
    3150      528459 :     for (i=0; i<f; i++) W[1+i] = Fl_add(w[1+i], w[1+i+lv], ell);
    3151             :   }
    3152      528374 :   for (i=i2=1; i<f; i++)
    3153             :   {
    3154      528289 :     W[i]=Fl_mul(W[1+i], vz_2f[1+i2], ell);
    3155      528289 :     if ((i2+=i+i+1)>=f2) i2%=f2;
    3156             :   }
    3157             :   /* W[r]=tau_{LL}^{sigma_r}, 1<= r <= f-1 */
    3158          85 :   return gerepilecopy(av, Flv_to_ZV(W));
    3159             : }
    3160             : 
    3161             : static GEN
    3162          90 : gauss_ZX_mul(ulong f, GEN elg, GEN ellg)
    3163             : {
    3164          90 :   pari_sp av = avma, av2;
    3165             :   ulong el, g_el, el_1, f2, lv, lu, m, i, i2, gi;
    3166          90 :   GEN  ell = gel(ellg, 1), g_ell, ell_1, z_2f, z_el, W, vz_2f, vz_el, u, v, w0;
    3167             : 
    3168          90 :   if (lgefint(ell) == 3) return gauss_Flx_mul(f, elg, ellg);
    3169           5 :   g_ell = gel(ellg, 2); ell_1 = subiu(ell, 1);
    3170           5 :   el = elg[1]; g_el = elg[2]; el_1 = el-1;
    3171           5 :   f2 = f<<1; lv=el_1; lu=f; m=el_1/f;
    3172           5 :   z_2f = Fp_pow(g_ell, diviuexact(ell_1, f2), ell);
    3173           5 :   vz_2f = Fp_powers(z_2f, f2-1, ell);
    3174           5 :   W = cgetg(f+1, t_VEC);
    3175           5 :   av2 = avma;
    3176           5 :   z_el = Fp_pow(g_ell, diviuexact(ell_1, el), ell);
    3177           5 :   vz_el = Fp_powers(z_el, el_1, ell);
    3178           5 :   u = cgetg(lu+2, t_POL); u[1] = evalsigne(1) | evalvarn(0);
    3179           5 :   v = cgetg(lv+2, t_POL); v[1] = evalsigne(1) | evalvarn(0);
    3180       35264 :   for (gi=1,i=i2=0; i<lu; i++)
    3181             :   {
    3182             :     long j2; /* i2=(i*i)%f2, gi=g_el^i */
    3183       35259 :     j2 = i2?f2-i2:i2;
    3184       35259 :     gel(u, 2+i) = gel(vz_2f, 1+j2);
    3185       35259 :     if ((i2+=i+i+1)>=f2) i2-=f2;
    3186             :   }
    3187       82787 :   for (gi=1,i=i2=0; i<lv; i++)
    3188             :   {
    3189       82782 :     gel(v, 2+i) = Fp_mul(gel(vz_2f, 1+i2), gel(vz_el, 1+gi), ell);
    3190       82782 :     gi = Fl_mul(gi, g_el, el);
    3191       82782 :     if ((i2+=i+i+1)>=f2) i2%=f2;
    3192             :   }
    3193           5 :   w0 = gerepileupto(av2, FpX_mul(u, v, ell)) + 1; av2 = avma;
    3194           5 :   if (m==1)
    3195             :   {
    3196           0 :     for (i=1; i < f; i++) gel(W,i) = Fp_add(gel(w0, i), gel(w0, i+lv), ell);
    3197           0 :     gel(W, f) = gel(w0, f);
    3198             :   }
    3199             :   else
    3200             :   {
    3201           5 :     ulong start = 1+f, end = f+el-1;
    3202           5 :     GEN w = cgetg(end+1, t_VEC);
    3203      118041 :     for (i=1; i<end; i++) gel(w, i) = gel(w0, i);
    3204           5 :     gel(w, end) = gen_0;
    3205          15 :     for (i=1; i<m; i++, start+=f)
    3206             :     {
    3207          13 :       w = both_odd(f,i)? FpV_shift_sub(w, w0, ell, start, end)
    3208          10 :                        : FpV_shift_add(w, w0, ell, start, end);
    3209          10 :       if ((i & 7) == 0) w = gerepilecopy(av2, w);
    3210             :     }
    3211       35264 :     for (i = 1; i <= f; i++) gel(W, i) = addii(gel(w, i), gel(w, i+lv));
    3212             :   }
    3213       35259 :   for (i = i2 = 1; i < f; i++)
    3214             :   {
    3215       35254 :     gel(W, i) = Fp_mul(gel(W, 1+i), gel(vz_2f, 1+i2), ell);
    3216       35254 :     if ((i2+=i+i+1) >= f2) i2 %= f2;
    3217             :   }
    3218           5 :   return gerepilecopy(av, W);  /* W[r]=tau_{LL}^{sigma_r}, 1<= r <= f-1 */
    3219             : }
    3220             : 
    3221             : /* fast but consumes memory */
    3222             : static GEN
    3223           4 : gauss_el_vell(ulong f, GEN elg, GEN vellg, GEN vz_2f)
    3224             : {
    3225           4 :   pari_sp av = avma, av2;
    3226           4 :   ulong el = elg[1], g_el = elg[2], el_1 = el-1;
    3227           4 :   ulong lv=el_1, f2=f<<1, lu=f, m=el_1/f;
    3228           4 :   GEN W = cgetg(f+1, t_VEC), vz_el, u, v, w0, M;
    3229             :   ulong i, i2, gi;
    3230             : 
    3231           4 :   av2 = avma;
    3232           4 :   vz_el = vz_vell(el, vellg, &M);
    3233           4 :   u = cgetg(lu+2, t_POL); u[1] = evalsigne(1) | evalvarn(0);
    3234           4 :   v = cgetg(lv+2, t_POL); v[1] = evalsigne(1) | evalvarn(0);
    3235       25554 :   for (i=i2=0; i<lu; i++)
    3236             :   {
    3237             :     long j2; /* i2=(i*i)%f2, gi=g_el^i */
    3238       25550 :     j2 = i2?f2-i2:i2;
    3239       25550 :     gel(u, 2+i) = gel(vz_2f, 1+j2);
    3240       25550 :     if ((i2+=i+i+1)>=f2) i2%=f2;
    3241             :   }
    3242       86874 :   for (gi=1,i=i2=0; i<lv; i++)
    3243             :   {
    3244       86870 :     gel(v, 2+i) = Fp_mul(gel(vz_2f, 1+i2), gel(vz_el, 1+gi), M);
    3245       86870 :     gi = Fl_mul(gi, g_el, el);
    3246       86870 :     if ((i2+=i+i+1)>=f2) i2%=f2;
    3247             :   }
    3248           4 :   M = gclone(M);
    3249           4 :   w0 = gerepileupto(av2, FpX_mul(u, v, M)) + 1;
    3250           4 :   u = M; M = icopy(M); gunclone(u);
    3251           4 :   av2 = avma;
    3252           4 :   if (m==1)
    3253             :   { /* el_1=f */
    3254           0 :     for (i=1; i < f; i++) gel(W,i) = Fp_add(gel(w0, i), gel(w0, i+lv), M);
    3255           0 :     gel(W, f) = gel(w0, f);
    3256             :   }
    3257             :   else
    3258             :   {
    3259           4 :     ulong start = 1+f, end = f+el-1;
    3260           4 :     GEN w = cgetg(end+1, t_VEC);
    3261      112420 :     for (i=1; i<end; i++) gel(w, i) = gel(w0, i);
    3262           4 :     gel(w, end) = gen_0;
    3263          19 :     for (i=1; i<m; i++, start+=f)
    3264             :     {
    3265          22 :       w = both_odd(f,i)? FpV_shift_sub(w, w0, M, start, end)
    3266          15 :                        : FpV_shift_add(w, w0, M, start, end);
    3267          15 :       if ((i & 7) == 0) w = gerepilecopy(av2, w);
    3268             :     }
    3269       25554 :     for (i = 1; i <= f; i++) gel(W, i) = Fp_add(gel(w, i), gel(w, i+lv), M);
    3270             :   }
    3271       25550 :   for (i = i2 = 1; i < f; i++)
    3272             :   {
    3273       25546 :     gel(W, i) = Fp_mul(gel(W, 1+i), gel(vz_2f, 1+i2), M);
    3274       25546 :     if ((i2+=i+i+1) >= f2) i2 %= f2;
    3275             :   }
    3276           4 :   return gerepilecopy(av, W);  /* W[r]=tau_{LL}^{sigma_r}, 1<= r <= f-1 */
    3277             : }
    3278             : 
    3279             : static GEN
    3280          94 : norm_chi(GEN K, GEN TAU, ulong p, long d_pow, GEN ell, long j0)
    3281             : {
    3282          94 :   pari_sp av = avma;
    3283          94 :   GEN H = K_get_H(K);
    3284          94 :   ulong d_K, f_K = K_get_f(K), h = K_get_h(K), g_K = K_get_g(K);
    3285          94 :   ulong i, j, gi, pd = upowuu(p, d_pow), d_chi = K_get_dchi(K);
    3286          94 :   GEN z = const_vec(d_chi, gen_1);
    3287          94 :   GEN e_chi = get_e_chi(K, j0, pd, &d_K);
    3288             : 
    3289        1420 :   for (gi=1, i=0; i<d_K; i++)
    3290             :   {
    3291        1326 :     GEN y = gen_1;
    3292      230862 :     for (j=1; j<=h; j++)
    3293      229536 :       y = Fp_mul(y, gel(TAU, Fl_mul(gi, H[j], f_K)), ell);
    3294        1326 :     gi = Fl_mul(gi, g_K, f_K);
    3295        2652 :     for (j=1; j<=d_chi; j++)
    3296             :     {
    3297        1326 :       GEN y2 = Fp_powu(y, e_chi[(i+j-1)%d_K], ell);
    3298        1326 :       gel(z, j) = Fp_mul(gel(z, j), y2, ell);
    3299             :     }
    3300             :   }
    3301          94 :   return gerepilecopy(av, z);
    3302             : }
    3303             : 
    3304             : static void
    3305           2 : imag_MLLn(long *y, GEN K, ulong p, long d_pow, long n,
    3306             :     GEN velg, GEN vellg, long j0)
    3307             : {
    3308           2 :   long f = K_get_f(K), d = upowuu(p, d_pow), row = lg(vellg)-1, i, j, k, nz;
    3309           2 :   GEN g, z, M, vz_2f = vz_vell(f << 1, vellg, &M);
    3310           6 :   for (i=1; i<=n; i++)
    3311             :   {
    3312           4 :     pari_sp av = avma;
    3313           4 :     GEN elg = gel(velg, i);
    3314           4 :     if (DEBUGLEVEL>1) err_printf("(f,el-1)=(%ld,%ld*%ld)\n", f,(elg[1]-1)/f,f);
    3315           4 :     g = gauss_el_vell(f, elg, vellg, vz_2f);
    3316           4 :     z = norm_chi(K, g, p, d_pow, M, j0);
    3317           4 :     nz = lg(z)-1;
    3318           8 :     for (k = 1; k <= nz; k++)
    3319          12 :       for (j = 1; j <= row; j++)
    3320           8 :         y[(j-1)*row+(i-1)*nz+k-1] = get_y(gel(z, k), gel(vellg, j), d);
    3321           4 :     set_avma(av);
    3322             :   }
    3323           2 : }
    3324             : 
    3325             : static void
    3326          42 : imag_MLL1(long *y, GEN K, ulong p, long d_pow, GEN velg, GEN vellg, long j0)
    3327             : {
    3328          42 :   long f = K_get_f(K), d = upowuu(p, d_pow);
    3329          42 :   GEN elg = gel(velg, 1), ellg = gel(vellg, 1), ell = gel(ellg, 1), g, z;
    3330             : 
    3331          42 :   if (DEBUGLEVEL>1) err_printf("(f,el-1)=(%ld,%ld*%ld)\n", f, (elg[1]-1)/f, f);
    3332          42 :   g = gauss_ZX_mul(f, elg, ellg);
    3333          42 :   z = norm_chi(K, g, p, d_pow, ell, j0);
    3334          42 :   y[0] = get_y(gel(z, 1), ellg, d);
    3335          42 : }
    3336             : 
    3337             : static void
    3338          12 : imag_MLL(long *y, GEN K, ulong p, long d_pow, long n, GEN velg, GEN vellg,
    3339             :     long j0)
    3340             : {
    3341          12 :   pari_sp av = avma;
    3342          12 :   long i, j, f = K_get_f(K), d = upowuu(p, d_pow), row = lg(vellg)-1;
    3343             : 
    3344          36 :   for (j=1; j<=row; j++)
    3345             :   {
    3346          24 :     GEN ellg = gel(vellg, j), ell = gel(ellg, 1);
    3347          72 :     for (i=1; i<=n; i++)
    3348             :     {
    3349          48 :       GEN elg = gel(velg, i), g, z;
    3350             :       ulong k, nz;
    3351          48 :       if (DEBUGLEVEL>1) err_printf("(f,el-1)=(%ld,%ld*%ld)\n",f,(elg[1]-1)/f,f);
    3352          48 :       g = gauss_ZX_mul(f, elg, ellg);
    3353          48 :       z = norm_chi(K, g, p, d_pow, ell, j0);
    3354          48 :       nz = lg(z)-1;
    3355          96 :       for (k = 1; k <= nz; k++)
    3356          48 :         y[(j-1)*row+(i-1)*nz+k-1] = get_y(gel(z, k), ellg, d);
    3357          48 :       set_avma(av);
    3358             :     }
    3359             :   }
    3360          12 : }
    3361             : 
    3362             : /* return an upper bound >= 0 if one was found, otherwise return -1.
    3363             :  * set chi-part to be (1) if chi is Teichmuller character.
    3364             :  * B_{1,omega^(-1)} is not p-adic integer. */
    3365             : static GEN
    3366          42 : cyc_imag_MLL(GEN K, ulong p, long d_pow, long j, long flag)
    3367             : {
    3368          42 :   long f = K_get_f(K), d_chi = K_get_dchi(K);
    3369          42 :   long n, n0 = 1, n_el = d_pow, d = upowuu(p, d_pow), rank = n_el*d_chi;
    3370          42 :   GEN df0, velg = const_vec(n_el, NULL), vellg = NULL;
    3371          42 :   GEN oldgr = mkvec2(gen_0, NULL), newgr = mkvec2(gen_0, NULL);
    3372          42 :   long *y0 = (long*)stack_calloc(sizeof(long)*rank*rank);
    3373             : 
    3374          42 :   if (DEBUGLEVEL>1)
    3375           0 :     err_printf("cyc_imag_MLL:p=%ld d_pow=%ld deg(K)=%ld cond(K)=%ld avma=%ld\n",
    3376             :         p, d_pow, K_get_d(K), f, avma);
    3377          42 :   df0 = muluu(d, f%p?f:f/p);
    3378          42 :   gel(velg, 1) = next_el_imag(mkvecsmall2(1, 1), f);
    3379          42 :   if (flag&USE_FULL_EL)
    3380             :   {
    3381           0 :     for (n=2; n<=n_el; n++) gel(velg, n) = next_el_imag(gel(velg, n-1), f);
    3382           0 :     n0 = n_el;
    3383             :   }
    3384          56 :   for (n=n0; n<=n_el; n++) /* loop while structure is unknown */
    3385             :   {
    3386          56 :     pari_sp av2 = avma;
    3387             :     pari_timer ti;
    3388             :     long n_ell, m, M;
    3389             :     GEN y;
    3390          56 :     vellg = set_ell_imag(velg, n, d_chi, df0);
    3391          56 :     n_ell = lg(vellg)-1;  /* equal to n*d_chi */
    3392          56 :     if (DEBUGLEVEL>2) err_printf("velg=%Ps\nvellg=%Ps\n", velg, vellg);
    3393          56 :     if (DEBUGLEVEL>2) timer_start(&ti);
    3394          56 :     if (n_ell==1)
    3395          42 :       imag_MLL1(y0, K, p, d_pow, velg, vellg, j);
    3396          14 :     else if (lgefint(gmael(vellg, n, 1))<=3 || (flag&SAVE_MEMORY))
    3397          12 :       imag_MLL(y0, K, p, d_pow, n, velg, vellg, j);
    3398             :     else
    3399           2 :       imag_MLLn(y0, K, p, d_pow, n, velg, vellg, j);
    3400          56 :     set_avma(av2);
    3401          56 :     if (DEBUGLEVEL>2) timer_printf(&ti, "gauss sum");
    3402          56 :     y = ary2mat(y0, n_ell);
    3403          56 :     if (DEBUGLEVEL>3) err_printf("y=%Ps\n", y);
    3404          56 :     y = ZM_snf(y);
    3405          56 :     if (DEBUGLEVEL>3) err_printf("y=%Ps\n", y);
    3406          56 :     y = make_p_part(y, p, d_pow);
    3407          56 :     if (DEBUGLEVEL>3) err_printf("y=%Ps\n", y);
    3408          56 :     newgr = structure_MLL(y, d_pow);
    3409          56 :     if (DEBUGLEVEL>3)
    3410           0 :       err_printf("d_pow=%ld d_chi=%ld old=%Ps new=%Ps\n",d_pow,d_chi,oldgr,newgr);
    3411          56 :     if (equalsi(d_pow*d_chi, gel(newgr, 1))) break;
    3412          14 :     if ((m = find_del_el(&oldgr, newgr, n, n_el, d_chi)))
    3413           0 :     { M = m = delete_el(velg, m); n--; }
    3414             :     else
    3415          14 :     { M = n+1; m = n; }
    3416          14 :     gel(velg, M) = next_el_imag(gel(velg, m), f);
    3417             :   }
    3418          42 :   return get_str(newgr);
    3419             : }
    3420             : 
    3421             : /* When |A_psi|=p^e, A_psi=(p^e1,...,p^er) (psi=chi^j),
    3422             :  *  return vec[e, [e1, ... ,er], 1].
    3423             :  * If gr str is not determined, return vec[e, [], 1].
    3424             :  * If |A_chi|=1, return vec[0, [], 1].
    3425             :  * If |A_chi|=p, return vec[1, [1], 1].
    3426             :  * If e is not determined, return vec[-1, [], 1].
    3427             :  * If psi is Teichmuller, return vec[0, [], 1].
    3428             :  * B_{1,omega^(-1)} is not p-adic integer. */
    3429             : static GEN
    3430       26334 : cyc_imag(GEN K, GEN B, GEN p, long j, GEN powp, long flag)
    3431             : {
    3432       26334 :   pari_sp av = avma;
    3433       26334 :   GEN MinPol = gel(K, 3), Chi = gel(K, 2), B1, B2, gr;
    3434       26334 :   long x, d_K = K_get_d(K), f_K = K_get_f(K), d_chi = K_get_dchi(K);
    3435             : 
    3436       26334 :   if (f_K == d_K+1 && equaliu(p, f_K) && j == 1) /* Teichmuller */
    3437          77 :     return mkvec3(gen_0, nullvec(), gen_1);
    3438       26257 :   B1 = FpX_rem(ZX_ber_conj(B, j, d_K), MinPol, powp);
    3439       26257 :   B2 = FpX_rem(ZX_ber_den(Chi, j, d_K), MinPol, powp);
    3440       26257 :   if (degpol(B1)<0 || degpol(B2)<0)
    3441             :   {
    3442           0 :     set_avma(av);
    3443           0 :     return mkvec3(gen_m1, nullvec(), gen_1); /* B=0(mod p^pow) */
    3444             :   }
    3445       26257 :   x = ZX_pval(B1, p) - ZX_pval(B2, p);
    3446       26257 :   set_avma(av);
    3447       26257 :   if (x<0) pari_err_BUG("subcyclopclgp [Bernoulli number]");
    3448       26257 :   if (DEBUGLEVEL && x) verbose_output(K, p, x, j);
    3449       26257 :   if (x==0) return mkvec3(gen_0, nullvec(), gen_1); /* trivial */
    3450         588 :   if (x==1) return mkvec3(utoi(d_chi), onevec(d_chi), gen_1);
    3451         140 :   if ((flag&USE_MLL)==0) return mkvec3(utoi(x*d_chi), nullvec(), gen_1);
    3452          42 :   gr = d_K == 2? cyc_buch(-f_K, p, x): cyc_imag_MLL(K, itou(p), x, j, flag);
    3453          42 :   return gerepilecopy(av, mkvec3(utoipos(d_chi * x), gr, gen_1));
    3454             : }
    3455             : 
    3456             : /* handle representatives of all injective characters, d_chi=[Q_p(zeta_d):Q_p],
    3457             :  * d=d_K */
    3458             : static GEN
    3459       10080 : pclgp_cyc_imag(GEN K, GEN p, long start_pow, long max_pow, long flag)
    3460             : {
    3461       10080 :   GEN C = gel(K, 5), Chi = gel(K, 2);
    3462       10080 :   long n_conj = K_get_nconj(K), d_K = K_get_d(K), f_K = K_get_f(K);
    3463       10080 :   long i, pow, n_done = 0;
    3464       10080 :   GEN gr = nullvec(), Done = const_vecsmall(n_conj, 0);
    3465       10080 :   GEN B = zx_ber_num(Chi, f_K, d_K), B_num;
    3466             : 
    3467       10080 :   if (lgefint(p)==3 && n_conj>10) /* mark trivial chi-part by pre-calculation */
    3468             :   {
    3469         595 :     ulong up = itou(p);
    3470         595 :     GEN minpol = ZX_to_Flx(gel(K, 3), up);
    3471        7350 :     for (i=1; i<=n_conj; i++)
    3472             :     {
    3473        7168 :       pari_sp av = avma;
    3474             :       long degB;
    3475        7168 :       B_num = Flx_rem(Flx_ber_conj(B, C[i], d_K, up), minpol, up);
    3476        7168 :       degB = degpol(B_num);
    3477        7168 :       set_avma(av);
    3478        7168 :       if (degB<0) continue;
    3479        6937 :       Done[i] = 1;
    3480        6937 :       if (++n_done == n_conj) return gr;
    3481             :     }
    3482             :   }
    3483        9667 :   for (pow = start_pow; pow<=max_pow; pow++)
    3484             :   {
    3485        9667 :     GEN powp = powiu(p, pow);
    3486       27503 :     for (i = 1; i <= n_conj; i++)
    3487             :     {
    3488             :       GEN z;
    3489       27503 :       if (Done[i]) continue;
    3490       26334 :       z = cyc_imag(K, B, p, C[i], powp, flag);
    3491       26334 :       if (equalim1(gel(z, 1))) continue;
    3492       26334 :       Done[i] = 1;
    3493       26334 :       if (!isintzero(gel(z, 1))) gr = vec_append(gr, z);
    3494       26334 :       if (++n_done == n_conj) return gr;
    3495             :     }
    3496             :   }
    3497           0 :   pari_err_BUG("pclgp_cyc_imag: max_pow is not enough");
    3498             :   return NULL; /*LCOV_EXCL_LINE*/
    3499             : }
    3500             : 
    3501             : static GEN
    3502         392 : gather_part(GEN g, long sgn)
    3503             : {
    3504         392 :   long i, j, l = lg(g), ord = 0, flag = 1;
    3505         392 :   GEN z2 = cgetg(l, t_VEC);
    3506        1778 :   for (i = j = 1; i < l; i++)
    3507             :   {
    3508        1386 :     GEN t = gel(g,i);
    3509        1386 :     if (equaliu(gel(t, 3), sgn))
    3510             :     {
    3511         693 :       ord += itou(gel(t, 1));
    3512         693 :       if (lg(gel(t, 2)) == 1) flag = 0;
    3513         693 :       gel(z2, j++) = gel(t, 2);
    3514             :     }
    3515             :   }
    3516         392 :   if (flag==0 || ord==0) z2 = nullvec();
    3517             :   else
    3518             :   {
    3519         126 :     setlg(z2, j); z2 = shallowconcat1(z2);
    3520         126 :     ZV_sort_inplace(z2); vecreverse_inplace(z2);
    3521             :   }
    3522         392 :   return mkvec2(utoi(ord), z2);
    3523             : }
    3524             : 
    3525             : #ifdef DEBUG
    3526             : static void
    3527             : handling(GEN K)
    3528             : {
    3529             :   long d_K = K_get_d(K), f_K = K_get_f(K), s_K = K_get_s(K), g_K = K_get_g(K);
    3530             :   long d_chi = K_get_dchi(K);
    3531             :   err_printf("  handling %s cyclic subfield K,\
    3532             :       deg(K)=%ld, cond(K)=%ld g_K=%ld d_chi=%ld H=%Ps\n",
    3533             :       s_K? "a real": "an imaginary",d_K,f_K,g_K,d_chi,zv_to_ZV(gmael3(K,1,1,1)));
    3534             : }
    3535             : #endif
    3536             : 
    3537             : /* HH a t_VECSMALL listing group generators
    3538             :  * Aoki and Fukuda, LNCS vol.4076 (2006), 56-74. */
    3539             : static GEN
    3540         161 : pclgp(GEN p0, long f, GEN HH, long degF, long flag)
    3541             : {
    3542             :   long start_pow, max_pow, ip, lp, i, n_f;
    3543         161 :   GEN vH1, z, vData, cycGH, vp = typ(p0) == t_INT? mkvec(p0): p0;
    3544             : 
    3545         161 :   vH1 = GHinit(f, HH, &cycGH); n_f = lg(vH1)-1;
    3546             : #ifdef DEBUG
    3547             :   err_printf("F is %s, deg(F)=%ld, ", srh_1(HH)? "real": "imaginary", degF);
    3548             :   err_printf("cond(F)=%ld, G(F/Q)=%Ps\n",f, cycGH);
    3549             :   err_printf("F has %ld cyclic subfield%s except for Q.\n", n_f,n_f>1?"s":"");
    3550             : #endif
    3551             : 
    3552         161 :   lp = lg(vp); z = cgetg(lp, t_MAT);
    3553         357 :   for (ip = 1; ip < lp; ip++)
    3554             :   {
    3555         196 :     pari_sp av = avma;
    3556         196 :     long n_sub=0, n_chi=0;
    3557         196 :     GEN gr=nullvec(), p = gel(vp, ip), zi;
    3558             :     /* find conductor e of cyclic subfield K and set the subgroup HE of (Z/eZ)^*
    3559             :      * corresponding to K */
    3560         196 :     set_p_f(p, f, &start_pow, &max_pow);
    3561         196 :     vData = const_vec(degF, NULL);
    3562             : 
    3563       16982 :     for (i=1; i<=n_f; i++) /* prescan. set Teichmuller */
    3564             :     {
    3565       16863 :       GEN H1 = gel(vH1, i);
    3566       16863 :       long d_K = _get_d(H1), f_K = _get_f(H1), g_K = _get_g(H1);
    3567             : 
    3568       16863 :       if (f_K == d_K+1 && equaliu(p, f_K)) /* found K=Q(zeta_p) */
    3569             :       {
    3570             :         pari_timer ti;
    3571          77 :         GEN pnmax = powiu(p, max_pow), vNewton, C, MinPol;
    3572          77 :         long d_chi = 1, n_conj = eulerphiu(d_K);
    3573          77 :         ulong pmodd = umodiu(p, d_K);
    3574             : 
    3575          77 :         C = set_C(pmodd, d_K, d_chi, n_conj);
    3576          77 :         MinPol = set_minpol_teich(g_K, p, max_pow);
    3577          77 :         if (DEBUGLEVEL>3) timer_start(&ti);
    3578          77 :         vNewton = FpX_Newton(MinPol, d_K+1, pnmax);
    3579          77 :         if (DEBUGLEVEL>3)
    3580           0 :           timer_printf(&ti, "FpX_Newton: teich: %ld %ld", degpol(MinPol), d_K);
    3581          77 :         gel(vData, d_K) = mkvec4(MinPol, vNewton, C,
    3582             :                                  mkvecsmall2(d_chi, n_conj));
    3583          77 :         break;
    3584             :       }
    3585             :     }
    3586             : 
    3587       20440 :     for (i=1; i<=n_f; i++) /* handle all cyclic K */
    3588             :     {
    3589       20244 :       GEN H1 = gel(vH1, i), K, z1, Chi;
    3590       20244 :       long d_K = _get_d(H1), s_K = _get_s(H1);
    3591             :       pari_sp av2;
    3592             : 
    3593       20244 :       if ((flag&SKIP_PROPER) && degF != d_K) continue;
    3594       20244 :       if (!gel(vData, d_K))
    3595             :       {
    3596             :         pari_timer ti;
    3597         819 :         GEN pnmax = powiu(p, max_pow), vNewton, C, MinPol;
    3598         819 :         ulong pmodd = umodiu(p, d_K);
    3599         819 :         long d_chi = order_f_x(d_K, pmodd), n_conj = eulerphiu(d_K)/d_chi;
    3600             : 
    3601         819 :         C = set_C(pmodd, d_K, d_chi, n_conj);
    3602         819 :         MinPol = set_minpol(d_K, p, max_pow, n_conj);
    3603         819 :         if (DEBUGLEVEL>3) timer_start(&ti);
    3604             :         /* vNewton[2+i] = vNewton[2+i+d_K]. We need vNewton[2+i] for
    3605             :          * 0 <= i < d_K. But vNewton[2+d_K-1] may be 0 and will be deleted.
    3606             :          * So we need vNewton[2+d_K] not to delete vNewton[2+d_K-1]. */
    3607         819 :         vNewton = FpX_Newton(MinPol, d_K+1, pnmax);
    3608         819 :         if (DEBUGLEVEL>3)
    3609           0 :           timer_printf(&ti, "FpX_Newton: %ld %ld", degpol(MinPol), d_K);
    3610         819 :         gel(vData, d_K) = mkvec4(MinPol, vNewton, C,
    3611             :                                  mkvecsmall2(d_chi, n_conj));
    3612             :       }
    3613       20244 :       av2 = avma;
    3614       20244 :       Chi = s_K? NULL: get_chi(H1);
    3615       20244 :       K = shallowconcat(mkvec2(H1, Chi), gel(vData, d_K));
    3616             : #ifdef DEBUG
    3617             :       handling(K);
    3618             : #endif
    3619       20244 :       if (s_K && !(flag&NO_PLUS_PART))
    3620       10164 :         z1 = pclgp_cyc_real(K, p, max_pow, flag);
    3621       10080 :       else if (!s_K && !(flag&NO_MINUS_PART))
    3622       10080 :         z1 = pclgp_cyc_imag(K, p, start_pow, max_pow, flag);
    3623           0 :       else { set_avma(av2); continue; }
    3624       20244 :       n_sub++; n_chi += gmael(vData, d_K, 4)[2]; /* += n_conj */
    3625       20244 :       if (lg(z1) == 1) set_avma(av2);
    3626         658 :       else gr = gerepilecopy(av2, shallowconcat(gr, z1));
    3627             :     }
    3628         196 :     zi = mkcol(p);
    3629         196 :     zi = vec_append(zi, (flag&NO_PLUS_PART)?nullvec():gather_part(gr, 0));
    3630         196 :     zi = vec_append(zi, (flag&NO_MINUS_PART)?nullvec():gather_part(gr, 1));
    3631         196 :     zi = shallowconcat(zi, mkcol3(cycGH, utoi(n_sub), utoi(n_chi)));
    3632         196 :     gel(z, ip) = gerepilecopy(av, zi);
    3633             :   }
    3634         161 :   return typ(p0) == t_INT? shallowtrans(gel(z,1)): shallowtrans(z);
    3635             : }
    3636             : 
    3637             : static GEN
    3638         413 : reduce_gcd(GEN x1, GEN x2)
    3639             : {
    3640         413 :   GEN d = gcdii(x1, x2);
    3641         413 :   x1 = diviiexact(x1, d);
    3642         413 :   x2 = diviiexact(x2, d);
    3643         413 :   return mkvec2(x1, x2);
    3644             : }
    3645             : 
    3646             : /* norm of x0 (= pol of zeta_d with deg <= d-1) by g of order n
    3647             :  * x0^{1+g+g^2+...+g^(n-1)} */
    3648             : static GEN
    3649          49 : ber_norm_cyc(GEN x0, long g, long n, long d)
    3650             : {
    3651          49 :   pari_sp av = avma;
    3652          49 :   long i, ei, di, fi = 0, l = ulogint(n, 2);
    3653          49 :   GEN xi = x0;
    3654          49 :   ei = 1L << l; di = n / ei;
    3655         203 :   for (i = 1; i <= l; i++)
    3656             :   {
    3657         154 :     if (odd(di)) fi += ei;
    3658         154 :     ei = 1L << (l-i); di = n / ei;
    3659         154 :     xi = ZX_mod_Xnm1(ZX_mul(xi, ber_conj(xi, Fl_powu(g, ei, d), d)), d);
    3660         154 :     if (odd(di))
    3661          42 :       xi = ZX_mod_Xnm1(ZX_mul(xi, ber_conj(x0, Fl_powu(g, fi, d), d)), d);
    3662             :   }
    3663          49 :   return gerepilecopy(av, xi);
    3664             : }
    3665             : 
    3666             : /* x0 a ZX of deg < d */
    3667             : static GEN
    3668          21 : ber_norm_by_cyc(GEN x0, long d, GEN MinPol)
    3669             : {
    3670          21 :   pari_sp av=avma;
    3671          21 :   GEN x = x0, z = znstar(utoi(d)), cyc = gel(z, 2), gen = gel(z, 3);
    3672          21 :   long i, l = lg(cyc);
    3673             :   pari_timer ti;
    3674             : 
    3675          21 :   if (DEBUGLEVEL>1) timer_start(&ti);
    3676          70 :   for (i = 1; i < l; i++)
    3677          49 :     x = ber_norm_cyc(x, itou(gmael(gen, i, 2)), itou(gel(cyc, i)), d);
    3678          21 :   if (DEBUGLEVEL>1) timer_printf(&ti, "ber_norm_by_cyc [ber_norm_cyc]");
    3679          21 :   x = ZX_rem(x, MinPol);  /* slow */
    3680          21 :   if (DEBUGLEVEL>1) timer_printf(&ti, "ber_norm_by_cyc [ZX_rem]");
    3681          21 :   if (lg(x) != 3) pari_err_BUG("subcyclohminus [norm of Bernoulli number]");
    3682          21 :   return gerepilecopy(av, gel(x, 2));
    3683             : }
    3684             : 
    3685             : /* MinPol = polcyclo(d_K, 0).
    3686             :  * MinPol = fac*cofac (mod p).
    3687             :  * B is zv.
    3688             :  * K : H1, MinPol, [fac, cofac], C, [d_chi, n_conj] */
    3689             : static long
    3690          98 : ber_norm_by_val(GEN K, GEN B, GEN p)
    3691             : {
    3692          98 :   pari_sp av = avma;
    3693          98 :   GEN MinPol = gel(K, 2), C = gel(K, 4);
    3694          98 :   GEN vfac = gel(K, 3), fac = gel(vfac, 1), cofac = gel(vfac, 2);
    3695          98 :   long d_chi = K_get_dchi(K), n_conj = K_get_nconj(K), d_K = K_get_d(K);
    3696          98 :   long i, r, n_done = 0, x = 0, dcofac = degpol(cofac);
    3697             :   GEN pr, Done;
    3698             : 
    3699          98 :   Done = const_vecsmall(n_conj, 0);
    3700          98 :   if (lgefint(p)==3)
    3701             :   { /* mark trivial chi-part by pre-calculation */
    3702          98 :     ulong up = itou(p);
    3703          98 :     GEN facs = ZX_to_Flx(fac, up);
    3704         196 :     for (i = 1; i <= n_conj; i++)
    3705             :     {
    3706          98 :       pari_sp av2 = avma;
    3707          98 :       GEN B_conj = Flx_rem(Flx_ber_conj(B, C[i], d_K, up), facs, up);
    3708          98 :       long degB = degpol(B_conj);
    3709          98 :       set_avma(av2); if (degB < 0) continue;
    3710           0 :       Done[i] = 1; if (++n_done == n_conj) return gc_long(av, x);
    3711             :     }
    3712             :   }
    3713             :   else
    3714             :   {
    3715           0 :     for (i = 1; i <= n_conj; i++)
    3716             :     {
    3717           0 :       pari_sp av2 = avma;
    3718           0 :       GEN B_conj = FpX_rem(FpX_ber_conj(B, C[i], d_K, p), fac, p);
    3719           0 :       long degB = degpol(B_conj);
    3720           0 :       set_avma(av2); if (degB < 0) continue;
    3721           0 :       Done[i] = 1; if (++n_done == n_conj) return gc_long(av, x);
    3722             :     }
    3723             :   }
    3724         252 :   for (pr = p, r = 2; r; r <<= 1)
    3725             :   {
    3726             :     GEN polr;
    3727         252 :     pr = sqri(pr); /* p^r */
    3728         252 :     polr = (dcofac==0)? FpX_red(MinPol, pr)
    3729         252 :                       : gel(ZpX_liftfact(MinPol, vfac, pr, p, r), 1);
    3730         406 :     for (i = 1; i <= n_conj; i++)
    3731             :     {
    3732         252 :       pari_sp av2 = avma;
    3733             :       GEN B_conj;
    3734             :       long degB;
    3735         252 :       if (Done[i]) continue;
    3736         252 :       B_conj = FpX_rem(FpX_ber_conj(B, C[i], d_K, pr), polr, pr);
    3737         252 :       degB = degpol(B_conj);
    3738         252 :       set_avma(av2); if (degB < 0) continue;
    3739          98 :       x += d_chi * ZX_pval(B_conj, p);
    3740          98 :       Done[i] = 1; if (++n_done == n_conj) return gc_long(av, x);
    3741             :     }
    3742             :   }
    3743             :   pari_err_BUG("ber_norm_by_val"); return 0;/*LCOV_EXCL_LINE*/
    3744             : }
    3745             : 
    3746             : /* n > 2, p = odd prime not dividing n, e > 0, pe = p^e; d = n*p^e
    3747             :  * return generators of the subgroup H of (Z/dZ)^* corresponding to
    3748             :  * Q(zeta_{p^e}): H = {1<=a<=d | gcd(a,n)=1, a=1(mod p^e)} */
    3749             : static GEN
    3750           0 : znstar_subgr(ulong n, ulong pe, ulong d)
    3751             : {
    3752           0 :   GEN z = znstar(utoi(n)), g = gel(z, 3), G;
    3753           0 :   long i, l = lg(g);
    3754           0 :   G = cgetg(l, t_VECSMALL);
    3755           0 :   for (i=1; i<l; i++) G[i] = u_chinese_coprime(itou(gmael(g,i,2)), 1, n, pe, d);
    3756           0 :   return mkvec2(gel(z,2), G);
    3757             : }
    3758             : 
    3759             : /* K is a cyclic extension of degree n*p^e (n>=4 is even).
    3760             :  * x a ZX of deg < n*p^e. */
    3761             : static long
    3762           0 : ber_norm_with_val(GEN x, long n, ulong p, ulong e)
    3763             : {
    3764           0 :   pari_sp av = avma;
    3765           0 :   long i, j, r, degx, pe = upowuu(p, e), d = n*pe;
    3766           0 :   GEN z, gr, gen, y = cgetg(pe+2, t_POL), MinPol = polcyclo(n, 0);
    3767           0 :   y[1] = evalsigne(1) | evalvarn(0);
    3768           0 :   z = znstar_subgr(n, pe, d);
    3769           0 :   gr = gel(z, 1); gen = gel(z, 2); r = lg(gr)-1;
    3770           0 :   for (i=1; i<=r; i++)
    3771           0 :     x = ber_norm_cyc(x, itou(gel(gen, i)), itou(gel(gr, i)), d);
    3772           0 :   degx = degpol(x);
    3773           0 :   for (j=0; j<pe; j++)
    3774             :   {
    3775           0 :     GEN t = pol_zero(n), z;
    3776           0 :     long a = j; /* a=i*pe+j */
    3777           0 :     for (i=0; i<n; i++)
    3778             :     {
    3779           0 :       if (a>degx) break;
    3780           0 :       gel(t, 2+a%n) = gel(x, 2+a);
    3781           0 :       a += pe;
    3782             :     }
    3783           0 :     z = ZX_rem(ZX_renormalize(t, 2+n), MinPol);
    3784           0 :     if (degpol(z)<0) gel(y, 2+j) = gen_0;
    3785           0 :     else if (degpol(z)==0) gel(y, 2+j) = gel(z, 2);
    3786           0 :     else pari_err_BUG("ber_norm_subgr");
    3787             :   }
    3788           0 :   y = ZX_renormalize(y, pe+2);
    3789           0 :   if (e>1) y = ZX_rem(y, polcyclo(pe, 0));
    3790           0 :   return gc_long(av,  ZX_p_val(y, p, e));
    3791             : }
    3792             : 
    3793             : /* K is a cyclic extension of degree 2*p^e. x a ZX of deg < 2*p^e. In most
    3794             :  * cases, deg(x)=2*p^e-1. But deg(x) can be any value in [0, 2*p^e-1]. */
    3795             : static long
    3796         301 : ber_norm_with_val2(GEN x, ulong p, ulong e)
    3797             : {
    3798         301 :   pari_sp av = avma;
    3799         301 :   long i, d = degpol(x), pe = upowuu(p, e);
    3800         301 :   GEN y = pol_zero(pe);
    3801         301 :   if (d == 2*pe-1)
    3802             :   {
    3803       38416 :     for (i = 0; i < pe; i++)
    3804       76230 :       gel(y, 2+i) = odd(i)? subii(gel(x, 2+i+pe), gel(x, 2+i))
    3805       38115 :                           : subii(gel(x, 2+i), gel(x, 2+i+pe));
    3806             :   }
    3807             :   else
    3808             :   {
    3809           0 :     for (i = 0; i < pe && i <= d; i++)
    3810           0 :       gel(y, 2+i) = odd(i)? negi(gel(x, 2+i)): gel(x, 2+i);
    3811           0 :     for (i = pe; i <= d; i++)
    3812           0 :       gel(y, 2+i-pe) = odd(i)? subii(gel(y, 2+i-pe), gel(x, 2+i))
    3813           0 :                              : addii(gel(y, 2+i-pe), gel(x, 2+i));
    3814             :   }
    3815         301 :   y = ZX_renormalize(y, 2+pe);
    3816         301 :   if (e > 1) y = ZX_rem(y, polcyclo(pe, 0));
    3817         301 :   return gc_long(av, ZX_p_val(y, p, e));
    3818             : }
    3819             : 
    3820             : /* K : H1, MinPol, [fac, cofac], C, [d_chi, n_conj] */
    3821             : static GEN
    3822         812 : ber_cyc5(GEN K, GEN p)
    3823             : {
    3824         812 :   pari_sp av = avma;
    3825         812 :   GEN MinPol = gel(K, 2), H = K_get_H(K);
    3826         812 :   long d = K_get_d(K), f = K_get_f(K), h = K_get_h(K), g = K_get_g(K);
    3827         812 :   GEN x, x1, x2, y, B = const_vecsmall(d+1, 0);
    3828         812 :   long i, j, gi, e, f2 = f>>1, dMinPol = degpol(MinPol), chi2 = -1, *B2 = B+2;
    3829             : 
    3830             :   /* get_chi inlined here to save memory */
    3831    18111989 :   for (j=1; j<=h; j++) /* i = 0 */
    3832             :   {
    3833    18111177 :     if (H[j] == 2) chi2 = 0;
    3834    18111177 :     if (H[j] <= f2) B2[0]++; /* Chi[H[j]] = 0 */
    3835             :   }
    3836       97314 :   for (i = 1, gi = g; i < d; i++)
    3837             :   {
    3838    93017085 :     for (j=1; j<=h; j++)
    3839             :     {
    3840    92920583 :       long t = Fl_mul(gi, H[j], f); /* Chi[t] = i */
    3841    92920583 :       if (t == 2) chi2 = i;
    3842    92920583 :       if (t <= f2) B2[i]++;
    3843             :     }
    3844       96502 :     gi = Fl_mul(gi, g, f);
    3845             :   }
    3846         812 :   y = zx_to_ZX(zx_renormalize(B, d+2));
    3847             : 
    3848         812 :   if (p)
    3849             :   {
    3850             :     ulong n;
    3851         399 :     e = u_pvalrem(d, p, &n);
    3852         399 :     if (e == 0)
    3853          98 :       x1 = utoi(ber_norm_by_val(K, B, p));
    3854         301 :     else if (n > 2)
    3855           0 :       x1 = utoi(ber_norm_with_val(y, n, itou(p), e));
    3856             :     else
    3857         301 :       x1 = utoi(ber_norm_with_val2(y, itou(p), e));
    3858             :   }
    3859             :   else
    3860             :   {
    3861         413 :     if (dMinPol > 100)
    3862          21 :       x1 = ber_norm_by_cyc(y, d, MinPol);
    3863             :     else
    3864         392 :       x1 = ZX_resultant(MinPol, ZX_rem(y, MinPol));
    3865             :   }
    3866             : 
    3867         812 :   if (chi2 < 0) /* chi2 = Chi[2] */
    3868           0 :     x2 = shifti(gen_1, 2*dMinPol);
    3869         812 :   else if (chi2 == 0)
    3870          21 :     x2 = shifti(gen_1, dMinPol);
    3871             :   else
    3872             :   {
    3873         791 :     long e = d/ugcd(chi2, d);
    3874         791 :     x2 = powiu(polcyclo_eval(e, gen_2), eulerphiu(d)/eulerphiu(e));
    3875         791 :     x2 = shifti(x2, dMinPol);
    3876             :   }
    3877         812 :   if (p) x = stoi(itou(x1)-Z_pval(x2, p)); else x = reduce_gcd(x1, x2);
    3878         812 :   return gerepilecopy(av, x);
    3879             : }
    3880             : 
    3881             : /*  Hirabayashi-Yoshino, Manuscripta Math. vol.60, 423-436 (1988), Theorem 1
    3882             :  *
    3883             :  *  F is a subfield of Q(zeta_f)
    3884             :  *  f=p^a => Q=1
    3885             :  *  If F=Q(zeta_f), Q=1 <=> f=p^a
    3886             :  *  If f=4*p^a, p^a*q^b (p,q are odd primes), Q=2 <=> [Q(zeta_f):F] is odd */
    3887             : static long
    3888          21 : unit_index(ulong d, ulong f)
    3889             : {
    3890             :   ulong r, d_f;
    3891          21 :   GEN fa = factoru(f), P = gel(fa, 1), E = gel(fa, 2); r = lg(P)-1;
    3892          21 :   if (r==1) return 1;  /* f=P^a */
    3893           7 :   d_f = eulerphiu_fact(fa);
    3894           7 :   if (d==d_f) return 2;  /* F=Q(zeta_f) */
    3895           0 :   if (r==2 && ((P[1]==2 && E[1]==2) || P[1]>2)) return odd(d_f/d)?2:1;
    3896           0 :   return 0;
    3897             : }
    3898             : 
    3899             : /* Compute relative class number h of the subfield K of Q(zeta_f)
    3900             :  * corresponding to the subgroup HH of (Z/fZ)^*.
    3901             :  * If p!=NULL, then return valuation(h,p). */
    3902             : static GEN
    3903         119 : rel_class_num(long f, GEN HH, long degF, GEN p)
    3904             : {
    3905             :   long i, n_f, W, Q;
    3906         119 :   GEN vH1, vData, x, z = gen_1, z1 = gen_0, z2 = mkvec2(gen_1, gen_1);
    3907             : 
    3908         119 :   vH1 = GHinit(f, HH, NULL); n_f = lg(vH1)-1;
    3909         119 :   vData = const_vec(degF, NULL);
    3910        1652 :   for (i=1; i<=n_f; i++)
    3911             :   {
    3912        1533 :     GEN H1 = gel(vH1, i), K;
    3913        1533 :     long d_K = _get_d(H1), s = _get_s(H1);
    3914             : 
    3915        1533 :     if (s) continue;  /* F is real */
    3916             : #ifdef DEBUG
    3917             :     err_printf("  handling %s cyclic subfield K, deg(K)=%ld, cond(K)=%ld\n",
    3918             :         s? "a real": "an imaginary", d_K, _get_f(H1));
    3919             : #endif
    3920         812 :     if (!gel(vData, d_K))
    3921             :     {
    3922             :       GEN C, MinPol, fac, cofac;
    3923             :       ulong d_chi, n_conj;
    3924         497 :       MinPol = polcyclo(d_K,0);
    3925         497 :       if (p && umodui(d_K, p))
    3926          98 :       {
    3927          98 :         ulong pmodd = umodiu(p, d_K);
    3928          98 :         GEN MinPol_p = FpX_red(MinPol, p);
    3929          98 :         d_chi = order_f_x(d_K, pmodd);
    3930          98 :         n_conj = eulerphiu(d_K)/d_chi;
    3931          98 :         if (n_conj==1) fac = MinPol_p;  /* polcyclo(d_K) is irred mod p */
    3932           0 :         else fac = FpX_one_cyclo(d_K, p);
    3933          98 :         cofac = FpX_div(MinPol_p, fac, p);
    3934          98 :         C = set_C(pmodd, d_K, d_chi, n_conj);
    3935             :       }
    3936             :       else
    3937             :       {
    3938         399 :         fac = cofac = C = NULL;
    3939         399 :         d_chi = n_conj = 0;
    3940             :       }
    3941         497 :       gel(vData, d_K) = mkvec5(MinPol, mkvec2(fac, cofac), C,
    3942             :                                NULL, mkvecsmall2(d_chi, n_conj));
    3943             :     }
    3944         812 :     K = vec_prepend(gel(vData, d_K), H1);
    3945         812 :     z = ber_cyc5(K, p);
    3946         812 :     if (p) z1 = addii(z1, z);
    3947             :     else
    3948             :     {
    3949         413 :       gel(z2, 1) = mulii(gel(z2, 1), gel(z, 1));
    3950         413 :       gel(z2, 2) = mulii(gel(z2, 2), gel(z, 2));
    3951             :     }
    3952             :   }
    3953         119 :   W = root_of_1(f, HH);
    3954         119 :   if (p) return addiu(z1, z_pval(W, p));
    3955          21 :   Q = unit_index(degF, f);
    3956          21 :   x = dvmdii(muliu(gel(z2,1), 2 * W), gel(z2,2), &z1);
    3957          21 :   if (signe(z1)) pari_err_BUG("subcyclohminus [norm of Bernoulli number]");
    3958          21 :   if (!Q && mpodd(x)) Q = 2; /* FIXME: can this happen ? */
    3959          21 :   if (Q == 1) x = shifti(x, -1);
    3960          21 :   return mkvec2(x, utoi(Q));
    3961             : }
    3962             : 
    3963             : static void
    3964         336 : checkp(const char *fun, long degF, GEN p)
    3965             : {
    3966         336 :   if (!BPSW_psp(p)) pari_err_PRIME(fun, p);
    3967         329 :   if (equaliu(p, 2)) pari_err_DOMAIN(fun,"p","=", gen_2, p);
    3968         315 :   if (degF && dvdsi(degF, p)) errpdiv(fun, p, degF);
    3969         301 : }
    3970             : 
    3971             : /* if flag is set, handle quadratic fields specially (don't set H) */
    3972             : static long
    3973         441 : subcyclo_init(const char *fun, GEN FH, long *pdegF, GEN *pH, long flag)
    3974             : {
    3975         441 :   long f = 0, degF = 2;
    3976         441 :   GEN F = NULL, H = NULL;
    3977         441 :   if (typ(FH) == t_POL)
    3978             :   {
    3979          77 :     degF = degpol(FH);
    3980          77 :     if (degF < 1 || !RgX_is_ZX(FH)) pari_err_TYPE(fun, FH);
    3981          77 :     if (flag && degF == 2)
    3982             :     {
    3983          56 :       F = coredisc(ZX_disc(FH));
    3984          56 :       if (is_bigint(F))
    3985           0 :         pari_err_IMPL(stack_sprintf("conductor f > %lu in %s", ULONG_MAX, fun));
    3986          56 :       f = itos(F); if (f == 1) degF = 1;
    3987             :     }
    3988             :     else
    3989             :     {
    3990          21 :       GEN z, bnf = Buchall(pol_x(fetch_var()), 0, DEFAULTPREC);
    3991          21 :       z = rnfconductor(bnf, FH); H = gel(z,3);
    3992          21 :       f = subcyclo_nH(fun, gel(z,2), &H);
    3993          21 :       delete_var();
    3994          21 :       H = znstar_generate(f, H); /* group elements */
    3995             :     }
    3996             :   }
    3997             :   else
    3998             :   {
    3999         364 :     long l = lg(FH), fH;
    4000         364 :     if (typ(FH) == t_INT) F = FH;
    4001         273 :     else if (typ(FH) == t_VEC && (l == 2 || l == 3))
    4002             :     {
    4003         273 :       F = gel(FH, 1);
    4004         273 :       if (l == 3) H = gel(FH, 2);
    4005             :     }
    4006           0 :     else pari_err_TYPE(fun, FH);
    4007         364 :     f = subcyclo_nH(fun, F, &H);
    4008         350 :     H = znstar_generate(f, H); /* group elements */
    4009         350 :     fH = znstar_conductor(H);
    4010         350 :     if (fH == 1) degF = 1;
    4011             :     else
    4012             :     {
    4013         350 :       if (fH != f) { H = znstar_reduce_modulus(H, fH); f = fH; }
    4014         350 :       degF = eulerphiu(f) / zv_prod(gel(H, 2));
    4015             :     }
    4016             :   }
    4017         427 :   *pH = H; *pdegF = degF; return f;
    4018             : }
    4019             : 
    4020             : GEN
    4021         210 : subcyclopclgp(GEN FH, GEN p, long flag)
    4022             : {
    4023         210 :   pari_sp av = avma;
    4024             :   GEN H;
    4025         210 :   long degF, f = subcyclo_init("subcyclopclgp", FH, &degF, &H, 0);
    4026         203 :   if (typ(p) == t_VEC)
    4027             :   {
    4028          28 :     long i, l = lg(p);
    4029          77 :     for (i = 1; i < l; i++) checkp("subcyclopclgp", degF, gel(p, i));
    4030          14 :     if (f == 1) { set_avma(av); return const_vec(l-1, nullvec()); }
    4031             :   }
    4032             :   else
    4033             :   {
    4034         175 :     checkp("subcyclopclgp", degF, p);
    4035         154 :     if (f == 1) { set_avma(av); return nullvec(); }
    4036             :   }
    4037         168 :   if (flag >= USE_BASIS) pari_err_FLAG("subcyclopclgp");
    4038         161 :   return gerepilecopy(av, pclgp(p, f, H, degF, flag));
    4039             : }
    4040             : 
    4041             : static GEN
    4042         126 : subcycloiwasawa_i(GEN FH, GEN P, long n)
    4043             : {
    4044             :   long B, p, f, degF;
    4045             :   GEN H;
    4046         126 :   const char *fun = "subcycloiwasawa";
    4047             : 
    4048         126 :   if (typ(P) != t_INT) pari_err_TYPE(fun, P);
    4049         126 :   if (n < 0) pari_err_DOMAIN(fun, "n", "<", gen_0, stoi(n));
    4050         126 :   B = 1L << (BITS_IN_LONG/4);
    4051         126 :   if (is_bigint(P) || cmpiu(P, B) > 0)
    4052           7 :     pari_err_IMPL(stack_sprintf("prime p > %ld in %s", B, fun));
    4053         119 :   p = itos(P);
    4054         119 :   if (p <= 1 || !uisprime(p)) pari_err_PRIME(fun, P);
    4055         119 :   if (!upowuu(p, n+1))
    4056           7 :     pari_err_IMPL(stack_sprintf("p^n > 2^%ld in %s", BITS_IN_LONG, fun));
    4057         112 :   f = subcyclo_init(fun, FH, &degF, &H, 1);
    4058         105 :   if (degF == 1) return NULL;
    4059         105 :   if (degF == 2)
    4060             :   {
    4061          56 :     long m = ((f & 3) == 0)? f / 4: f;
    4062          56 :     if (H && !srh_1(H)) m = -m;
    4063          56 :     if (!n) return quadlambda(p, m);
    4064          28 :     return m < 0? imagquadstkpol(p, m, n): realquadstkpol(p, m, n);
    4065             :   }
    4066          49 :   if (p == 2) pari_err_DOMAIN(fun, "p", "=", gen_2, gen_2);
    4067          49 :   if (srh_1(H)) return NULL;
    4068          49 :   if (degF % p == 0) errpdiv("abeliwasawa", P, degF);
    4069          49 :   return abeliwasawa(p, f, H, degF, n);
    4070             : }
    4071             : GEN
    4072         126 : subcycloiwasawa(GEN FH, GEN P, long n)
    4073             : {
    4074         126 :   pari_sp av = avma;
    4075         126 :   GEN z = subcycloiwasawa_i(FH, P, n);
    4076         105 :   if (!z) { set_avma(av); return n? nullvec(): mkvec(gen_0); }
    4077         105 :   return gerepilecopy(av, z);
    4078             : }
    4079             : 
    4080             : GEN
    4081         119 : subcyclohminus(GEN FH, GEN P)
    4082             : {
    4083         119 :   const char *fun = "subcyclohminus";
    4084         119 :   pari_sp av = avma;
    4085             :   GEN H;
    4086         119 :   long degF, f = subcyclo_init(fun, FH, &degF, &H, 0);
    4087         119 :   if (P)
    4088             :   {
    4089          98 :     if (typ(P) != t_INT) pari_err_TYPE(fun, P);
    4090          98 :     if (isintzero(P)) P = NULL; else checkp(fun, 0, P);
    4091             :   }
    4092         119 :   if (degF == 1 ||  srh_1(H) == 1) return gen_1;
    4093         119 :   return gerepilecopy(av, rel_class_num(f, H, degF, P));
    4094             : }

Generated by: LCOV version 1.16