Code coverage tests

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

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

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

LCOV - code coverage report
Current view: top level - basemath - subcyclo.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.8.0 lcov report (development 19611-73a567d) Lines: 524 556 94.2 %
Date: 2016-09-27 05:54:40 Functions: 38 39 97.4 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2000  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             : /*************************************************************************/
      18             : /**                                                                     **/
      19             : /**              Routines for handling subgroups of (Z/nZ)^*            **/
      20             : /**              without requiring discrete logarithms.                 **/
      21             : /**                                                                     **/
      22             : /*************************************************************************/
      23             : /* Subgroups are [gen,ord,bits] where
      24             :  * gen is a vecsmall of generators
      25             :  * ord is theirs relative orders
      26             :  * bits is a bit vector of the elements, of length(n). */
      27             : 
      28             : /*The algorithm is similar to testpermutation*/
      29             : static void
      30        4074 : znstar_partial_coset_func(long n, GEN H, void (*func)(void *data,long c)
      31             :     , void *data, long d, long c)
      32             : {
      33             :   GEN gen, ord, cache;
      34             :   long i, j, card;
      35             : 
      36        8148 :   if (!d) { (*func)(data,c); return; }
      37             : 
      38        2086 :   cache = const_vecsmall(d,c);
      39        2086 :   (*func)(data,c);  /* AFTER cache: may contain gerepileupto statement */
      40        2086 :   gen = gel(H,1);
      41        2086 :   ord = gel(H,2);
      42        2086 :   card = ord[1]; for (i = 2; i <= d; i++) card *= ord[i];
      43       66661 :   for(i=1; i<card; i++)
      44             :   {
      45       64575 :     long k, m = i;
      46       64575 :     for(j=1; j<d && m%ord[j]==0 ;j++) m /= ord[j];
      47       64575 :     cache[j] = Fl_mul(cache[j],gen[j],n);
      48       64575 :     for (k=1; k<j; k++) cache[k] = cache[j];
      49       64575 :     (*func)(data, cache[j]);
      50             :   }
      51             : }
      52             : 
      53             : static void
      54         728 : znstar_coset_func(long n, GEN H, void (*func)(void *data,long c)
      55             :     , void *data, long c)
      56             : {
      57         728 :   znstar_partial_coset_func(n, H, func,data, lg(gel(H,1))-1, c);
      58         728 : }
      59             : 
      60             : /* Add the element of the bitvec of the coset c modulo the subgroup of H
      61             :  * generated by the first d generators to the bitvec bits.*/
      62             : 
      63             : static void
      64        3346 : znstar_partial_coset_bits_inplace(long n, GEN H, GEN bits, long d, long c)
      65             : {
      66        3346 :   pari_sp av = avma;
      67        3346 :   znstar_partial_coset_func(n,H, (void (*)(void *,long)) &F2v_set,
      68             :       (void *) bits, d, c);
      69        3346 :   avma = av;
      70        3346 : }
      71             : 
      72             : static void
      73         364 : znstar_coset_bits_inplace(long n, GEN H, GEN bits, long c)
      74             : {
      75         364 :   znstar_partial_coset_bits_inplace(n, H, bits, lg(gel(H,1))-1, c);
      76         364 : }
      77             : 
      78             : static GEN
      79        2982 : znstar_partial_coset_bits(long n, GEN H, long d, long c)
      80             : {
      81        2982 :   GEN bits = zero_F2v(n);
      82        2982 :   znstar_partial_coset_bits_inplace(n,H,bits,d,c);
      83        2982 :   return bits;
      84             : }
      85             : 
      86             : /* Compute the bitvec of the elements of the subgroup of H generated by the
      87             :  * first d generators.*/
      88             : static GEN
      89        2982 : znstar_partial_bits(long n, GEN H, long d)
      90             : {
      91        2982 :   return znstar_partial_coset_bits(n, H, d, 1);
      92             : }
      93             : 
      94             : /* Compute the bitvec of the elements of H. */
      95             : GEN
      96           0 : znstar_bits(long n, GEN H)
      97             : {
      98           0 :   return znstar_partial_bits(n,H,lg(gel(H,1))-1);
      99             : }
     100             : 
     101             : /* Compute the subgroup of (Z/nZ)^* generated by the elements of
     102             :  * the vecsmall V */
     103             : GEN
     104        1778 : znstar_generate(long n, GEN V)
     105             : {
     106        1778 :   pari_sp av = avma;
     107        1778 :   GEN gen = cgetg(lg(V),t_VECSMALL);
     108        1778 :   GEN ord = cgetg(lg(V),t_VECSMALL), res = mkvec2(gen,ord);
     109        1778 :   GEN bits = znstar_partial_bits(n,NULL,0);
     110        1778 :   long i, r = 0;
     111        4235 :   for(i=1; i<lg(V); i++)
     112             :   {
     113        2457 :     ulong v = uel(V,i), g = v;
     114        2457 :     long o = 0;
     115        2457 :     while (!F2v_coeff(bits, (long)g)) { g = Fl_mul(g, v, (ulong)n); o++; }
     116        2457 :     if (!o) continue;
     117        1204 :     r++;
     118        1204 :     gen[r] = v;
     119        1204 :     ord[r] = o+1;
     120        1204 :     cgiv(bits); bits = znstar_partial_bits(n,res,r);
     121             :   }
     122        1778 :   setlg(gen,r+1);
     123        1778 :   setlg(ord,r+1); return gerepilecopy(av, mkvec3(gen,ord,bits));
     124             : }
     125             : 
     126             : static ulong
     127        1834 : znstar_order(GEN H) { return zv_prod(gel(H,2)); }
     128             : 
     129             : /* Return the lists of element of H.
     130             :  * This can be implemented with znstar_coset_func instead. */
     131             : GEN
     132        1596 : znstar_elts(long n, GEN H)
     133             : {
     134        1596 :   long card = znstar_order(H);
     135        1596 :   GEN gen = gel(H,1), ord = gel(H,2);
     136        1596 :   GEN sg = cgetg(1 + card, t_VECSMALL);
     137             :   long k, j, l;
     138        1596 :   sg[1] = 1;
     139        2576 :   for (j = 1, l = 1; j < lg(gen); j++)
     140             :   {
     141         980 :     long c = l * (ord[j]-1);
     142         980 :     for (k = 1; k <= c; k++) sg[++l] = Fl_mul(sg[k], gen[j], n);
     143             :   }
     144        1596 :   vecsmall_sort(sg); return sg;
     145             : }
     146             : 
     147             : /* Take a znstar H and n dividing the modulus of H.
     148             :  * Output H reduced to modulus n */
     149             : GEN
     150          35 : znstar_reduce_modulus(GEN H, long n)
     151             : {
     152          35 :   pari_sp ltop=avma;
     153          35 :   GEN gen=cgetg(lgcols(H),t_VECSMALL);
     154             :   long i;
     155         119 :   for(i=1; i < lg(gen); i++)
     156          84 :     gen[i] = mael(H,1,i)%n;
     157          35 :   return gerepileupto(ltop, znstar_generate(n,gen));
     158             : }
     159             : 
     160             : /* Compute conductor of H */
     161             : long
     162         147 : znstar_conductor(long n, GEN H)
     163             : {
     164         147 :   pari_sp av = avma;
     165         147 :   long i, j, cnd = n;
     166         147 :   GEN F = factoru(n), P = gel(F,1), E = gel(F,2);
     167         357 :   for (i = lg(P)-1; i > 0; i--)
     168             :   {
     169         210 :     long p = P[i], e = E[i], q = n;
     170         210 :     if (DEBUGLEVEL>=4) err_printf("SubCyclo: testing %ld^%ld\n",p,e);
     171         308 :     for (  ; e >= 1; e--)
     172             :     {
     173         252 :       long z = 1;
     174         252 :       q /= p;
     175        3115 :       for (j = 1; j < p; j++)
     176             :       {
     177        3017 :         z += q;
     178        3017 :         if (!F2v_coeff(gel(H,3),z) && ugcd(z,n)==1) break;
     179             :       }
     180         252 :       if (j < p)
     181             :       {
     182         154 :         if (DEBUGLEVEL>=4) err_printf("SubCyclo: %ld not found\n",z);
     183         154 :         break;
     184             :       }
     185          98 :       cnd /= p;
     186          98 :       if (DEBUGLEVEL>=4) err_printf("SubCyclo: new conductor:%ld\n",cnd);
     187             :     }
     188             :   }
     189         147 :   if (DEBUGLEVEL>=6) err_printf("SubCyclo: conductor:%ld\n",cnd);
     190         147 :   avma = av; return cnd;
     191             : }
     192             : 
     193             : /* Compute the orbits of a subgroups of Z/nZ given by a generator
     194             :  * or a set of generators given as a vector.
     195             :  */
     196             : GEN
     197         119 : znstar_cosets(long n, long phi_n, GEN H)
     198             : {
     199             :   long    k;
     200         119 :   long    c = 0;
     201         119 :   long    card   = znstar_order(H);
     202         119 :   long    index  = phi_n/card;
     203         119 :   GEN     cosets = cgetg(index+1,t_VECSMALL);
     204         119 :   pari_sp ltop = avma;
     205         119 :   GEN     bits   = zero_F2v(n);
     206         483 :   for (k = 1; k <= index; k++)
     207             :   {
     208         364 :     for (c++ ; F2v_coeff(bits,c) || ugcd(c,n)!=1; c++);
     209         364 :     cosets[k]=c;
     210         364 :     znstar_coset_bits_inplace(n, H, bits, c);
     211             :   }
     212         119 :   avma=ltop;
     213         119 :   return cosets;
     214             : }
     215             : 
     216             : 
     217             : /*************************************************************************/
     218             : /**                                                                     **/
     219             : /**                     znstar/HNF interface                            **/
     220             : /**                                                                     **/
     221             : /*************************************************************************/
     222             : static GEN
     223         749 : vecmod_to_vecsmall(GEN z)
     224             : {
     225         749 :   long i, l = lg(z);
     226         749 :   GEN x = cgetg(l, t_VECSMALL);
     227        1729 :   for (i=1; i<l; i++) {
     228         980 :     GEN c = gel(z,i);
     229         980 :     if (typ(c) == t_INTMOD) c = gel(c,2);
     230         980 :     x[i] = itos(c);
     231             :   }
     232         749 :   return x;
     233             : }
     234             : /* Convert a true znstar output by znstar to a `small znstar' */
     235             : GEN
     236         749 : znstar_small(GEN zn)
     237             : {
     238         749 :   GEN Z = cgetg(4,t_VEC);
     239         749 :   gel(Z,1) = icopy(gmael3(zn,3,1,1));
     240         749 :   gel(Z,2) = vec_to_vecsmall(gel(zn,2));
     241         749 :   gel(Z,3) = vecmod_to_vecsmall(gel(zn,3)); return Z;
     242             : }
     243             : 
     244             : /* Compute generators for the subgroup of (Z/nZ)* given in HNF. */
     245             : GEN
     246        1659 : znstar_hnf_generators(GEN Z, GEN M)
     247             : {
     248        1659 :   long j, h, l = lg(M);
     249        1659 :   GEN gen = cgetg(l, t_VECSMALL);
     250        1659 :   pari_sp ltop = avma;
     251        1659 :   GEN zgen = gel(Z,3);
     252        1659 :   ulong n = itou(gel(Z,1));
     253        3948 :   for (j = 1; j < l; j++)
     254             :   {
     255        2289 :     GEN Mj = gel(M,j);
     256        2289 :     gen[j] = 1;
     257        5978 :     for (h = 1; h < l; h++)
     258             :     {
     259        3689 :       ulong u = itou(gel(Mj,h));
     260        3689 :       if (!u) continue;
     261        2457 :       gen[j] = Fl_mul(uel(gen,j), Fl_powu(uel(zgen,h), u, n), n);
     262             :     }
     263             :   }
     264        1659 :   avma = ltop; return gen;
     265             : }
     266             : 
     267             : GEN
     268        1596 : znstar_hnf(GEN Z, GEN M)
     269             : {
     270        1596 :   return znstar_generate(itos(gel(Z,1)),znstar_hnf_generators(Z,M));
     271             : }
     272             : 
     273             : GEN
     274        1596 : znstar_hnf_elts(GEN Z, GEN H)
     275             : {
     276        1596 :   pari_sp ltop = avma;
     277        1596 :   GEN G = znstar_hnf(Z,H);
     278        1596 :   long n = itos(gel(Z,1));
     279        1596 :   return gerepileupto(ltop, znstar_elts(n,G));
     280             : }
     281             : 
     282             : /*************************************************************************/
     283             : /**                                                                     **/
     284             : /**                     polsubcyclo                                     **/
     285             : /**                                                                     **/
     286             : /*************************************************************************/
     287             : 
     288             : static GEN
     289         147 : gscycloconductor(GEN g, long n, long flag)
     290             : {
     291         147 :   if (flag==2) retmkvec2(gcopy(g), stoi(n));
     292         147 :   return g;
     293             : }
     294             : 
     295             : static long
     296          91 : lift_check_modulus(GEN H, long n)
     297             : {
     298             :   long h;
     299          91 :   switch(typ(H))
     300             :   {
     301             :     case t_INTMOD:
     302           7 :       if (!equalsi(n, gel(H,1)))
     303           7 :         pari_err_MODULUS("galoissubcyclo", stoi(n), gel(H,1));
     304           0 :       H = gel(H,2);
     305             :     case t_INT:
     306          84 :       h = smodis(H,n);
     307          84 :       if (ugcd(h,n) != 1) pari_err_COPRIME("galoissubcyclo", H,stoi(n));
     308          84 :       return h;
     309             :   }
     310           0 :   pari_err_TYPE("galoissubcyclo [subgroup]", H);
     311           0 :   return 0;/*not reached*/
     312             : }
     313             : 
     314             : /* Compute z^ex using the baby-step/giant-step table powz
     315             :  * with only one multiply.
     316             :  * In the modular case, the result is not reduced. */
     317             : static GEN
     318      215348 : polsubcyclo_powz(GEN powz, long ex)
     319             : {
     320      215348 :   long m = lg(gel(powz,1))-1, q = ex/m, r = ex%m; /*ex=m*q+r*/
     321      215348 :   GEN g = gmael(powz,1,r+1), G = gmael(powz,2,q+1);
     322      215348 :   return (lg(powz)==4)? mulreal(g,G): gmul(g,G);
     323             : }
     324             : 
     325             : static GEN
     326        3059 : polsubcyclo_complex_bound(pari_sp av, GEN V, long prec)
     327             : {
     328        3059 :   GEN pol = real_i(roots_to_pol(V,0));
     329        3059 :   return gerepileuptoint(av, ceil_safe(gsupnorm(pol,prec)));
     330             : }
     331             : 
     332             : /* Newton sums mod le. if le==NULL, works with complex instead */
     333             : static GEN
     334        5880 : polsubcyclo_cyclic(long n, long d, long m ,long z, long g, GEN powz, GEN le)
     335             : {
     336        5880 :   GEN V = cgetg(d+1,t_VEC);
     337        5880 :   ulong base = 1;
     338             :   long i,k;
     339             :   pari_timer ti;
     340        5880 :   if (DEBUGLEVEL >= 6) timer_start(&ti);
     341       52136 :   for (i=1; i<=d; i++, base = Fl_mul(base,z,n))
     342             :   {
     343       46256 :     pari_sp av = avma;
     344       46256 :     long ex = base;
     345       46256 :     GEN s = gen_0;
     346      233968 :     for (k=0; k<m; k++, ex = Fl_mul(ex,g,n))
     347             :     {
     348      187712 :       s = gadd(s, polsubcyclo_powz(powz,ex));
     349      187712 :       if ((k&0xff)==0) s = gerepileupto(av,s);
     350             :     }
     351       46256 :     if (le) s = modii(s, le);
     352       46256 :     gel(V,i) = gerepileupto(av, s);
     353             :   }
     354        5880 :   if (DEBUGLEVEL >= 6) timer_printf(&ti, "polsubcyclo_cyclic");
     355        5880 :   return V;
     356             : }
     357             : 
     358             : struct _subcyclo_orbits_s
     359             : {
     360             :   GEN powz;
     361             :   GEN *s;
     362             :   ulong count;
     363             :   pari_sp ltop;
     364             : };
     365             : 
     366             : static void
     367       27636 : _subcyclo_orbits(struct _subcyclo_orbits_s *data, long k)
     368             : {
     369       27636 :   GEN powz = data->powz;
     370       27636 :   GEN *s = data->s;
     371             : 
     372       27636 :   if (!data->count) data->ltop = avma;
     373       27636 :   *s = gadd(*s, polsubcyclo_powz(powz,k));
     374       27636 :   data->count++;
     375       27636 :   if ((data->count & 0xffUL) == 0) *s = gerepileupto(data->ltop, *s);
     376       27636 : }
     377             : 
     378             : /* Newton sums mod le. if le==NULL, works with complex instead */
     379             : static GEN
     380         238 : polsubcyclo_orbits(long n, GEN H, GEN O, GEN powz, GEN le)
     381             : {
     382         238 :   long i, d = lg(O);
     383         238 :   GEN V = cgetg(d,t_VEC);
     384             :   struct _subcyclo_orbits_s data;
     385         238 :   long lle = le?lg(le)*2+1: 2*lg(gmael(powz,1,2))+3;/*dvmdii uses lx+ly space*/
     386         238 :   data.powz = powz;
     387         966 :   for(i=1; i<d; i++)
     388             :   {
     389         728 :     GEN s = gen_0;
     390         728 :     pari_sp av = avma;
     391         728 :     (void)new_chunk(lle);
     392         728 :     data.count = 0;
     393         728 :     data.s     = &s;
     394         728 :     znstar_coset_func(n, H, (void (*)(void *,long)) _subcyclo_orbits,
     395         728 :       (void *) &data, O[i]);
     396         728 :     avma = av; /* HACK */
     397         728 :     gel(V,i) = le? modii(s,le): gcopy(s);
     398             :   }
     399         238 :   return V;
     400             : }
     401             : 
     402             : static GEN
     403        8491 : polsubcyclo_start(long n, long d, long o, GEN borne, long *ptr_val,long *ptr_l)
     404             : {
     405             :   pari_sp av;
     406             :   GEN le, z, gl;
     407             :   long i, l, e, val;
     408        8491 :   l = n+1; e = 1;
     409        8491 :   while(!uisprime(l)) { l += n; e++; }
     410        8491 :   if (DEBUGLEVEL >= 4) err_printf("Subcyclo: prime l=%ld\n",l);
     411        8491 :   gl = utoipos(l); av = avma;
     412        8491 :   if (!borne)
     413             :   { /* Use vecmax(Vec((x+o)^d)) = max{binomial(d,i)*o^i ;1<=i<=d} */
     414         154 :     i = d-(1+d)/(1+o);
     415         154 :     borne = mulii(binomial(utoipos(d),i),powuu(o,i));
     416             :   }
     417        8491 :   if (DEBUGLEVEL >= 4) err_printf("Subcyclo: bound=2^%ld\n",expi(borne));
     418        8491 :   val = logint(shifti(borne,2), gl) + 1;
     419        8491 :   avma = av;
     420        8491 :   if (DEBUGLEVEL >= 4) err_printf("Subcyclo: val=%ld\n",val);
     421        8491 :   le = powiu(gl,val);
     422        8491 :   z = utoipos( Fl_powu(pgener_Fl(l), e, l) );
     423        8491 :   z = Zp_sqrtnlift(gen_1,utoipos(n),z,gl,val);
     424        8491 :   *ptr_val = val;
     425        8491 :   *ptr_l = l;
     426        8491 :   return gmodulo(z,le);
     427             : }
     428             : 
     429             : /*Fill in the powz table:
     430             :  *  powz[1]: baby-step
     431             :  *  powz[2]: giant-step
     432             :  *  powz[3] exists only if the field is real (value is ignored). */
     433             : static GEN
     434        3059 : polsubcyclo_complex_roots(long n, long real, long prec)
     435             : {
     436        3059 :   long i, m = (long)(1+sqrt((double) n));
     437        3059 :   GEN bab, gig, powz = cgetg(real?4:3, t_VEC);
     438             : 
     439        3059 :   bab = cgetg(m+1,t_VEC);
     440        3059 :   gel(bab,1) = gen_1;
     441        3059 :   gel(bab,2) = rootsof1u_cx(n, prec); /* = e_n(1) */
     442        3059 :   for (i=3; i<=m; i++) gel(bab,i) = gmul(gel(bab,2),gel(bab,i-1));
     443        3059 :   gig = cgetg(m+1,t_VEC);
     444        3059 :   gel(gig,1) = gen_1;
     445        3059 :   gel(gig,2) = gmul(gel(bab,2),gel(bab,m));;
     446        3059 :   for (i=3; i<=m; i++) gel(gig,i) = gmul(gel(gig,2),gel(gig,i-1));
     447        3059 :   gel(powz,1) = bab;
     448        3059 :   gel(powz,2) = gig;
     449        3059 :   if (real) gel(powz,3) = gen_0;
     450        3059 :   return powz;
     451             : }
     452             : 
     453             : static GEN
     454       22477 : muliimod_sz(GEN x, GEN y, GEN l, long siz)
     455             : {
     456       22477 :   pari_sp av = avma;
     457             :   GEN p1;
     458       22477 :   (void)new_chunk(siz); /* HACK */
     459       22477 :   p1 = mulii(x,y);
     460       22477 :   avma = av; return modii(p1,l);
     461             : }
     462             : 
     463             : static GEN
     464        3059 : polsubcyclo_roots(long n, GEN zl)
     465             : {
     466        3059 :   GEN le = gel(zl,1), z = gel(zl,2);
     467        3059 :   long i, lle = lg(le)*3; /*Assume dvmdii use lx+ly space*/
     468        3059 :   long m = (long)(1+sqrt((double) n));
     469        3059 :   GEN bab, gig, powz = cgetg(3,t_VEC);
     470             :   pari_timer ti;
     471        3059 :   if (DEBUGLEVEL >= 6) timer_start(&ti);
     472        3059 :   bab = cgetg(m+1,t_VEC);
     473        3059 :   gel(bab,1) = gen_1;
     474        3059 :   gel(bab,2) = icopy(z);
     475        3059 :   for (i=3; i<=m; i++) gel(bab,i) = muliimod_sz(z,gel(bab,i-1),le,lle);
     476        3059 :   gig = cgetg(m+1,t_VEC);
     477        3059 :   gel(gig,1) = gen_1;
     478        3059 :   gel(gig,2) = muliimod_sz(z,gel(bab,m),le,lle);;
     479        3059 :   for (i=3; i<=m; i++) gel(gig,i) = muliimod_sz(gel(gig,2),gel(gig,i-1),le,lle);
     480        3059 :   if (DEBUGLEVEL >= 6) timer_printf(&ti, "polsubcyclo_roots");
     481        3059 :   gel(powz,1) = bab;
     482        3059 :   gel(powz,2) = gig; return powz;
     483             : }
     484             : 
     485             : GEN
     486         154 : galoiscyclo(long n, long v)
     487             : {
     488         154 :   ulong av = avma;
     489             :   GEN grp, G, z, le, L, elts;
     490             :   long val, l, i, j, k;
     491         154 :   GEN zn = znstar(stoi(n));
     492         154 :   long card = itos(gel(zn,1));
     493         154 :   GEN gen = vec_to_vecsmall(lift_shallow(gel(zn,3)));
     494         154 :   GEN ord = gtovecsmall(gel(zn,2));
     495             : 
     496         154 :   z = polsubcyclo_start(n,card/2,2,NULL,&val,&l);
     497         154 :   le = gel(z,1); z = gel(z,2);
     498         154 :   L = cgetg(1+card,t_VEC);
     499         154 :   gel(L,1) = z;
     500         336 :   for (j = 1, i = 1; j < lg(gen); j++)
     501             :   {
     502         182 :     long c = i * (ord[j]-1);
     503         182 :     for (k = 1; k <= c; k++) gel(L,++i) = Fp_powu(gel(L,k), gen[j], le);
     504             :   }
     505         154 :   G = abelian_group(ord);
     506         154 :   elts = group_elts(G, card); /*not stack clean*/
     507         154 :   grp = cgetg(9, t_VEC);
     508         154 :   gel(grp,1) = polcyclo(n,v);
     509         154 :   gel(grp,2) = mkvec3(stoi(l), stoi(val), icopy(le));
     510         154 :   gel(grp,3) = gcopy(L);
     511         154 :   gel(grp,4) = FpV_invVandermonde(L,  NULL, le);
     512         154 :   gel(grp,5) = gen_1;
     513         154 :   gel(grp,6) = gcopy(elts);
     514         154 :   gel(grp,7) = gcopy(gel(G,1));
     515         154 :   gel(grp,8) = gcopy(gel(G,2));
     516         154 :   return gerepileupto(av, grp);
     517             : }
     518             : 
     519             : /* Convert a bnrinit(Q,n) to a znstar(n)
     520             :  * complex is set to 0 if the bnr is real and to 1 if it is complex.
     521             :  * Not stack clean */
     522             : GEN
     523          14 : bnr_to_znstar(GEN bnr, long *complex)
     524             : {
     525             :   GEN gen, cond, v, bid;
     526             :   long l2, i;
     527          14 :   checkbnr(bnr);
     528          14 :   bid = bnr_get_bid(bnr);
     529          14 :   gen = bnr_get_gen(bnr);
     530          14 :   if (nf_get_degree(bnr_get_nf(bnr)) != 1)
     531           7 :     pari_err_DOMAIN("bnr_to_znstar", "bnr", "!=", strtoGENstr("Q"), bnr);
     532             :   /* cond is the finite part of the conductor,
     533             :    * complex is the infinite part*/
     534           7 :   cond = gcoeff(bid_get_ideal(bid), 1, 1);
     535           7 :   *complex = signe(gel(bid_get_arch(bid), 1));
     536           7 :   l2 = lg(gen);
     537           7 :   v = cgetg(l2, t_VEC);
     538          35 :   for (i = 1; i < l2; ++i)
     539             :   {
     540          28 :     GEN x = gel(gen,i);
     541          28 :     switch(typ(x))
     542             :     {
     543           0 :       case t_MAT: x = gcoeff(x,1,1); break;
     544           0 :       case t_COL: x = gel(x,1); break;
     545             :     }
     546          28 :     gel(v,i) = gmodulo(absi(x), cond);
     547             :   }
     548           7 :   return mkvec3(bnr_get_no(bnr), bnr_get_cyc(bnr), v);
     549             : }
     550             : 
     551             : GEN
     552         196 : galoissubcyclo(GEN N, GEN sg, long flag, long v)
     553             : {
     554         196 :   pari_sp ltop= avma, av;
     555         196 :   GEN H, V, B, zl, L, T, le, powz, O, Z = NULL;
     556         196 :   long i, card, phi_n, val,l, n, cnd, complex=1;
     557             :   pari_timer ti;
     558             : 
     559         196 :   if (flag<0 || flag>2) pari_err_FLAG("galoissubcyclo");
     560         196 :   if (v < 0) v = 0;
     561         196 :   if (!sg) sg = gen_1;
     562         196 :   switch(typ(N))
     563             :   {
     564             :     case t_INT:
     565         112 :       n = itos(N);
     566         112 :       if (n < 1)
     567           7 :         pari_err_DOMAIN("galoissubcyclo", "degree", "<=", gen_0, stoi(n));
     568         105 :       break;
     569             :     case t_VEC:
     570          84 :       if (lg(N)==7) N = bnr_to_znstar(N,&complex);
     571          77 :       if (lg(N)==4)
     572             :       { /* znstar */
     573          77 :         GEN gen = gel(N,3);
     574          77 :         Z = N;
     575          77 :         if (typ(gen)!=t_VEC) pari_err_TYPE("galoissubcyclo",gen);
     576          77 :         if (lg(gen) == 1) n = 1;
     577          77 :         else if (typ(gel(gen,1)) == t_INTMOD)
     578             :         {
     579          70 :           GEN z = gel(gen,1);
     580          70 :           n = itos(gel(z,1));
     581             :         } else
     582             :         {
     583           7 :           pari_err_TYPE("galoissubcyclo",N);
     584           0 :           return NULL;/*Not reached*/
     585             :         }
     586          70 :         break;
     587             :       }
     588             :     default: /*fall through*/
     589           0 :       pari_err_TYPE("galoissubcyclo",N);
     590           0 :       return NULL;/*Not reached*/
     591             :   }
     592         175 :   if (n==1) { avma = ltop; return deg1pol_shallow(gen_1,gen_m1,v); }
     593             : 
     594         175 :   switch(typ(sg))
     595             :   {
     596             :      case t_INTMOD: case t_INT:
     597          91 :       V = mkvecsmall( lift_check_modulus(sg,n) );
     598          84 :       break;
     599             :     case t_VECSMALL:
     600           0 :       V = gcopy(sg);
     601           0 :       for (i=1; i<lg(V); i++) { V[i] %= n; if (V[i] < 0) V[i] += n; }
     602           0 :       break;
     603             :     case t_VEC:
     604             :     case t_COL:
     605           0 :       V = cgetg(lg(sg),t_VECSMALL);
     606           0 :       for(i=1;i<lg(sg);i++) V[i] = lift_check_modulus(gel(sg,i),n);
     607           0 :       break;
     608             :     case t_MAT:
     609          84 :       if (lg(sg) == 1 || lg(sg) != lgcols(sg))
     610           7 :         pari_err_TYPE("galoissubcyclo [H not in HNF]", sg);
     611          77 :       if (!Z) pari_err_TYPE("galoissubcyclo [N not a bnrinit or znstar]", sg);
     612          70 :       if ( lg(gel(Z,2)) != lg(sg) ) pari_err_DIM("galoissubcyclo");
     613          63 :       V = znstar_hnf_generators(znstar_small(Z),sg);
     614          63 :       break;
     615             :     default:
     616           0 :       pari_err_TYPE("galoissubcyclo",sg);
     617           0 :       return NULL;/*Not reached*/
     618             :   }
     619         147 :   if (!complex) V = vecsmall_append(V,n-1); /*add complex conjugation*/
     620         147 :   H = znstar_generate(n,V);
     621         147 :   if (DEBUGLEVEL >= 6)
     622             :   {
     623           0 :     err_printf("Subcyclo: elements:");
     624           0 :     for (i=1;i<n;i++)
     625           0 :       if (F2v_coeff(gel(H,3),i)) err_printf(" %ld",i);
     626           0 :     err_printf("\n");
     627             :   }
     628             :   /* field is real iff z -> conj(z) = z^-1 = z^(n-1) is in H */
     629         147 :   complex = !F2v_coeff(gel(H,3),n-1);
     630         147 :   if (DEBUGLEVEL >= 6) err_printf("Subcyclo: complex=%ld\n",complex);
     631         147 :   if (DEBUGLEVEL >= 1) timer_start(&ti);
     632         147 :   cnd = znstar_conductor(n,H);
     633         147 :   if (DEBUGLEVEL >= 1) timer_printf(&ti, "znstar_conductor");
     634         147 :   if (flag == 1)  { avma=ltop; return stoi(cnd); }
     635         147 :   if (cnd == 1)
     636             :   {
     637          28 :     avma = ltop;
     638          28 :     return gscycloconductor(deg1pol_shallow(gen_1,gen_m1,v),1,flag);
     639             :   }
     640         119 :   if (n != cnd)
     641             :   {
     642          35 :     H = znstar_reduce_modulus(H, cnd);
     643          35 :     n = cnd;
     644             :   }
     645         119 :   card = znstar_order(H);
     646         119 :   phi_n = eulerphiu(n);
     647         119 :   if (card == phi_n)
     648             :   {
     649           0 :     avma = ltop;
     650           0 :     return gscycloconductor(polcyclo(n,v),n,flag);
     651             :   }
     652         119 :   O = znstar_cosets(n, phi_n, H);
     653         119 :   if (DEBUGLEVEL >= 1) timer_printf(&ti, "znstar_cosets");
     654         119 :   if (DEBUGLEVEL >= 6) err_printf("Subcyclo: orbits=%Ps\n",O);
     655         119 :   if (DEBUGLEVEL >= 4)
     656           0 :     err_printf("Subcyclo: %ld orbits with %ld elements each\n",phi_n/card,card);
     657         119 :   av = avma;
     658         119 :   powz = polsubcyclo_complex_roots(n,!complex,LOWDEFAULTPREC);
     659         119 :   L = polsubcyclo_orbits(n,H,O,powz,NULL);
     660         119 :   B = polsubcyclo_complex_bound(av,L,LOWDEFAULTPREC);
     661         119 :   zl = polsubcyclo_start(n,phi_n/card,card,B,&val,&l);
     662         119 :   powz = polsubcyclo_roots(n,zl);
     663         119 :   le = gel(zl,1);
     664         119 :   L = polsubcyclo_orbits(n,H,O,powz,le);
     665         119 :   if (DEBUGLEVEL >= 6) timer_start(&ti);
     666         119 :   T = FpV_roots_to_pol(L,le,v);
     667         119 :   if (DEBUGLEVEL >= 6) timer_printf(&ti, "roots_to_pol");
     668         119 :   T = FpX_center(T,le,shifti(le,-1));
     669         119 :   return gerepileupto(ltop, gscycloconductor(T,n,flag));
     670             : }
     671             : 
     672             : /* Z = znstar(n) cyclic. n = 1,2,4,p^a or 2p^a,
     673             :  * and d | phi(n) = 1,1,2,(p-1)p^(a-1) */
     674             : static GEN
     675        3003 : polsubcyclo_g(long n, long d, GEN Z, long v)
     676             : {
     677        3003 :   pari_sp ltop = avma;
     678             :   long o, p, r, g, gd, l , val;
     679             :   GEN zl, L, T, le, B, powz;
     680             :   pari_timer ti;
     681        3003 :   if (d==1) return deg1pol_shallow(gen_1,gen_m1,v); /* get rid of n=1,2 */
     682        3003 :   if ((n & 3) == 2) n >>= 1;
     683             :   /* n = 4 or p^a, p odd */
     684        3003 :   o = itos(gel(Z,1));
     685        3003 :   g = itos(gmael3(Z,3,1,2));
     686        3003 :   p = n / ugcd(n,o); /* p^a / gcd(p^a,phi(p^a)) = p*/
     687        3003 :   r = ugcd(d,n); /* = p^(v_p(d)) < n */
     688        3003 :   n = r*p; /* n is now the conductor */
     689        3003 :   o = n-r; /* = phi(n) */
     690        3003 :   if (o == d) return polcyclo(n,v);
     691        2940 :   o /= d;
     692        2940 :   gd = Fl_powu(g%n, d, n);
     693             :   /*FIXME: If degree is small, the computation of B is a waste of time*/
     694        2940 :   powz = polsubcyclo_complex_roots(n,(o&1)==0,LOWDEFAULTPREC);
     695        2940 :   L = polsubcyclo_cyclic(n,d,o,g,gd,powz,NULL);
     696        2940 :   B = polsubcyclo_complex_bound(ltop,L,LOWDEFAULTPREC);
     697        2940 :   zl = polsubcyclo_start(n,d,o,B,&val,&l);
     698        2940 :   le = gel(zl,1);
     699        2940 :   powz = polsubcyclo_roots(n,zl);
     700        2940 :   L = polsubcyclo_cyclic(n,d,o,g,gd,powz,le);
     701        2940 :   if (DEBUGLEVEL >= 6) timer_start(&ti);
     702        2940 :   T = FpV_roots_to_pol(L,le,v);
     703        2940 :   if (DEBUGLEVEL >= 6) timer_printf(&ti, "roots_to_pol");
     704        2940 :   return gerepileupto(ltop, FpX_center(T,le,shifti(le,-1)));
     705             : }
     706             : 
     707             : GEN
     708        3031 : polsubcyclo(long n, long d, long v)
     709             : {
     710        3031 :   pari_sp ltop = avma;
     711             :   GEN L, Z;
     712        3031 :   if (v<0) v = 0;
     713        3031 :   if (d<=0) pari_err_DOMAIN("polsubcyclo","d","<=",gen_0,stoi(d));
     714        3024 :   if (n<=0) pari_err_DOMAIN("polsubcyclo","n","<=",gen_0,stoi(n));
     715        3017 :   Z = znstar(stoi(n));
     716        3017 :   if (!dvdis(gel(Z,1), d)) { avma = ltop; return cgetg(1, t_VEC); }
     717        3017 :   if (lg(gel(Z,2)) == 2)
     718             :   { /* faster but Z must be cyclic */
     719        3003 :     avma = ltop;
     720        3003 :     return polsubcyclo_g(n, d, Z, v);
     721             :   }
     722          14 :   L = subgrouplist(gel(Z,2), mkvec(stoi(d)));
     723          14 :   if (lg(L) == 2)
     724           7 :     return gerepileupto(ltop, galoissubcyclo(Z, gel(L,1), 0, v));
     725             :   else
     726             :   {
     727           7 :     GEN V = cgetg(lg(L),t_VEC);
     728             :     long i;
     729           7 :     for (i=1; i< lg(V); i++) gel(V,i) = galoissubcyclo(Z, gel(L,i), 0, v);
     730           7 :     return gerepileupto(ltop, V);
     731             :   }
     732             : }
     733             : 
     734             : struct aurifeuille_t {
     735             :   GEN z, le;
     736             :   ulong l;
     737             :   long e;
     738             : };
     739             : 
     740             : /* Let z a primitive n-th root of 1, n > 1, A an integer such that
     741             :  * Aurifeuillian factorization of Phi_n(A) exists ( z.A is a square in Q(z) ).
     742             :  * Let G(p) the Gauss sum mod p prime:
     743             :  *      sum_x (x|p) z^(xn/p) for p odd,  i - 1 for p = 2 [ i := z^(n/4) ]
     744             :  * We have N(-1) = Nz = 1 (n != 1,2), and
     745             :  *      G^2 = (-1|p) p for p odd,  G^2 = -2i for p = 2
     746             :  * In particular, for odd A, (-1|A) A = g^2 is a square. If A = prod p^{e_p},
     747             :  * sigma_j(g) = \prod_p (sigma_j G(p)))^e_p = \prod_p (j|p)^e_p g = (j|A) g
     748             :  * n odd  : z^2 is a primitive root, A = g^2
     749             :  *   Phi_n(A) = N(A - z^2) = N(g - z) N(g + z)
     750             :  *
     751             :  * n = 2 (4) : -z^2 is a primitive root, -A = g^2
     752             :  *   Phi_n(A) = N(A - (-z^2)) = N(g^2 - z^2)  [ N(-1) = 1 ]
     753             :  *                            = N(g - z) N(g + z)
     754             :  *
     755             :  * n = 4 (8) : i z^2 primitive root, -Ai = g^2
     756             :  *   Phi_n(A) = N(A - i z^2) = N(-Ai -  z^2) = N(g - z) N(g + z)
     757             :  * sigma_j(g) / g =  (j|A)  if j = 1 (4)
     758             :  *                  (-j|A)i if j = 3 (4)
     759             :  *   */
     760             : /* factor Phi_n(A), Astar: A* = squarefree kernel of A, P = odd prime divisors
     761             :  * of n */
     762             : static GEN
     763        5278 : factor_Aurifeuille_aux(GEN A, long Astar, long n, GEN P,
     764             :                        struct aurifeuille_t *S)
     765             : {
     766             :   pari_sp av;
     767        5278 :   GEN f, a, b, s, powers, z = S->z, le = S->le;
     768        5278 :   long j, k, maxjump, lastj, e = S->e;
     769        5278 :   ulong l = S->l;
     770             :   char *invertible;
     771             : 
     772        5278 :   if ((n & 7) == 4)
     773             :   { /* A^* even */
     774        5215 :     GEN i = Fp_powu(z, n>>2, le), z2 = Fp_sqr(z, le);
     775             : 
     776        5215 :     invertible = stack_malloc(n); /* even indices unused */
     777        5215 :     for (j = 1; j < n; j+=2) invertible[j] = 1;
     778        5271 :     for (k = 1; k < lg(P); k++)
     779             :     {
     780          56 :       long p = P[k];
     781          56 :       for (j = p; j < n; j += 2*p) invertible[j] = 0;
     782             :     }
     783        5215 :     lastj = 1; maxjump = 2;
     784       10990 :     for (j= 3; j < n; j+=2)
     785        5775 :       if (invertible[j]) {
     786        5663 :         long jump = j - lastj;
     787        5663 :         if (jump > maxjump) maxjump = jump;
     788        5663 :         lastj = j;
     789             :       }
     790        5215 :     powers = cgetg(maxjump+1, t_VEC); /* powers[k] = z^k, odd indices unused */
     791        5215 :     gel(powers,2) = z2;
     792        5271 :     for (k = 4; k <= maxjump; k+=2)
     793         112 :       gel(powers,k) = odd(k>>1)? Fp_mul(gel(powers, k-2), z2, le)
     794          56 :                                : Fp_sqr(gel(powers, k>>1), le);
     795             : 
     796        5215 :     if (Astar == 2)
     797             :     { /* important special case (includes A=2), split for efficiency */
     798        5194 :       if (!equalis(A, 2))
     799             :       {
     800          14 :         GEN f = sqrti(shifti(A,-1)), mf = Fp_neg(f,le), fi = Fp_mul(f,i,le);
     801          14 :         a = Fp_add(mf, fi, le);
     802          14 :         b = Fp_sub(mf, fi, le);
     803             :       }
     804             :       else
     805             :       {
     806        5180 :         a = addsi(-1,i);
     807        5180 :         b = subsi(-1,i);
     808             :       }
     809        5194 :       av = avma;
     810        5194 :       s = z; f = subii(a, s); lastj = 1;
     811       10640 :       for (j = 3, k = 0; j < n; j+=2)
     812        5446 :         if (invertible[j])
     813             :         {
     814        5376 :           s = Fp_mul(gel(powers, j-lastj), s, le); /* z^j */
     815        5376 :           lastj = j;
     816        5376 :           f = Fp_mul(f, subii((j & 3) == 1? a: b, s), le);
     817        5376 :           if (++k == 0x1ff) { gerepileall(av, 2, &s, &f); k = 0; }
     818             :         }
     819             :     }
     820             :     else
     821             :     {
     822          21 :       GEN ma, mb, B = Fp_mul(A, i, le), gl = utoipos(l);
     823             :       long t;
     824          21 :       Astar >>= 1;
     825          21 :       t = Astar & 3; if (Astar < 0) t = 4-t; /* t = 1 or 3 */
     826          21 :       if (t == 1) B = Fp_neg(B, le);
     827          21 :       a = Zp_sqrtlift(B, Fp_sqrt(B, gl), gl, e);
     828          21 :       b = Fp_mul(a, i, le);
     829          21 :       ma = Fp_neg(a, le);
     830          21 :       mb = Fp_neg(b, le);
     831          21 :       av = avma;
     832          21 :       s = z; f = subii(a, s); lastj = 1;
     833         350 :       for (j = 3, k = 0; j<n; j+=2)
     834         329 :         if (invertible[j])
     835             :         {
     836             :           GEN t;
     837         287 :           if ((j & 3) == 1) t = (kross(j, Astar) < 0)? ma: a;
     838         154 :           else              t = (kross(j, Astar) < 0)? mb: b;
     839         287 :           s = Fp_mul(gel(powers, j-lastj), s, le); /* z^j */
     840         287 :           lastj = j;
     841         287 :           f = Fp_mul(f, subii(t, s), le);
     842         287 :           if (++k == 0x1ff) { gerepileall(av, 2, &s, &f); k = 0; }
     843             :         }
     844             :     }
     845             :   }
     846             :   else /* A^* odd */
     847             :   {
     848             :     ulong g;
     849          63 :     if ((n & 3) == 2)
     850             :     { /* A^* = 3 (mod 4) */
     851           0 :       A = negi(A); Astar = -Astar;
     852           0 :       z = Fp_neg(z, le);
     853           0 :       n >>= 1;
     854             :     }
     855             :     /* A^* = 1 (mod 4) */
     856          63 :     g = Fl_sqrt(umodiu(A,l), l);
     857          63 :     a = Zp_sqrtlift(A, utoipos(g), utoipos(l), e);
     858          63 :     b = negi(a);
     859             : 
     860          63 :     invertible = stack_malloc(n);
     861          63 :     for (j = 1; j < n; j++) invertible[j] = 1;
     862         168 :     for (k = 1; k < lg(P); k++)
     863             :     {
     864         105 :       long p = P[k];
     865         105 :       for (j = p; j < n; j += p) invertible[j] = 0;
     866             :     }
     867          63 :     lastj = 2; maxjump = 1;
     868        1141 :     for (j= 3; j < n; j++)
     869        1078 :       if (invertible[j]) {
     870         714 :         long jump = j - lastj;
     871         714 :         if (jump > maxjump) maxjump = jump;
     872         714 :         lastj = j;
     873             :       }
     874          63 :     powers = cgetg(maxjump+1, t_VEC); /* powers[k] = z^k */
     875          63 :     gel(powers,1) = z;
     876         147 :     for (k = 2; k <= maxjump; k++)
     877         210 :       gel(powers,k) = odd(k)? Fp_mul(gel(powers, k-1), z, le)
     878         126 :                             : Fp_sqr(gel(powers, k>>1), le);
     879          63 :     av = avma;
     880          63 :     s = z; f = subii(a, s); lastj = 1;
     881        1204 :     for(j = 2, k = 0; j < n; j++)
     882        1141 :       if (invertible[j])
     883             :       {
     884         777 :         s = Fp_mul(gel(powers, j-lastj), s, le);
     885         777 :         lastj = j;
     886         777 :         f = Fp_mul(f, subii(kross(j,Astar)==1? a: b, s), le);
     887         777 :         if (++k == 0x1ff) { gerepileall(av, 2, &s, &f); k = 0; }
     888             :       }
     889             :   }
     890        5278 :   return f;
     891             : }
     892             : 
     893             : /* fd = factoru(odd part of d = d or d/4). Return eulerphi(d) */
     894             : static ulong
     895        5278 : phi(long d, GEN fd)
     896             : {
     897        5278 :   GEN P = gel(fd,1), E = gel(fd,2);
     898        5278 :   long i, l = lg(P);
     899        5278 :   ulong phi = 1;
     900        5439 :   for (i = 1; i < l; i++)
     901             :   {
     902         161 :     ulong p = P[i], e = E[i];
     903         161 :     phi *= upowuu(p, e-1)*(p-1);
     904             :   }
     905        5278 :   if (!odd(d)) phi <<= 1;
     906        5278 :   return phi;
     907             : }
     908             : 
     909             : static void
     910        5278 : Aurifeuille_init(GEN a, long d, GEN fd, struct aurifeuille_t *S)
     911             : {
     912        5278 :   GEN sqrta = sqrtr_abs(itor(a, LOWDEFAULTPREC));
     913        5278 :   GEN bound = ceil_safe(powru(addrs(sqrta,1), phi(d, fd)));
     914        5278 :   GEN zl = polsubcyclo_start(d, 0, 0, bound, &(S->e), (long*)&(S->l));
     915        5278 :   S->le = gel(zl,1);
     916        5278 :   S->z  = gel(zl,2);
     917        5278 : }
     918             : 
     919             : GEN
     920        5215 : factor_Aurifeuille_prime(GEN p, long d)
     921             : {
     922        5215 :   pari_sp av = avma;
     923             :   struct aurifeuille_t S;
     924             :   GEN fd;
     925             :   long pp;
     926        5215 :   if ((d & 3) == 2) { d >>= 1; p = negi(p); }
     927        5215 :   fd = factoru(odd(d)? d: d>>2);
     928        5215 :   pp = itos(p);
     929        5215 :   Aurifeuille_init(p, d, fd, &S);
     930        5215 :   return gerepileuptoint(av, factor_Aurifeuille_aux(p, pp, d, gel(fd,1), &S));
     931             : }
     932             : 
     933             : /* an algebraic factor of Phi_d(a), a != 0 */
     934             : GEN
     935          63 : factor_Aurifeuille(GEN a, long d)
     936             : {
     937          63 :   pari_sp av = avma;
     938             :   GEN fd, P, A;
     939          63 :   long i, lP, va = vali(a), sa, astar, D;
     940             :   struct aurifeuille_t S;
     941             : 
     942          63 :   if (d <= 0)
     943           0 :     pari_err_DOMAIN("factor_Aurifeuille", "degre", "<=",gen_0,stoi(d));
     944          63 :   if ((d & 3) == 2) { d >>= 1; a = negi(a); }
     945          63 :   if ((va & 1) == (d & 1)) { avma = av; return gen_1; }
     946          63 :   sa = signe(a);
     947          63 :   if (odd(d))
     948             :   {
     949             :     long a4;
     950          28 :     if (d == 1)
     951             :     {
     952           0 :       if (!Z_issquareall(a, &A)) return gen_1;
     953           0 :       return gerepileuptoint(av, addis(A,1));
     954             :     }
     955          28 :     A = va? shifti(a, -va): a;
     956          28 :     a4 = mod4(A); if (sa < 0) a4 = 4 - a4;
     957          28 :     if (a4 != 1) { avma = av; return gen_1; }
     958             :   }
     959          35 :   else if ((d & 7) == 4)
     960          35 :     A = shifti(a, -va);
     961             :   else
     962             :   {
     963           0 :     avma = av; return gen_1;
     964             :   }
     965             :   /* v_2(d) = 0 or 2. Kill 2 from factorization (minor efficiency gain) */
     966          63 :   fd = factoru(odd(d)? d: d>>2); P = gel(fd,1); lP = lg(P);
     967          63 :   astar = sa;
     968          63 :   if (odd(va)) astar <<= 1;
     969         147 :   for (i = 1; i < lP; i++)
     970          84 :     if (odd( (Z_lvalrem(A, P[i], &A)) ) ) astar *= P[i];
     971          63 :   if (sa < 0)
     972             :   { /* negate in place if possible */
     973          14 :     if (A == a) A = icopy(A);
     974          14 :     setabssign(A);
     975             :   }
     976          63 :   if (!Z_issquare(A)) { avma = av; return gen_1; }
     977             : 
     978          63 :   D = odd(d)? 1: 4;
     979          63 :   for (i = 1; i < lP; i++) D *= P[i];
     980          63 :   if (D != d) { a = powiu(a, d/D); d = D; }
     981             : 
     982          63 :   Aurifeuille_init(a, d, fd, &S);
     983          63 :   return gerepileuptoint(av, factor_Aurifeuille_aux(a, astar, d, P, &S));
     984             : }

Generated by: LCOV version 1.11