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 - hnf_snf.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.13.0 lcov report (development 25825-edc109b529) Lines: 1512 1674 90.3 %
Date: 2020-09-21 06:08:33 Functions: 86 95 90.5 %
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             : /**                HERMITE NORMAL FORM REDUCTION             **/
      19             : /**                                                          **/
      20             : /**************************************************************/
      21             : static GEN ZV_hnfgcdext(GEN A);
      22             : static GEN
      23          14 : hnfallgen(GEN x)
      24             : {
      25          14 :   GEN z = cgetg(3, t_VEC);
      26          14 :   gel(z,1) = RgM_hnfall(x, (GEN*)(z+2), 1);
      27          14 :   return z;
      28             : }
      29             : GEN
      30         252 : mathnf0(GEN x, long flag)
      31             : {
      32         252 :   switch(typ(x))
      33             :   {
      34          70 :     case t_VEC:
      35          70 :       if (RgV_is_ZV(x))
      36             :         switch (flag)
      37             :         {
      38          14 :           case 0:
      39          14 :             if (lg(x) == 1) return cgetg(1, t_MAT);
      40           7 :             retmkmat(mkcol(ZV_content(x)));
      41          21 :           case 1:
      42             :           case 4:
      43          21 :             return ZV_hnfgcdext(x);
      44             :         }
      45          35 :       x = gtomat(x); break;
      46         182 :     case t_MAT: break;
      47           0 :     default: pari_err_TYPE("mathnf0",x);
      48             :   }
      49             : 
      50         217 :   switch(flag)
      51             :   {
      52         168 :     case 0: case 2: return RgM_is_ZM(x)? ZM_hnf(x): RgM_hnfall(x,NULL,1);
      53          28 :     case 1: case 3: return RgM_is_ZM(x)? hnfall(x): hnfallgen(x);
      54           7 :     case 4: RgM_check_ZM(x, "mathnf0"); return hnflll(x);
      55          14 :     case 5: RgM_check_ZM(x, "mathnf0"); return hnfperm(x);
      56           0 :     default: pari_err_FLAG("mathnf");
      57             :   }
      58             :   return NULL; /* LCOV_EXCL_LINE */
      59             : }
      60             : 
      61             : /*******************************************************************/
      62             : /*                                                                 */
      63             : /*                SPECIAL HNF (FOR INTERNAL USE !!!)               */
      64             : /*                                                                 */
      65             : /*******************************************************************/
      66             : static int
      67     4359783 : count(GEN mat, long row, long len, long *firstnonzero)
      68             : {
      69     4359783 :   long j, n = 0;
      70             : 
      71   374976850 :   for (j=1; j<=len; j++)
      72             :   {
      73   371863080 :     long p = mael(mat,j,row);
      74   371863080 :     if (p)
      75             :     {
      76    11832772 :       if (labs(p)!=1) return -1;
      77    10586759 :       n++; *firstnonzero=j;
      78             :     }
      79             :   }
      80     3113770 :   return n;
      81             : }
      82             : 
      83             : static int
      84      164002 : count2(GEN mat, long row, long len)
      85             : {
      86             :   long j;
      87     4627949 :   for (j=len; j; j--)
      88     4544328 :     if (labs(mael(mat,j,row)) == 1) return j;
      89       83621 :   return 0;
      90             : }
      91             : 
      92             : static GEN
      93       55211 : hnffinal(GEN matgen,GEN perm,GEN* ptdep,GEN* ptB,GEN* ptC)
      94             : {
      95             :   GEN p1,p2,U,H,Hnew,Bnew,Cnew,diagH1;
      96       55211 :   GEN B = *ptB, C = *ptC, dep = *ptdep, depnew;
      97             :   pari_sp av;
      98             :   long i,j,k,s,i1,j1,zc;
      99       55211 :   long co = lg(C);
     100       55211 :   long col = lg(matgen)-1;
     101             :   long lnz, nlze, lig;
     102             : 
     103       55211 :   if (col == 0) return matgen;
     104       54875 :   lnz = nbrows(matgen); /* was called lnz-1 - nr in hnfspec */
     105       54875 :   nlze = nbrows(dep);
     106       54875 :   lig = nlze + lnz;
     107             :   /* H: lnz x lnz [disregarding initial 0 cols], U: col x col */
     108       54875 :   H = ZM_hnflll(matgen, &U, 0);
     109       54875 :   H += (lg(H)-1 - lnz); H[0] = evaltyp(t_MAT) | evallg(lnz+1);
     110             :   /* Only keep the part above the H (above the 0s is 0 since the dep rows
     111             :    * are dependent from the ones in matgen) */
     112       54875 :   zc = col - lnz; /* # of 0 columns, correspond to units */
     113       54875 :   if (nlze) { dep = ZM_mul(dep,U); dep += zc; }
     114             : 
     115       54875 :   diagH1 = new_chunk(lnz+1); /* diagH1[i] = 0 iff H[i,i] != 1 (set later) */
     116             : 
     117       54875 :   av = avma;
     118       54875 :   Cnew = cgetg(co, typ(C));
     119       54875 :   setlg(C, col+1); p1 = gmul(C,U);
     120      470522 :   for (j=1; j<=col; j++) gel(Cnew,j) = gel(p1,j);
     121     2296832 :   for (   ; j<co ; j++)  gel(Cnew,j) = gel(C,j);
     122             : 
     123             :   /* Clean up B using new H */
     124      261441 :   for (s=0,i=lnz; i; i--)
     125             :   {
     126      206566 :     GEN Di = gel(dep,i), Hi = gel(H,i);
     127      206566 :     GEN h = gel(Hi,i); /* H[i,i] */
     128      206566 :     if ( (diagH1[i] = is_pm1(h)) ) { h = NULL; s++; }
     129    15185242 :     for (j=col+1; j<co; j++)
     130             :     {
     131    14978676 :       GEN z = gel(B,j-col);
     132    14978676 :       p1 = gel(z,i+nlze);
     133    14978676 :       if (h) p1 = truedivii(p1,h);
     134    14978676 :       if (!signe(p1)) continue;
     135   150601711 :       for (k=1; k<=nlze; k++) gel(z,k) = subii(gel(z,k), mulii(p1, gel(Di,k)));
     136   116835394 :       for (   ; k<=lig;  k++) gel(z,k) = subii(gel(z,k), mulii(p1, gel(Hi,k-nlze)));
     137     8434263 :       gel(Cnew,j) = gsub(gel(Cnew,j), gmul(p1, gel(Cnew,i+zc)));
     138             :     }
     139      206566 :     if (gc_needed(av,2))
     140             :     {
     141         782 :       if(DEBUGMEM>1) pari_warn(warnmem,"hnffinal, i = %ld",i);
     142         782 :       gerepileall(av, 2, &Cnew, &B);
     143             :     }
     144             :   }
     145       54875 :   p1 = cgetg(lnz+1,t_VEC); p2 = perm + nlze;
     146      261441 :   for (i1=0, j1=lnz-s, i=1; i<=lnz; i++) /* push the 1 rows down */
     147      206566 :     if (diagH1[i])
     148       96166 :       gel(p1,++j1) = gel(p2,i);
     149             :     else
     150      110400 :       gel(p2,++i1) = gel(p2,i);
     151      151041 :   for (i=i1+1; i<=lnz; i++) gel(p2,i) = gel(p1,i);
     152             : 
     153             :   /* s = # extra redundant generators taken from H
     154             :    *          zc  col-s  co   zc = col - lnz
     155             :    *       [ 0 |dep |     ]    i = nlze + lnz - s = lig - s
     156             :    *  nlze [--------|  B' ]
     157             :    *       [ 0 | H' |     ]    H' = H minus the s rows with a 1 on diagonal
     158             :    *     i [--------|-----] lig-s           (= "1-rows")
     159             :    *       [   0    | Id  ]
     160             :    *       [        |     ] li */
     161       54875 :   lig -= s; col -= s; lnz -= s;
     162       54875 :   Hnew = cgetg(lnz+1,t_MAT);
     163       54875 :   depnew = cgetg(lnz+1,t_MAT); /* only used if nlze > 0 */
     164       54875 :   Bnew = cgetg(co-col,t_MAT);
     165       54875 :   C = shallowcopy(Cnew);
     166      261441 :   for (j=1,i1=j1=0; j<=lnz+s; j++)
     167             :   {
     168      206566 :     GEN z = gel(H,j);
     169      206566 :     if (diagH1[j])
     170             :     { /* hit exactly s times */
     171       96166 :       i1++; C[i1+col] = Cnew[j+zc];
     172       96166 :       p1 = cgetg(lig+1,t_COL); gel(Bnew,i1) = p1;
     173     1572729 :       for (i=1; i<=nlze; i++) gel(p1,i) = gcoeff(dep,i,j);
     174       96166 :       p1 += nlze;
     175             :     }
     176             :     else
     177             :     {
     178      110400 :       j1++; C[j1+zc] = Cnew[j+zc];
     179      110400 :       p1 = cgetg(lnz+1,t_COL); gel(Hnew,j1) = p1;
     180      110400 :       depnew[j1] = dep[j];
     181             :     }
     182     1344960 :     for (i=k=1; k<=lnz; i++)
     183     1138394 :       if (!diagH1[i]) p1[k++] = z[i];
     184             :   }
     185     2296832 :   for (j=s+1; j<co-col; j++)
     186             :   {
     187     2241957 :     GEN z = gel(B,j-s);
     188     2241957 :     p1 = cgetg(lig+1,t_COL); gel(Bnew,j) = p1;
     189    33921670 :     for (i=1; i<=nlze; i++) gel(p1,i) = gel(z,i);
     190     2241957 :     z += nlze; p1 += nlze;
     191    12400434 :     for (i=k=1; k<=lnz; i++)
     192    10158477 :       if (!diagH1[i]) gel(p1,k++) = gel(z,i);
     193             :   }
     194       54875 :   *ptdep = depnew;
     195       54875 :   *ptC = C;
     196       54875 :   *ptB = Bnew; return Hnew;
     197             : }
     198             : 
     199             : /* for debugging */
     200             : static void
     201           0 : p_mat(GEN mat, GEN perm, long k)
     202             : {
     203           0 :   pari_sp av = avma;
     204           0 :   perm = vecslice(perm, k+1, lg(perm)-1);
     205           0 :   err_printf("Permutation: %Ps\n",perm);
     206           0 :   if (DEBUGLEVEL > 6)
     207           0 :     err_printf("matgen = %Ps\n", zm_to_ZM( rowpermute(mat, perm) ));
     208           0 :   set_avma(av);
     209           0 : }
     210             : 
     211             : static GEN
     212      703618 : col_dup(long l, GEN col)
     213             : {
     214      703618 :   GEN c = new_chunk(l);
     215      703618 :   memcpy(c,col,l * sizeof(long)); return c;
     216             : }
     217             : 
     218             : /* permutation giving imagecompl(x') | image(x'), x' = transpose of x */
     219             : static GEN
     220       54875 : ZM_rowrankprofile(GEN x, long *nlze)
     221             : {
     222       54875 :   pari_sp av = avma;
     223             :   GEN d, y;
     224             :   long i, j, k, l, r;
     225             : 
     226       54875 :   x = shallowtrans(x); l = lg(x);
     227       54875 :   (void)new_chunk(l); /* HACK */
     228       54875 :   d = ZM_pivots(x,&r); set_avma(av);
     229       54875 :   *nlze = r;
     230       54875 :   if (!d) return identity_perm(l-1);
     231       51211 :   y = cgetg(l,t_VECSMALL);
     232      422186 :   for (i = j = 1, k = r+1; i<l; i++)
     233      370975 :     if (d[i]) y[k++] = i; else y[j++] = i;
     234       51211 :   return y;
     235             : }
     236             : 
     237             : /* HNF reduce a relation matrix (column operations + row permutation)
     238             : ** Input:
     239             : **   mat = (li-1) x (co-1) matrix of long
     240             : **   C   = r x (co-1) matrix of GEN
     241             : **   perm= permutation vector (length li-1), indexing the rows of mat: easier
     242             : **     to maintain perm than to copy rows. For columns we can do it directly
     243             : **     using e.g. swap(mat[i], mat[j])
     244             : **   k0 = integer. The k0 first lines of mat are dense, the others are sparse.
     245             : ** Output: cf ASCII art in the function body
     246             : **
     247             : ** row permutations applied to perm
     248             : ** column operations applied to C. IN PLACE
     249             : **/
     250             : GEN
     251       28543 : hnfspec_i(GEN mat0, GEN perm, GEN* ptdep, GEN* ptB, GEN* ptC, long k0)
     252             : {
     253             :   pari_sp av;
     254             :   long co, n, s, nlze, lnz, nr, i, j, k, lk0, col, lig, *p;
     255             :   GEN mat;
     256             :   GEN p1, p2, matb, matbnew, vmax, matt, T, extramat, B, C, H, dep, permpro;
     257       28543 :   const long li = lg(perm); /* = lgcols(mat0) */
     258       28543 :   const long CO = lg(mat0);
     259             : 
     260       28543 :   n = 0; /* -Wall */
     261             : 
     262       28543 :   C = *ptC; co = CO;
     263       28543 :   if (co > 300 && co > 1.5 * li)
     264             :   { /* treat the rest at the end */
     265           0 :     co = (long)(1.2 * li);
     266           0 :     setlg(C, co);
     267             :   }
     268             : 
     269       28543 :   if (DEBUGLEVEL>5)
     270             :   {
     271           0 :     err_printf("Entering hnfspec\n");
     272           0 :     p_mat(mat0,perm,0);
     273             :   }
     274       28543 :   matt = cgetg(co, t_MAT); /* dense part of mat (top) */
     275       28543 :   mat = cgetg(co, t_MAT);
     276      732161 :   for (j = 1; j < co; j++)
     277             :   {
     278      703618 :     GEN matj = col_dup(li, gel(mat0,j));
     279      703618 :     p1 = cgetg(k0+1,t_COL); gel(matt,j) = p1; gel(mat,j) = matj;
     280     4103933 :     for (i=1; i<=k0; i++) gel(p1,i) = stoi(matj[perm[i]]);
     281             :   }
     282       28543 :   av = avma;
     283             : 
     284       28543 :   i = lig = li-1; col = co-1; lk0 = k0;
     285       28543 :   T = (k0 || (lg(C) > 1 && lgcols(C) > 1))? matid(col): NULL;
     286             :   /* Look for lines with a single non-0 entry, equal to 1 in absolute value */
     287     3316830 :   while (i > lk0 && col)
     288     3288287 :     switch( count(mat,perm[i],col,&n) )
     289             :     {
     290       29995 :       case 0: /* move zero lines between k0+1 and lk0 */
     291       29995 :         lk0++; lswap(perm[i], perm[lk0]);
     292       29995 :         i = lig; continue;
     293             : 
     294      146991 :       case 1: /* move trivial generator between lig+1 and li */
     295      146991 :         lswap(perm[i], perm[lig]);
     296      146991 :         if (T) swap(gel(T,n), gel(T,col));
     297      146991 :         swap(gel(mat,n), gel(mat,col)); p = gel(mat,col);
     298      146991 :         if (p[perm[lig]] < 0) /* = -1 */
     299             :         { /* convert relation -g = 0 to g = 0 */
     300     5997406 :           for (i=lk0+1; i<lig; i++) p[perm[i]] = -p[perm[i]];
     301      113404 :           if (T)
     302             :           {
     303      113404 :             p1 = gel(T,col);
     304      113404 :             for (i=1; ; i++) /* T is a permuted identity: single non-0 entry */
     305     5854882 :               if (signe(gel(p1,i))) { togglesign_safe(&gel(p1,i)); break; }
     306             :           }
     307             :         }
     308      146991 :         lig--; col--; i = lig; continue;
     309             : 
     310     3111301 :       default: i--;
     311             :     }
     312       28543 :   if (DEBUGLEVEL>5) { err_printf("    after phase1:\n"); p_mat(mat,perm,0); }
     313             : 
     314             : #define absmax(s,z) {long _z; _z = labs(z); if (_z > s) s = _z;}
     315             :   /* Get rid of all lines containing only 0 and +/- 1, keeping track of column
     316             :    * operations in T. Leave the rows 1..lk0 alone [up to k0, coefficient
     317             :    * explosion, between k0+1 and lk0, row is 0] */
     318       28543 :   s = 0;
     319      260629 :   while (lig > lk0 && col && s < (long)(HIGHBIT>>1))
     320             :   {
     321     1096720 :     for (i=lig; i>lk0; i--)
     322     1071496 :       if (count(mat,perm[i],col,&n) > 0) break;
     323      257310 :     if (i == lk0) break;
     324             : 
     325             :     /* only 0, +/- 1 entries, at least 2 of them non-zero */
     326      232086 :     lswap(perm[i], perm[lig]);
     327      232086 :     swap(gel(mat,n), gel(mat,col)); p = gel(mat,col);
     328      232086 :     if (T) swap(gel(T,n), gel(T,col));
     329      232086 :     if (p[perm[lig]] < 0)
     330             :     {
     331     5090502 :       for (i=lk0+1; i<=lig; i++) p[perm[i]] = -p[perm[i]];
     332      158109 :       if (T) ZV_togglesign(gel(T,col));
     333             :     }
     334    10757471 :     for (j=1; j<col; j++)
     335             :     {
     336    10525385 :       GEN matj = gel(mat,j);
     337             :       long t;
     338    10525385 :       if (! (t = matj[perm[lig]]) ) continue;
     339      633624 :       if (t == 1) {
     340    11721267 :         for (i=lk0+1; i<=lig; i++) absmax(s, matj[perm[i]] -= p[perm[i]]);
     341             :       }
     342             :       else { /* t = -1 */
     343    10608219 :         for (i=lk0+1; i<=lig; i++) absmax(s, matj[perm[i]] += p[perm[i]]);
     344             :       }
     345      633624 :       if (T) ZC_lincomb1_inplace(gel(T,j), gel(T,col), stoi(-t));
     346             :     }
     347      232086 :     lig--; col--;
     348      232086 :     if (gc_needed(av,3))
     349             :     {
     350           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"hnfspec[1]");
     351           0 :       if (T) T = gerepilecopy(av, T); else set_avma(av);
     352             :     }
     353             :   }
     354             :   /* As above with lines containing a +/- 1 (no other assumption).
     355             :    * Stop when single precision becomes dangerous */
     356       28543 :   vmax = cgetg(co,t_VECSMALL);
     357      353084 :   for (j=1; j<=col; j++)
     358             :   {
     359      324541 :     GEN matj = gel(mat,j);
     360     2635710 :     for (s=0, i=lk0+1; i<=lig; i++) absmax(s, matj[i]);
     361      324541 :     vmax[j] = s;
     362             :   }
     363      108919 :   while (lig > lk0 && col)
     364             :   {
     365      170946 :     for (i=lig; i>lk0; i--)
     366      164002 :       if ( (n = count2(mat,perm[i],col)) ) break;
     367       87325 :     if (i == lk0) break;
     368             : 
     369       80381 :     lswap(vmax[n], vmax[col]);
     370       80381 :     lswap(perm[i], perm[lig]);
     371       80381 :     swap(gel(mat,n), gel(mat,col)); p = gel(mat,col);
     372       80381 :     if (T) swap(gel(T,n), gel(T,col));
     373       80381 :     if (p[perm[lig]] < 0)
     374             :     {
     375      297682 :       for (i=lk0+1; i<=lig; i++) p[perm[i]] = -p[perm[i]];
     376       33762 :       if (T) ZV_togglesign(gel(T,col));
     377             :     }
     378     1440097 :     for (j=1; j<col; j++)
     379             :     {
     380     1359721 :       GEN matj = gel(mat,j);
     381             :       long t;
     382     1359721 :       if (! (t = matj[perm[lig]]) ) continue;
     383      728083 :       if (vmax[col] && (ulong)labs(t) >= (HIGHBIT-vmax[j]) / vmax[col])
     384           5 :         goto END2;
     385             : 
     386    14154178 :       for (s=0, i=lk0+1; i<=lig; i++) absmax(s, matj[perm[i]] -= t*p[perm[i]]);
     387      728078 :       vmax[j] = s;
     388      728078 :       if (T) ZC_lincomb1_inplace(gel(T,j), gel(T,col), stoi(-t));
     389             :     }
     390       80376 :     lig--; col--;
     391       80376 :     if (gc_needed(av,3))
     392             :     {
     393           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"hnfspec[2]");
     394           0 :       gerepileall(av, T? 2: 1, &vmax, &T);
     395             :     }
     396             :   }
     397             : 
     398       21594 : END2: /* clean up mat: remove everything to the right of the 1s on diagonal */
     399             :   /* go multiprecision first */
     400       28543 :   matb = cgetg(co,t_MAT); /* bottom part (complement of matt) */
     401      732161 :   for (j=1; j<co; j++)
     402             :   {
     403      703618 :     GEN matj = gel(mat,j);
     404      703618 :     p1 = cgetg(li-k0,t_COL); gel(matb,j) = p1;
     405      703618 :     p1 -= k0;
     406    46009556 :     for (i=k0+1; i<li; i++) gel(p1,i) = stoi(matj[perm[i]]);
     407             :   }
     408       28543 :   if (DEBUGLEVEL>5)
     409             :   {
     410           0 :     err_printf("    after phase2:\n");
     411           0 :     p_mat(mat,perm,lk0);
     412             :   }
     413      459585 :   for (i=li-2; i>lig; i--)
     414             :   {
     415      431042 :     long h, i0 = i - k0, k = i + co-li;
     416      431042 :     GEN Bk = gel(matb,k);
     417    16788523 :     for (j=k+1; j<co; j++)
     418             :     {
     419    16357481 :       GEN Bj = gel(matb,j), v = gel(Bj,i0);
     420    16357481 :       s = signe(v); if (!s) continue;
     421             : 
     422     3281756 :       gel(Bj,i0) = gen_0;
     423     3281756 :       if (is_pm1(v))
     424             :       {
     425     1650202 :         if (s > 0) /* v = 1 */
     426    27975802 :         { for (h=1; h<i0; h++) gel(Bj,h) = subii(gel(Bj,h), gel(Bk,h)); }
     427             :         else /* v = -1 */
     428    24946171 :         { for (h=1; h<i0; h++) gel(Bj,h) = addii(gel(Bj,h), gel(Bk,h)); }
     429             :       }
     430             :       else {
     431    57296319 :         for (h=1; h<i0; h++) gel(Bj,h) = subii(gel(Bj,h), mulii(v,gel(Bk,h)));
     432             :       }
     433     3281756 :       if (T) ZC_lincomb1_inplace(gel(T,j), gel(T,k), negi(v));
     434             :     }
     435      431042 :     if (gc_needed(av,2))
     436             :     {
     437          36 :       if(DEBUGMEM>1) pari_warn(warnmem,"hnfspec[3], i = %ld", i);
     438       17592 :       for (h=1; h<co; h++) setlg(matb[h], i0+1); /* bottom can be forgotten */
     439          36 :       gerepileall(av, T? 2: 1, &matb, &T);
     440             :     }
     441             :   }
     442      732161 :   for (j=1; j<co; j++) setlg(matb[j], lig-k0+1); /* bottom can be forgotten */
     443       28543 :   gerepileall(av, T? 2: 1, &matb, &T);
     444       28543 :   if (DEBUGLEVEL>5) err_printf("    matb cleaned up (using Id block)\n");
     445             : 
     446       28543 :   nlze = lk0 - k0;  /* # of 0 rows */
     447       28543 :   lnz = lig-nlze+1; /* 1 + # of non-0 rows (!= 0...0 1 0 ... 0) */
     448       28543 :   if (T) matt = ZM_mul(matt,T); /* update top rows */
     449       28543 :   extramat = cgetg(col+1,t_MAT); /* = new C minus the 0 rows */
     450      272708 :   for (j=1; j<=col; j++)
     451             :   {
     452      244165 :     GEN z = gel(matt,j);
     453      244165 :     GEN t = (gel(matb,j)) + nlze - k0;
     454      244165 :     p2=cgetg(lnz,t_COL); gel(extramat,j) = p2;
     455     1171701 :     for (i=1; i<=k0; i++) gel(p2,i) = gel(z,i); /* top k0 rows */
     456      619650 :     for (   ; i<lnz; i++) gel(p2,i) = gel(t,i); /* other non-0 rows */
     457             :   }
     458       28543 :   if (!col) {
     459         336 :     permpro = identity_perm(lnz);
     460         336 :     nr = lnz;
     461             :   }
     462             :   else
     463       28207 :     permpro = ZM_rowrankprofile(extramat, &nr);
     464             :   /* lnz = lg(permpro) */
     465       28543 :   if (nlze)
     466             :   { /* put the nlze 0 rows (trivial generators) at the top */
     467        4116 :     p1 = new_chunk(lk0+1);
     468       34111 :     for (i=1; i<=nlze; i++) p1[i] = perm[i + k0];
     469       27015 :     for (   ; i<=lk0; i++)  p1[i] = perm[i - nlze];
     470       57010 :     for (i=1; i<=lk0; i++)  perm[i] = p1[i];
     471             :   }
     472             :   /* sort other rows according to permpro (nr redundant generators first) */
     473       28543 :   p1 = new_chunk(lnz); p2 = perm + nlze;
     474      137112 :   for (i=1; i<lnz; i++) p1[i] = p2[permpro[i]];
     475      137112 :   for (i=1; i<lnz; i++) p2[i] = p1[i];
     476             :   /* perm indexes the rows of mat
     477             :    *   |_0__|__redund__|__dense__|__too big__|_____done______|
     478             :    *   0  nlze                              lig             li
     479             :    *         \___nr___/ \___k0__/
     480             :    *         \____________lnz ______________/
     481             :    *
     482             :    *               col   co
     483             :    *       [dep     |     ]
     484             :    *    i0 [--------|  B  ] (i0 = nlze + nr)
     485             :    *       [matbnew |     ] matbnew has maximal rank = lnz-1 - nr
     486             :    * mat = [--------|-----] lig
     487             :    *       [   0    | Id  ]
     488             :    *       [        |     ] li */
     489             : 
     490       28543 :   matbnew = cgetg(col+1,t_MAT); /* dense+toobig, maximal rank. For hnffinal */
     491       28543 :   dep    = cgetg(col+1,t_MAT); /* rows dependent from the ones in matbnew */
     492      272708 :   for (j=1; j<=col; j++)
     493             :   {
     494      244165 :     GEN z = gel(extramat,j);
     495      244165 :     p1 = cgetg(nlze+nr+1,t_COL); gel(dep,j) = p1;
     496      244165 :     p2 = cgetg(lnz-nr,t_COL); gel(matbnew,j) = p2;
     497      577667 :     for (i=1; i<=nlze; i++) gel(p1,i) = gen_0;
     498      312300 :     p1 += nlze; for (i=1; i<=nr; i++) p1[i] = z[permpro[i]];
     499     1479051 :     p2 -= nr;   for (   ; i<lnz; i++) p2[i] = z[permpro[i]];
     500             :   }
     501             : 
     502             :   /* redundant generators in terms of the genuine generators
     503             :    * (x_i) = - (g_i) B */
     504       28543 :   B = cgetg(co-col,t_MAT);
     505      487996 :   for (j=col+1; j<co; j++)
     506             :   {
     507      459453 :     GEN y = gel(matt,j);
     508      459453 :     GEN z = gel(matb,j);
     509      459453 :     p1=cgetg(lig+1,t_COL); gel(B,j-col) = p1;
     510     3066557 :     for (i=1; i<=nlze; i++) gel(p1,i) = gel(z,i);
     511      459453 :     p1 += nlze; z += nlze-k0;
     512     4751224 :     for (k=1; k<lnz; k++)
     513             :     {
     514     4291771 :       i = permpro[k];
     515     4291771 :       gel(p1,k) = (i <= k0)? gel(y,i): gel(z,i);
     516             :     }
     517             :   }
     518       28543 :   if (T) C = typ(C)==t_MAT? RgM_ZM_mul(C,T): RgV_RgM_mul(C,T);
     519       28543 :   gerepileall(av, 4, &matbnew, &B, &dep, &C);
     520       28543 :   *ptdep = dep;
     521       28543 :   *ptB = B;
     522       28543 :   H = hnffinal(matbnew, perm, ptdep, ptB, &C);
     523       28543 :   if (CO > co)
     524             :   { /* treat the rest, N cols at a time (hnflll slow otherwise) */
     525           0 :     const long N = 300;
     526           0 :     long a, L = CO - co, l = minss(L, N); /* L columns to add */
     527           0 :     GEN CC = *ptC, m0 = mat0;
     528           0 :     setlg(CC, CO); /* restore */
     529           0 :     CC += co-1;
     530           0 :     m0 += co-1;
     531           0 :     for (a = l;;)
     532           0 :     {
     533             :       GEN MAT, emb;
     534           0 :       gerepileall(av, 4, &H,&C,ptB,ptdep);
     535           0 :       MAT = cgetg(l + 1, t_MAT);
     536           0 :       emb = cgetg(l + 1, typ(C));
     537           0 :       for (j = 1 ; j <= l; j++)
     538             :       {
     539           0 :         gel(MAT,j) = gel(m0,j);
     540           0 :         emb[j] = CC[j];
     541             :       }
     542           0 :       H = hnfadd_i(H, perm, ptdep, ptB, &C, MAT, emb);
     543           0 :       if (a == L) break;
     544           0 :       CC += l;
     545           0 :       m0 += l;
     546           0 :       a += l; if (a > L) { l = L - (a - l); a = L; }
     547             :     }
     548             :   }
     549       28543 :   *ptC = C; return H;
     550             : }
     551             : 
     552             : GEN
     553           0 : hnfspec(GEN mat, GEN perm, GEN* ptdep, GEN* ptB, GEN* ptC, long k0)
     554             : {
     555           0 :   pari_sp av = avma;
     556           0 :   GEN H = hnfspec_i(mat, perm, ptdep, ptB, ptC, k0);
     557           0 :   gerepileall(av, 4, ptC, ptdep, ptB, &H); return H;
     558             : }
     559             : 
     560             : /* HNF reduce x, apply same transforms to C */
     561             : GEN
     562           0 : mathnfspec(GEN x, GEN *pperm, GEN *pdep, GEN *pB, GEN *pC)
     563             : {
     564           0 :   long i, j, k, l, n, ly, lx = lg(x);
     565             :   GEN z, v1, perm;
     566           0 :   if (lx == 1) return cgetg(1, t_MAT);
     567           0 :   ly = lgcols(x);
     568           0 :   *pperm = perm = identity_perm(ly-1);
     569           0 :   z = cgetg(lx,t_MAT);
     570           0 :   for (i=1; i<lx; i++)
     571             :   {
     572           0 :     GEN C = cgetg(ly,t_COL), D = gel(x,i);
     573           0 :     gel(z,i) = C;
     574           0 :     for (j=1; j<ly; j++)
     575             :     {
     576           0 :       GEN d = gel(D,j);
     577           0 :       if (is_bigint(d)) goto TOOLARGE;
     578           0 :       C[j] = itos(d);
     579             :     }
     580             :   }
     581             :   /*  [ dep |     ]
     582             :    *  [-----|  B  ]
     583             :    *  [  H  |     ]
     584             :    *  [-----|-----]
     585             :    *  [  0  | Id  ] */
     586           0 :   return hnfspec(z,perm, pdep, pB, pC, 0);
     587             : 
     588           0 : TOOLARGE:
     589           0 :   if (lg(*pC) > 1 && lgcols(*pC) > 1)
     590           0 :     pari_err_IMPL("mathnfspec with large entries");
     591           0 :   x = ZM_hnf(x); lx = lg(x);
     592           0 :   v1 = cgetg(ly, t_VECSMALL);
     593           0 :   n = lx - ly;
     594           0 :   for (i = k = l = 1; i < ly; i++)
     595           0 :     if (equali1(gcoeff(x,i,i + n))) v1[l++] = i; else perm[k++] = i;
     596           0 :   setlg(perm, k);
     597           0 :   setlg(v1, l);
     598           0 :   x = rowpermute(x, perm); /* upper part */
     599           0 :   *pperm = vecsmall_concat(perm, v1);
     600           0 :   *pB = vecslice(x, k+n, lx-1);
     601           0 :   setlg(x, k);
     602           0 :   *pdep = rowslice(x, 1, n);
     603           0 :   return n? rowslice(x, n+1, k-1): x; /* H */
     604             : }
     605             : 
     606             : /* add new relations to a matrix treated by hnfspec (extramat / extraC) */
     607             : GEN
     608       26668 : hnfadd_i(GEN H, GEN perm, GEN* ptdep, GEN* ptB, GEN* ptC, /* cf hnfspec */
     609             :        GEN extramat,GEN extraC)
     610             : {
     611       26668 :   GEN matb, extratop, Cnew, permpro, B = *ptB, C = *ptC, dep = *ptdep;
     612             :   long i, lH, lB, li, lig, co, col, nlze;
     613             : 
     614       26668 :   if (lg(extramat) == 1) return H;
     615       26668 :   co = lg(C)-1;
     616       26668 :   lH = lg(H)-1;
     617       26668 :   lB = lg(B)-1;
     618       26668 :   li = lg(perm)-1;
     619       26668 :   lig = li - lB;
     620       26668 :   col = co - lB;
     621       26668 :   nlze = lig - lH;
     622             : 
     623             :  /*               col    co
     624             :   *       [ 0 |dep |     ]
     625             :   *  nlze [--------|  B  ]
     626             :   *       [ 0 | H  |     ]
     627             :   *       [--------|-----] lig
     628             :   *       [   0    | Id  ]
     629             :   *       [        |     ] li */
     630       26668 :   extratop = zm_to_ZM( rowslicepermute(extramat, perm, 1, lig) );
     631       26668 :   if (li != lig)
     632             :   { /* zero out bottom part, using the Id block */
     633       26668 :     GEN A = vecslice(C, col+1, co);
     634       26668 :     GEN c = rowslicepermute(extramat, perm, lig+1, li);
     635       26668 :     extraC   = gsub(extraC, typ(A)==t_MAT? RgM_zm_mul(A, c): RgV_zm_mul(A,c));
     636       26668 :     extratop = ZM_sub(extratop, ZM_zm_mul(B, c));
     637             :   }
     638             : 
     639       26668 :   extramat = shallowconcat(extratop, vconcat(dep, H));
     640       26668 :   Cnew     = shallowconcat(extraC, vecslice(C, col-lH+1, co));
     641       26668 :   if (DEBUGLEVEL>5) err_printf("    1st phase done\n");
     642       26668 :   permpro = ZM_rowrankprofile(extramat, &nlze);
     643       26668 :   extramat = rowpermute(extramat, permpro);
     644       26668 :   *ptB     = rowpermute(B,        permpro);
     645       26668 :   permpro = vecsmallpermute(perm, permpro);
     646      289816 :   for (i=1; i<=lig; i++) perm[i] = permpro[i]; /* perm o= permpro */
     647             : 
     648       26668 :   *ptdep  = rowslice(extramat, 1, nlze);
     649       26668 :   matb    = rowslice(extramat, nlze+1, lig);
     650       26668 :   if (DEBUGLEVEL>5) err_printf("    2nd phase done\n");
     651       26668 :   H = hnffinal(matb,perm,ptdep,ptB,&Cnew);
     652       26668 :   *ptC = shallowconcat(vecslice(C, 1, col-lH), Cnew);
     653       26668 :   return H;
     654             : }
     655             : 
     656             : GEN
     657           0 : hnfadd(GEN H, GEN perm, GEN* ptdep, GEN* ptB, GEN* ptC, /* cf hnfspec */
     658             :        GEN extramat,GEN extraC)
     659             : {
     660           0 :   pari_sp av = avma;
     661           0 :   H = hnfadd_i(H, perm, ptdep, ptB, ptC, ZM_to_zm(extramat), extraC);
     662           0 :   gerepileall(av, 4, ptC, ptdep, ptB, &H); return H;
     663             : }
     664             : 
     665             : /* zero aj = Aij (!= 0)  using  ak = Aik (maybe 0), via linear combination of
     666             :  * A[j] and A[k] of determinant 1. If U != NULL, likewise update its columns */
     667             : static void
     668    28867012 : ZC_elem(GEN aj, GEN ak, GEN A, GEN U, long j, long k)
     669             : {
     670             :   GEN p1,u,v,d;
     671             : 
     672    28867012 :   if (!signe(ak)) {
     673       19150 :     swap(gel(A,j), gel(A,k));
     674       19150 :     if (U) swap(gel(U,j), gel(U,k));
     675    25481826 :     return;
     676             :   }
     677    28847862 :   d = bezout(aj,ak,&u,&v);
     678             :   /* frequent special case (u,v) = (1,0) or (0,1) */
     679    28846971 :   if (!signe(u))
     680             :   { /* ak | aj */
     681    11861465 :     p1 = diviiexact(aj,ak); togglesign(p1);
     682    11861748 :     ZC_lincomb1_inplace(gel(A,j), gel(A,k), p1);
     683    11862090 :     if (U)
     684      540389 :       ZC_lincomb1_inplace(gel(U,j), gel(U,k), p1);
     685    11862126 :     return;
     686             :   }
     687    16985506 :   if (!signe(v))
     688             :   { /* aj | ak */
     689    13600616 :     p1 = diviiexact(ak,aj); togglesign(p1);
     690    13600509 :     ZC_lincomb1_inplace(gel(A,k), gel(A,j), p1);
     691    13600550 :     swap(gel(A,j), gel(A,k));
     692    13600550 :     if (U) {
     693       62567 :       ZC_lincomb1_inplace(gel(U,k), gel(U,j), p1);
     694       62567 :       swap(gel(U,j), gel(U,k));
     695             :     }
     696    13600550 :     return;
     697             :   }
     698             : 
     699     3384890 :   if (!is_pm1(d)) { aj = diviiexact(aj, d); ak = diviiexact(ak, d); }
     700     3384885 :   p1 = gel(A,k); aj = negi(aj); /* NOT togglesign */
     701     3384916 :   gel(A,k) = ZC_lincomb(u,v, gel(A,j),p1);
     702     3384876 :   gel(A,j) = ZC_lincomb(aj,ak, p1,gel(A,j));
     703     3384889 :   if (U)
     704             :   {
     705      186703 :     p1 = gel(U,k);
     706      186703 :     gel(U,k) = ZC_lincomb(u,v, gel(U,j),p1);
     707      186695 :     gel(U,j) = ZC_lincomb(aj,ak, p1,gel(U,j));
     708             :   }
     709             : }
     710             : 
     711             : INLINE int
     712        1001 : is_RgX(GEN a, long v) { return typ(a) == t_POL && varn(a)==v; }
     713             : /* set u,v such that au + bv = gcd(a,b), divide a,b by the gcd */
     714             : static GEN
     715         266 : gbezout_step(GEN *pa, GEN *pb, GEN *pu, GEN *pv, long vx)
     716             : {
     717         266 :   GEN a = *pa, b = *pb, d;
     718         266 :   if (gequal0(a))
     719             :   {
     720           0 :     *pa = gen_0; *pu = gen_0;
     721           0 :     *pb = gen_1; *pv = gen_1; return b;
     722             :   }
     723         266 :   a = is_RgX(a,vx)? RgX_renormalize(a): scalarpol(a, vx);
     724         266 :   b = is_RgX(b,vx)? RgX_renormalize(b): scalarpol(b, vx);
     725         266 :   d = RgX_extgcd(a,b, pu,pv);
     726         266 :   if (degpol(d)) { a = RgX_div(a, d); b = RgX_div(b, d); }
     727         154 :   else if (typ(gel(d,2)) == t_REAL && lg(gel(d,2)) <= 3)
     728             : #if 1
     729             :   { /* possible accuracy problem */
     730           0 :     GEN D = RgX_gcd_simple(a,b);
     731           0 :     if (degpol(D)) {
     732           0 :       D = RgX_normalize(D);
     733           0 :       a = RgX_div(a, D);
     734           0 :       b = RgX_div(b, D);
     735           0 :       d = RgX_extgcd(a,b, pu,pv); /* retry now */
     736           0 :       d = RgX_mul(d, D);
     737             :     }
     738             :   }
     739             : #else
     740             :   { /* less stable */
     741             :     d = RgX_extgcd_simple(a,b, pu,pv);
     742             :     if (degpol(d)) { a = RgX_div(a, d); b = RgX_div(b, d); }
     743             :   }
     744             : #endif
     745         266 :   *pa = a;
     746         266 :   *pb = b; return d;
     747             : }
     748             : static GEN
     749      465392 : col_mul(GEN x, GEN c)
     750             : {
     751      465392 :   if (typ(x) == t_INT)
     752             :   {
     753      465392 :     long s = signe(x);
     754      465392 :     if (!s) return NULL;
     755      353255 :     if (is_pm1(x)) return (s > 0)? c: RgC_neg(c);
     756             :   }
     757       75407 :   return RgC_Rg_mul(c, x);
     758             : }
     759             : static void
     760           0 : do_zero(GEN x)
     761             : {
     762           0 :   long i, lx = lg(x);
     763           0 :   for (i=1; i<lx; i++) gel(x,i) = gen_0;
     764           0 : }
     765             : 
     766             : /* (c1, c2) *= [u,-b; v,a] */
     767             : static void
     768      116348 : update(GEN u, GEN v, GEN a, GEN b, GEN *c1, GEN *c2)
     769             : {
     770             :   GEN p1,p2;
     771             : 
     772      116348 :   u = col_mul(u,*c1);
     773      116348 :   v = col_mul(v,*c2);
     774      116348 :   if (u) p1 = v? gadd(u,v): u;
     775        1541 :   else   p1 = v? v: NULL;
     776             : 
     777      116348 :   a = col_mul(a,*c2);
     778      116348 :   b = col_mul(gneg_i(b),*c1);
     779      116348 :   if (a) p2 = b? RgC_add(a,b): a;
     780           0 :   else   p2 = b? b: NULL;
     781             : 
     782      116348 :   if (!p1) do_zero(*c1); else *c1 = p1;
     783      116348 :   if (!p2) do_zero(*c2); else *c2 = p2;
     784      116348 : }
     785             : 
     786             : /* zero aj = Aij (!= 0)  using  ak = Aik (maybe 0), via linear combination of
     787             :  * A[j] and A[k] of determinant 1. If U != NULL, likewise update its columns */
     788             : static void
     789           7 : RgC_elem(GEN aj, GEN ak, GEN A, GEN V, long j, long k, long li, long vx)
     790             : {
     791           7 :   GEN u,v, d = gbezout_step(&aj, &ak, &u, &v, vx);
     792             :   long l;
     793             :   /* (A[,k], A[,j]) *= [v, -aj; u, ak ] */
     794          14 :   for (l = 1; l < li; l++)
     795             :   {
     796           7 :     GEN t = gadd(gmul(u,gcoeff(A,l,j)), gmul(v,gcoeff(A,l,k)));
     797           7 :     gcoeff(A,l,j) = gsub(gmul(ak,gcoeff(A,l,j)), gmul(aj,gcoeff(A,l,k)));
     798           7 :     gcoeff(A,l,k) = t;
     799             :   }
     800           7 :   gcoeff(A,li,j) = gen_0;
     801           7 :   gcoeff(A,li,k) = d;
     802           7 :   if (V) update(v,u,ak,aj,(GEN*)(V+k),(GEN*)(V+j));
     803           7 : }
     804             : 
     805             : /* reduce A[i,j] mod A[i,j0] for j=j0+1... via column operations */
     806             : static void
     807     1418559 : ZM_reduce(GEN A, GEN U, long i, long j0)
     808             : {
     809     1418559 :   long j, lA = lg(A);
     810     1418559 :   GEN d = gcoeff(A,i,j0);
     811     1418559 :   if (signe(d) < 0)
     812             :   {
     813        5444 :     ZV_neg_inplace(gel(A,j0));
     814        5444 :     if (U) ZV_togglesign(gel(U,j0));
     815        5444 :     d = gcoeff(A,i,j0);
     816             :   }
     817     3348496 :   for (j=j0+1; j<lA; j++)
     818             :   {
     819     1929937 :     GEN q = truedivii(gcoeff(A,i,j), d);
     820     1929937 :     if (!signe(q)) continue;
     821             : 
     822      103809 :     togglesign(q);
     823      103809 :     ZC_lincomb1_inplace(gel(A,j), gel(A,j0), q);
     824      103809 :     if (U) ZC_lincomb1_inplace(gel(U,j), gel(U,j0), q);
     825             :   }
     826     1418559 : }
     827             : 
     828             : /* normalize T as if it were a t_POL in variable v */
     829             : static GEN
     830         210 : normalize_as_RgX(GEN T, long v, GEN *pd)
     831             : {
     832             :   GEN d;
     833         210 :   if (!is_RgX(T,v)) { *pd = T; return gen_1; }
     834         175 :   d = leading_coeff(T);
     835         175 :   while (gequal0(d) || (typ(d) == t_REAL && lg(d) == 3
     836           0 :                         && gexpo(T) - expo(d) > (long)BITS_IN_LONG)) {
     837          14 :     T = normalizepol_lg(T, lg(T)-1);
     838          14 :     if (!signe(T)) { *pd = gen_1; return T; }
     839           0 :     d = leading_coeff(T);
     840             :   }
     841         161 :   if (degpol(T)) T = RgX_Rg_div(T,d); else { d = gel(T,2); T = gen_1; }
     842         161 :   *pd = d; return T;
     843             : }
     844             : /* reduce A[i,j] mod A[i,j0] for j=j0+1... via column operations */
     845             : static void
     846          42 : RgM_reduce(GEN A, GEN U, long i, long j0, long vx)
     847             : {
     848          42 :   long j, lA = lg(A);
     849          42 :   GEN d, T = normalize_as_RgX(gcoeff(A,i,j0), vx, &d);
     850          42 :   if (U && !gequal1(d)) gel(U,j0) = RgC_Rg_div(gel(U,j0), d);
     851          42 :   gcoeff(A,i,j0) = T;
     852             : 
     853          49 :   for (j=j0+1; j<lA; j++)
     854             :   {
     855           7 :     GEN t = gcoeff(A,i,j), q;
     856           7 :     if (gequal0(t)) continue;
     857           7 :     if (T == gen_1)
     858           0 :       q = t;
     859           7 :     else if (is_RgX(t,vx))
     860           7 :       q = RgX_div(t, T);
     861           0 :     else continue;
     862             : 
     863           7 :     if (gequal0(q)) continue;
     864           0 :     gel(A,j) = RgC_sub(gel(A,j), RgC_Rg_mul(gel(A,j0), q));
     865           0 :     if (U) gel(U,j) = RgC_sub(gel(U,j), RgC_Rg_mul(gel(U,j0), q));
     866             :   }
     867          42 : }
     868             : 
     869             : /* A,B square integral in upper HNF, of the same dimension > 0. Return Au
     870             :  * in Z^n (v in Z^n not computed), such that Au + Bv = [1, 0, ..., 0] */
     871             : GEN
     872      172638 : hnfmerge_get_1(GEN A, GEN B)
     873             : {
     874      172638 :   pari_sp av = avma;
     875      172638 :   long j, k, l = lg(A), lb;
     876      172638 :   GEN b, U = cgetg(l + 1, t_MAT), C = cgetg(l + 1, t_VEC);
     877             : 
     878      172643 :   b = gcoeff(B,1,1); lb = lgefint(b);
     879      334575 :   for (j = 1; j < l; j++)
     880             :   {
     881             :     GEN t;
     882      334575 :     long c = j+1;
     883      334575 :     gel(U,j) = col_ei(l-1, j);
     884      334661 :     gel(U,c) = zerocol(l-1); /* dummy */
     885      334681 :     gel(C,j) = vecslice(gel(A,j), 1,j);
     886      334665 :     gel(C,c) = vecslice(gel(B,j), 1,j);
     887      942705 :     for (k = j; k > 0; k--)
     888             :     {
     889      608240 :       t = gcoeff(C,k,c);
     890      608240 :       if (gequal0(t)) continue;
     891      536727 :       setlg(C[c], k+1);
     892      536722 :       ZC_elem(t, gcoeff(C,k,k), C, U, c, k);
     893      536286 :       if (lgefint(gcoeff(C,k,k)) > lb) gel(C,k) = FpC_red(gel(C,k), b);
     894      536286 :       if (j > 4)
     895             :       {
     896       52905 :         GEN u = gel(U,k);
     897             :         long h;
     898      707177 :         for (h=1; h<l; h++)
     899      654272 :           if (lgefint(gel(u,h)) > lb) gel(u,h) = remii(gel(u,h), b);
     900             :       }
     901             :     }
     902      334465 :     if (j == 1)
     903      172587 :       t = gcoeff(C,1,1);
     904             :     else
     905             :     {
     906             :       GEN u;
     907      161878 :       t = bezout(gcoeff(C,1,1), b, &u, NULL); /* >= 0 */
     908      162055 :       if (signe(u) && !equali1(u)) gel(U,1) = ZC_Z_mul(gel(U,1), u);
     909      162049 :       gcoeff(C,1,1) = t;
     910             :     }
     911      334636 :     if (equali1(t)) break;
     912             :   }
     913      172637 :   if (j >= l) return NULL;
     914      172637 :   b = lcmii(gcoeff(A,1,1),b);
     915      172605 :   A = FpC_red(ZM_ZC_mul(A,gel(U,1)), b);
     916      172615 :   return gerepileupto(av, FpC_center(A, b, shifti(b,-1)));
     917             : }
     918             : 
     919             : /* remove the first r columns */
     920             : static void
     921       63498 : remove_0cols(long r, GEN *pA, GEN *pB, long remove)
     922             : {
     923       63498 :   GEN A = *pA, B = *pB;
     924       63498 :   long l = lg(A);
     925       63498 :   A += r; A[0] = evaltyp(t_MAT) | evallg(l-r);
     926       63498 :   if (B && remove == 2) { B += r; B[0] = A[0]; }
     927       63498 :   *pA = A; *pB = B;
     928       63498 : }
     929             : 
     930             : /* Inefficient compared to hnfall. 'remove' = throw away lin.dep columns */
     931             : static GEN
     932       18193 : hnf_i(GEN A, int remove)
     933             : {
     934       18193 :   pari_sp av0 = avma, av;
     935             :   long s, n, m, j, k, li, def, ldef;
     936             : 
     937       18193 :   RgM_dimensions(A, &m, &n);
     938       18193 :   if (!n) return cgetg(1,t_MAT);
     939       18137 :   av = avma;
     940       18137 :   A = RgM_shallowcopy(A);
     941       18137 :   def = n; ldef = (m>n)? m-n: 0;
     942       62167 :   for (li=m; li>ldef; li--)
     943             :   {
     944      103770 :     for (j=def-1; j; j--)
     945             :     {
     946       59740 :       GEN a = gcoeff(A,li,j);
     947       59740 :       if (!signe(a)) continue;
     948             : 
     949             :       /* zero a = Aij  using  b = Aik */
     950       31362 :       k = (j==1)? def: j-1;
     951       31362 :       ZC_elem(a,gcoeff(A,li,k), A,NULL, j,k);
     952       31362 :       if (gc_needed(av,1))
     953             :       {
     954           0 :         if (DEBUGMEM>1) pari_warn(warnmem,"ZM_hnf[1]. li=%ld",li);
     955           0 :         A = gerepilecopy(av, A);
     956             :       }
     957             :     }
     958       44030 :     s = signe(gcoeff(A,li,def));
     959       44030 :     if (s)
     960             :     {
     961       43652 :       if (s < 0) ZV_neg_inplace(gel(A,def));
     962       43652 :       ZM_reduce(A, NULL, li,def);
     963       43652 :       def--;
     964             :     }
     965             :     else
     966         378 :       if (ldef) ldef--;
     967       44030 :     if (gc_needed(av,1))
     968             :     {
     969           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"ZM_hnf[2]. li=%ld",li);
     970           0 :       A = gerepilecopy(av, A);
     971             :     }
     972             :   }
     973             :   /* rank A = n - def */
     974       18137 :   if (remove) { GEN B = NULL; remove_0cols(def, &A, &B, remove); }
     975       18137 :   return gerepileupto(av0, ZM_copy(A));
     976             : }
     977             : 
     978             : GEN
     979       22758 : ZM_hnf(GEN x) { return lg(x) > 8? ZM_hnfall(x, NULL, 1): hnf_i(x, 1); }
     980             : 
     981             : /* u*z[1..k] mod p, in place */
     982             : static void
     983     3743387 : FpV_Fp_mul_part_ip(GEN z, GEN u, GEN p, long k)
     984             : {
     985             :   long i;
     986     3743387 :   if (is_pm1(u)) {
     987       62685 :     if (signe(u) > 0) {
     988      136422 :       for (i = 1; i <= k; i++)
     989      104563 :         if (signe(gel(z,i))) gel(z,i) = modii(gel(z,i), p);
     990             :     } else {
     991       78107 :       for (i = 1; i <= k; i++)
     992       47281 :         if (signe(gel(z,i))) gel(z,i) = modii(negi(gel(z,i)), p);
     993             :     }
     994             :   }
     995             :   else {
     996    14041635 :     for (i = 1; i <= k; i++)
     997    10361105 :       if (signe(gel(z,i))) gel(z,i) = Fp_mul(u,gel(z,i), p);
     998             :   }
     999     3743215 : }
    1000             : static void
    1001    24991373 : FpV_red_part_ipvec(GEN z, GEN p, long k)
    1002             : {
    1003             :   long i;
    1004    75889970 :   for (i = 1; i <= k; i++) gel(z,i) = modii(gel(z,i), gel(p,i));
    1005    24991362 : }
    1006             : 
    1007             : /* return x * U, in echelon form (mod p^m), where (det(U),p) = 1.
    1008             :  * If early_abort is set, return NULL as soon as one pivot is 0 (mod p^m) */
    1009             : GEN
    1010       60595 : ZpM_echelon(GEN x, long early_abort, GEN p, GEN pm)
    1011             : {
    1012       60595 :   pari_sp av0 = avma, av;
    1013             :   long m, li, co, i, j, k, def, ldef;
    1014             : 
    1015       60595 :   co = lg(x); if (co == 1) return cgetg(1,t_MAT);
    1016       60595 :   li = lgcols(x);
    1017       60595 :   av = avma;
    1018       60595 :   x = RgM_shallowcopy(x);
    1019       60595 :   m = Z_pval(pm, p);
    1020             : 
    1021       60595 :   ldef = (li > co)? li - co: 0;
    1022      545575 :   for (def = co-1,i = li-1; i > ldef; i--)
    1023             :   {
    1024      485631 :     long vmin = LONG_MAX, kmin = 0;
    1025      485631 :     GEN umin = gen_0, pvmin, q;
    1026     4302309 :     for (k = 1; k <= def; k++)
    1027             :     {
    1028     3912298 :       GEN u = gcoeff(x,i,k);
    1029             :       long v;
    1030     3912298 :       if (!signe(u)) continue;
    1031     1821623 :       v = Z_pvalrem(u, p, &u);
    1032     1821623 :       if (v >= m) gcoeff(x,i,k) = gen_0;
    1033     1186117 :       else if (v < vmin) {
    1034      456577 :         vmin = v; kmin = k; umin = u;
    1035      456577 :         if (!vmin) break;
    1036             :       }
    1037             :     }
    1038      485631 :     if (!kmin)
    1039             :     {
    1040      200922 :       if (early_abort) return NULL;
    1041      200271 :       gcoeff(x,i,def) = gen_0;
    1042      200271 :       ldef--;
    1043      200271 :       if (ldef < 0) ldef = 0;
    1044      200271 :       continue;
    1045             :     }
    1046      284709 :     if (kmin != def) swap(gel(x,def), gel(x,kmin));
    1047      284709 :     q = vmin? powiu(p, m-vmin): pm;
    1048             :     /* pivot has valuation vmin */
    1049      284709 :     umin = modii(umin, q);
    1050      284709 :     if (!equali1(umin))
    1051      218281 :       FpV_Fp_mul_part_ip(gel(x,def), Fp_inv(umin,q), pm, i-1);
    1052      284709 :     gcoeff(x, i, def) = pvmin = powiu(p, vmin);
    1053     3100806 :     for (j = def-1; j; j--)
    1054             :     { /* zero x[i, 1..def-1] using x[i,def] = pvmin */
    1055     2816097 :       GEN t, a = gcoeff(x,i,j) = modii(gcoeff(x,i,j), pm);
    1056     2816097 :       if (!signe(a)) continue;
    1057             : 
    1058     1382898 :       t = diviiexact(a, pvmin); togglesign(t);
    1059     1382898 :       ZC_lincomb1_inplace(gel(x,j), gel(x,def), t);
    1060     1382898 :       if (gc_needed(av,1))
    1061             :       {
    1062          14 :         if (DEBUGMEM>1) pari_warn(warnmem,"ZpM_echelon. i=%ld",i);
    1063          14 :         x = gerepilecopy(av, x); pvmin = gcoeff(x,i,def);
    1064             :       }
    1065             :     }
    1066      284709 :     def--;
    1067             :   }
    1068       59944 :   if (co > li)
    1069             :   {
    1070           0 :     x += co - li;
    1071           0 :     x[0] = evaltyp(t_MAT) | evallg(li);
    1072             :   }
    1073       59944 :   return gerepilecopy(av0, x);
    1074             : }
    1075             : GEN
    1076     1026783 : zlm_echelon(GEN x, long early_abort, ulong p, ulong pm)
    1077             : {
    1078     1026783 :   pari_sp av0 = avma;
    1079             :   long li, co, i, j, k, def, ldef;
    1080             :   ulong m;
    1081             : 
    1082     1026783 :   co = lg(x); if (co == 1) return cgetg(1,t_MAT);
    1083     1026783 :   li = lgcols(x);
    1084     1026791 :   x = Flm_copy(x);
    1085     1026796 :   m = u_lval(pm, p);
    1086             : 
    1087     1026801 :   ldef = (li > co)? li - co: 0;
    1088     5590634 :   for (def = co-1,i = li-1; i > ldef; i--)
    1089             :   {
    1090     4638460 :     long vmin = LONG_MAX, kmin = 0;
    1091     4638460 :     ulong umin = 0, pvmin, q;
    1092    18322451 :     for (k = 1; k <= def; k++)
    1093             :     {
    1094    16430773 :       ulong u = ucoeff(x,i,k);
    1095             :       long v;
    1096    16430773 :       if (!u) continue;
    1097     7933995 :       v = u_lvalrem(u, p, &u);
    1098     7933992 :       if (v >= (long) m) ucoeff(x,i,k) = 0;
    1099     7933992 :       else if (v < vmin) {
    1100     4765038 :         vmin = v; kmin = k; umin = u;
    1101     4765038 :         if (!vmin) break;
    1102             :       }
    1103             :     }
    1104     4638457 :     if (!kmin)
    1105             :     {
    1106      311786 :       if (early_abort) return NULL;
    1107      237165 :       ucoeff(x,i,def) = 0;
    1108      237165 :       ldef--;
    1109      237165 :       if (ldef < 0) ldef = 0;
    1110      237165 :       continue;
    1111             :     }
    1112     4326671 :     if (kmin != def) swap(gel(x,def), gel(x,kmin));
    1113     4326671 :     q = vmin? upowuu(p, m-vmin): pm;
    1114             :     /* pivot has valuation vmin */
    1115     4326679 :     umin %= q;
    1116     4326679 :     if (umin != 1)
    1117     2276479 :       Flv_Fl_mul_part_inplace(gel(x,def), Fl_inv(umin,q), pm, i-1);
    1118     4326673 :     ucoeff(x, i, def) = pvmin = upowuu(p, vmin);
    1119    26201734 :     for (j = def-1; j; j--)
    1120             :     { /* zero x[i, 1..def-1] using x[i,def] = pvmin */
    1121    21875066 :       ulong t, a = ucoeff(x,i,j);
    1122    21875066 :       if (!a) continue;
    1123             : 
    1124    12775684 :       t = Fl_neg(a / pvmin, q);
    1125    12775682 :       Flc_lincomb1_inplace(gel(x,j), gel(x,def), t, pm);
    1126             :     }
    1127     4326668 :     def--;
    1128             :   }
    1129      952174 :   if (co > li)
    1130             :   {
    1131           0 :     x += co - li;
    1132           0 :     x[0] = evaltyp(t_MAT) | evallg(li);
    1133             :   }
    1134      952174 :   return gerepilecopy(av0, x);
    1135             : }
    1136             : 
    1137             : static int
    1138       15260 : ZV_allequal(GEN v)
    1139             : {
    1140       15260 :   long i, l = lg(v);
    1141       15260 :   if (l > 1)
    1142             :   {
    1143       15260 :     GEN x = gel(v,1);
    1144       27510 :     for (i = 2; i < l; i++) if (!equalii(x,gel(v,i))) return 0;
    1145             :   }
    1146        9926 :   return 1;
    1147             : }
    1148             : /* compute optimal D for hnfmod: x upper triangular */
    1149             : static GEN
    1150     3477472 : optimal_D(GEN x, GEN D)
    1151             : {
    1152     3477472 :   long i, n = nbrows(x);
    1153     3477473 :   GEN C = shallowcopy(D);
    1154     3477488 :   gel(C,1) = gcoeff(x,1,1);
    1155     3791342 :   for (i = 2; i < n; i++)
    1156             :   {
    1157     1246459 :     GEN c = mulii(gel(C,i-1), gcoeff(x,i,i));
    1158     1246441 :     if (signe(c) < 0) togglesign(c);
    1159     1246442 :     if (cmpii(c, gel(D,i)) >= 0) break;
    1160      313854 :     gel(C,i) = c;
    1161             :   }
    1162     3477487 :   return C;
    1163             : }
    1164             : 
    1165             : /* D = multiple of det x (usually detint(x)) or vector of positive moduli
    1166             :  * (compute hnf(x | D))
    1167             :  * flag & hnf_MODID: reduce mod D * matid [ otherwise as above ].
    1168             :  * flag & hnf_PART: don't reduce once diagonal is known
    1169             :  * flag & hnf_CENTER: centermod HNF (2|x[i,j]| <] x[i,i]) */
    1170             : GEN
    1171     3478000 : ZM_hnfmodall_i(GEN x, GEN D, long flag)
    1172             : {
    1173             :   pari_sp av;
    1174     3478000 :   const long center = (flag & hnf_CENTER);
    1175             :   long moddiag, modsame, nli, li, co, i, j, k, def, ldef;
    1176             :   GEN u, LDM;
    1177             : 
    1178     3478000 :   co = lg(x);
    1179     3478000 :   if (co == 1)
    1180             :   {
    1181         469 :     if (typ(D) == t_INT || lg(D) == 1) return cgetg(1,t_MAT);
    1182          14 :     x = diagonal_shallow(D); /* handle flags properly */
    1183          14 :     co = lg(x);
    1184             :   }
    1185     3477545 :   li = lgcols(x);
    1186     3477547 :   if (li == 1)
    1187             :   {
    1188          21 :     if (typ(D) != t_INT && lg(D) != li) pari_err_DIM("ZM_hnfmod");
    1189          21 :     return cgetg(1,t_MAT);
    1190             :   }
    1191     3477526 :   nli = li - 1;
    1192     3477526 :   modsame = typ(D)==t_INT;
    1193     3477526 :   if (!modsame)
    1194             :   {
    1195       15274 :     if (lg(D) != li) pari_err_DIM("ZM_hnfmod");
    1196       15260 :     if (ZV_allequal(D)) { modsame = 1; D = gel(D,1); }
    1197             :   }
    1198     3477512 :   moddiag = (flag & hnf_MODID) || !modsame;
    1199             :   /* modsame: triangularize mod fixed d*Id;
    1200             :    * moddiag: modulo diagonal matrix, else modulo multiple of determinant */
    1201             : 
    1202     3477512 :   if (modsame)
    1203             :   {
    1204     3472180 :     LDM = const_vecsmall(nli, 2*lgefint(D)-2);
    1205     3472194 :     D = const_vec(nli,D);
    1206             :   }
    1207             :   else
    1208             :   {
    1209        5332 :     LDM = cgetg(li, t_VECSMALL);
    1210       20412 :     for (i=1; i<li; i++) LDM[i] = lgefint(gel(D,i));
    1211             :   }
    1212     3477531 :   av = avma;
    1213     3477531 :   x = RgM_shallowcopy(x);
    1214             : 
    1215     3477523 :   ldef = 0;
    1216     3477523 :   if (li > co)
    1217             :   {
    1218       50057 :     ldef = li - co;
    1219       50057 :     if (!moddiag)
    1220           7 :       pari_err_DOMAIN("ZM_hnfmod","nb lines",">", strtoGENstr("nb columns"), x);
    1221             :   }
    1222    13462864 :   for (def = co-1,i = nli; i > ldef; i--,def--)
    1223             :   {
    1224     9985402 :     GEN d = gel(D,i);
    1225     9985402 :     long add_N = modsame;
    1226    54362311 :     for (j = 1; j < def; j++)
    1227             :     {
    1228    44376963 :       GEN p1, p2, b, a = gcoeff(x,i,j) = remii(gcoeff(x,i,j), d);
    1229    57282724 :       if (!signe(a)) continue;
    1230             : 
    1231    30434105 :       k = j+1;
    1232    30434105 :       b = gcoeff(x,i,k) = remii(gcoeff(x,i,k), d);
    1233    30434165 :       if (!signe(b)) { swap(gel(x,j), gel(x,k)); continue; }
    1234    17528312 :       if (add_N)
    1235             :       { /* ensure the moving pivot on row i divides d from now on */
    1236     5436328 :         add_N = 0;
    1237     5436328 :         if (!equali1(a))
    1238             :         { /* x[j] *= u; after this, a = x[i,j] | d */
    1239     4503320 :           GEN u = Fp_invgen(a, d, &a);
    1240             :           long t;
    1241     4503320 :           p1 = gel(x,j);
    1242    15379056 :           for (t = 1; t < i; t++) gel(p1,t) = mulii(gel(p1,t), u);
    1243     4503319 :           FpV_red_part_ipvec(p1, D, i-1);
    1244     4503320 :           gel(p1,i) = a;
    1245     4503320 :           if (2*lg(a) < lg(b))
    1246             :           { /* reduce x[i,k] mod x[i,j]: helps ZC_elem */
    1247        9986 :             GEN r, q = dvmdii(b, a, &r);
    1248        9986 :             togglesign(q);
    1249        9986 :             ZC_lincomb1_inplace_i(gel(x,k), gel(x,j), q, i-1);
    1250        9986 :             FpV_red_part_ipvec(gel(x,k), D, i-1);
    1251        9986 :             gcoeff(x,i,k) = b = r;
    1252             :           }
    1253             :         }
    1254             :       }
    1255    17528312 :       ZC_elem(a,b, x, NULL, j,k);
    1256    17528290 :       p1 = gel(x,j);
    1257    17528290 :       p2 = gel(x,k);
    1258             :       /* prevent coeffs explosion */
    1259    83094976 :       for (k = 1; k < i; k++)
    1260             :       {
    1261    65566686 :         if (lgefint(gel(p1,k)) > LDM[k]) gel(p1,k) = remii(gel(p1,k),gel(D,k));
    1262    65566686 :         if (lgefint(gel(p2,k)) > LDM[k]) gel(p2,k) = remii(gel(p2,k),gel(D,k));
    1263             :       }
    1264             :     }
    1265     9985348 :     if (gc_needed(av,2))
    1266             :     {
    1267           7 :       if (DEBUGMEM>1) pari_warn(warnmem,"ZM_hnfmod[1]. i=%ld",i);
    1268           7 :       x = gerepilecopy(av, x);
    1269             :     }
    1270     9985348 :     if (moddiag && !signe(gcoeff(x,i,def)))
    1271             :     { /* missing pivot on line i, insert column */
    1272      425673 :       GEN a = cgetg(co + 1, t_MAT);
    1273     2169925 :       for (k = 1; k <= def; k++) gel(a,k) = gel(x,k);
    1274      425673 :       gel(a,k++) = Rg_col_ei(gel(D,i), nli, i);
    1275     2106240 :       for (     ; k <= co;  k++) gel(a,k) = gel(x,k-1);
    1276      425673 :       ldef--; if (ldef < 0) ldef = 0;
    1277      425673 :       co++; def++; x = a;
    1278             :     }
    1279             :   }
    1280     3477462 :   if (co < li)
    1281             :   { /* implies moddiag, add missing diag(D) components */
    1282       29015 :     GEN a = cgetg(li+1, t_MAT);
    1283       92694 :     for (k = 1; k <= li-co; k++) gel(a,k) = Rg_col_ei(gel(D,k), nli, k);
    1284      112364 :     for (i = 1; i < co; i++) gel(a,k-1+i) = gel(x,i);
    1285       29015 :     gel(a,li) = zerocol(nli); x = a;
    1286             :   }
    1287             :   else
    1288             :   {
    1289     3448447 :     x += co - li;
    1290     3448447 :     x[0] = evaltyp(t_MAT) | evallg(li); /* kill 0 columns */
    1291     3448468 :     if (moddiag) x = shallowconcat(x, zerocol(nli));
    1292             :   }
    1293     3477497 :   if (moddiag)
    1294             :   { /* x[li]: extra column, an accumulator discarded at the end */
    1295             :     GEN D2;
    1296     3418319 :     gcoeff(x,1,1) = gcdii(gcoeff(x,1,1), gel(D,1));
    1297     3418295 :     D2 = optimal_D(x,D);
    1298             :     /* add up missing diag(D) components */
    1299    13329683 :     for (i = nli; i > 0; i--)
    1300             :     {
    1301     9911365 :       gcoeff(x, i, li) = gel(D,i);
    1302    38750342 :       for (j = i; j > 0; j--)
    1303             :       {
    1304    28839021 :         GEN a = gcoeff(x, j, li);
    1305    28839021 :         if (!signe(a)) continue;
    1306    10239246 :         ZC_elem(a, gcoeff(x,j,j), x, NULL, li,j);
    1307    10239184 :         FpV_red_part_ipvec(gel(x,li), D, j-1);
    1308    10239137 :         FpV_red_part_ipvec(gel(x,j),  D, j-1);
    1309             :       }
    1310     9911321 :       if (gc_needed(av,1))
    1311             :       {
    1312           0 :         if (DEBUGMEM>1) pari_warn(warnmem,"ZM_hnfmod[2]. i=%ld", i);
    1313           0 :         gerepileall(av, 2, &x, &D2);
    1314             :       }
    1315             :     }
    1316     3418318 :     D = D2;
    1317             :   }
    1318             :   else
    1319             :   {
    1320       59178 :     GEN b = gel(D,1);
    1321      196861 :     for (i = nli; i > 0; i--)
    1322             :     {
    1323      137683 :       GEN d = bezout(gcoeff(x,i,i),b, &u,NULL);
    1324      137683 :       gcoeff(x,i,i) = d;
    1325      137683 :       FpV_Fp_mul_part_ip(gel(x,i), u, b, i-1);
    1326      137683 :       if (i > 1) b = diviiexact(b,d);
    1327             :     }
    1328       59178 :     D = optimal_D(x,D);
    1329             :   }
    1330     3477496 :   x[0] = evaltyp(t_MAT) | evallg(li); /* kill 0 columns / discard accumulator */
    1331     3477492 :   if (flag & hnf_PART) return x;
    1332             : 
    1333    13526482 :   for (i = nli; i > 0; i--)
    1334             :   {
    1335    10049056 :     GEN diag = gcoeff(x,i,i);
    1336    10049056 :     if (signe(diag) < 0) { gel(x,i) = ZC_neg(gel(x,i)); diag = gcoeff(x,i,i); }
    1337    10049056 :     if (i != nli)
    1338    19098507 :       for (j = 1; j < i; j++) gcoeff(x,j,i) = remii(gcoeff(x,j,i), gel(D,j));
    1339    29147493 :     for (j = i+1; j < li; j++)
    1340             :     {
    1341    19098483 :       GEN b = gcoeff(x,i,j) = remii(gcoeff(x,i,j), gel(D,i));
    1342    19098475 :       GEN r, q = truedvmdii(b, diag, &r);
    1343             :       /* ensure -diag/2 <= r < diag/2 */
    1344    19098442 :       if (center && signe(r) && abscmpii(shifti(r,1),diag) >= 0)
    1345      453724 :       { r = subii(r,diag); q = addiu(q,1); }
    1346    19098442 :       if (!signe(q)) continue;
    1347     8254981 :       togglesign(q);
    1348     8254981 :       ZC_lincomb1_inplace_i(gel(x,j), gel(x,i), q, i-1);
    1349     8254975 :       gcoeff(x,i,j) = r;
    1350             :     }
    1351    10049010 :     if (gc_needed(av,1))
    1352             :     {
    1353          14 :       if (DEBUGMEM>1) pari_warn(warnmem,"ZM_hnfmod[3]. i=%ld", i);
    1354          14 :       gerepileall(av, 2, &x, &D);
    1355             :     }
    1356             :   }
    1357     3477426 :   return x;
    1358             : }
    1359             : GEN
    1360     3467603 : ZM_hnfmodall(GEN x, GEN dm, long flag)
    1361             : {
    1362     3467603 :   pari_sp av = avma;
    1363     3467603 :   return gerepilecopy(av, ZM_hnfmodall_i(x, dm, flag));
    1364             : }
    1365             : GEN
    1366       59171 : ZM_hnfmod(GEN x, GEN d) { return ZM_hnfmodall(x,d,0); }
    1367             : GEN
    1368     3172435 : ZM_hnfmodid(GEN x, GEN d)
    1369     3172435 : { return ZM_hnfmodall(x,d,hnf_MODID); }
    1370             : /* return the column echelon form of x with 1's as pivots,
    1371             :  * P contains the row indices containing the pivots in increasing order */
    1372             : static GEN
    1373     2129346 : FpM_echelon(GEN x, GEN *pP, GEN p)
    1374             : {
    1375             :   pari_sp av;
    1376             :   long iP, li, co, i, j, k, def, ldef;
    1377             :   GEN P;
    1378             : 
    1379     2129346 :   co = lg(x); if (co == 1) { *pP = cgetg(1,t_VECSMALL); return cgetg(1,t_MAT); }
    1380     2129346 :   li = lgcols(x);
    1381     2129340 :   iP = 1;
    1382     2129340 :   *pP = P = cgetg(li, t_VECSMALL);
    1383     2129340 :   av = avma;
    1384     2129340 :   x = FpM_red(x, p);
    1385             : 
    1386     2129309 :   ldef = (li > co)? li - co: 0;
    1387     8028000 :   for (def = co-1,i = li-1; i > ldef; i--)
    1388             :   {
    1389     5898412 :     GEN u = NULL;
    1390     8694178 :     for (k = def; k; k--)
    1391             :     {
    1392     6467005 :       u = gcoeff(x,i,k);
    1393     6467005 :       if (signe(u)) break;
    1394             :     }
    1395     5898412 :     if (!k)
    1396             :     {
    1397     2227418 :       if (--ldef < 0) ldef = 0;
    1398     2227418 :       continue;
    1399             :     }
    1400     3670994 :     P[iP++] = i;
    1401     3670994 :     if (k != def) swap(gel(x,def), gel(x,k));
    1402     3670994 :     if (!equali1(u))
    1403     3387596 :       FpV_Fp_mul_part_ip(gel(x,def), Fp_inv(u,p), p, i-1);
    1404     3671108 :     gcoeff(x, i, def) = gen_1;
    1405    12582857 :     for (j = def-1; j; j--)
    1406             :     { /* zero x[i, 1..def-1] using x[i,def] = 1*/
    1407     8911624 :       GEN xj = gel(x,j), u = gel(xj,i);
    1408     8911624 :       if (!signe(u)) continue;
    1409             : 
    1410     5188199 :       ZC_lincomb1_inplace(xj, gel(x,def), negi(u));
    1411    27080324 :       for (k = 1; k < i; k++) gel(xj,k) = modii(gel(xj,k), p);
    1412             :     }
    1413     3671233 :     if (gc_needed(av,2))
    1414             :     {
    1415           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"FpM_echelon. i=%ld",i);
    1416           0 :       x = gerepilecopy(av, x);
    1417             :     }
    1418     3671273 :     def--;
    1419             :   }
    1420             :   /* rank = iP - 1 */
    1421     2129588 :   setlg(P, iP); vecsmall_sort(P);
    1422     2129339 :   if (co > iP) x += co - iP;
    1423     2129339 :   x[0] = evaltyp(t_MAT) | evallg(iP);
    1424     2129336 :   return x;
    1425             : }
    1426             : /* given x square of maximal rank with 1 or p on diagonal from hnfmodid
    1427             :  * (=> a column containing p has its other entries at 0 ), return the HNF */
    1428             : static GEN
    1429     2129375 : FpM_hnfend(pari_sp av, GEN x, GEN p)
    1430             : {
    1431     2129375 :   long i, l = lgcols(x);
    1432     8027649 :   for (i = l-1; i > 0; i--)
    1433             :   {
    1434     5898648 :     GEN diag = gcoeff(x,i,i);
    1435             :     long j;
    1436     5898648 :     if (is_pm1(diag))
    1437     8550600 :       for (j = i+1; j < l; j++)
    1438             :       {
    1439     4879272 :         GEN xj = gel(x,j), b = gel(xj,i);
    1440             :         long k;
    1441     4879272 :         if (!signe(b)) continue;
    1442     3603628 :         ZC_lincomb1_inplace(xj, gel(x,i), negi(b));
    1443    15244461 :         for (k=1; k<i; k++)
    1444    11640923 :           if (lgefint(gel(xj,k)) > 3) gel(xj,k) = remii(gel(xj,k), p);
    1445             :       }
    1446             :     else
    1447     6423138 :       for (j = i+1; j < l; j++) gcoeff(x,i,j) = modii(gcoeff(x,i,j), p);
    1448     5898503 :     if (gc_needed(av,2))
    1449             :     {
    1450           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"FpM_hnfend. i=%ld",i);
    1451           0 :       x = gerepilecopy(av, x);
    1452             :     }
    1453             :   }
    1454     2129001 :   return x;
    1455             : }
    1456             : GEN
    1457     2129345 : ZM_hnfmodprime(GEN x, GEN p)
    1458             : {
    1459     2129345 :   pari_sp av = avma;
    1460             :   GEN P, y;
    1461             :   long l, lP, i;
    1462     2129345 :   if (lg(x) == 1) return cgetg(1, t_MAT);
    1463     2129345 :   l = lgcols(x);
    1464     2129345 :   x = FpM_echelon(x, &P, p);
    1465     2129334 :   lP = lg(P); /* rank = lP-1 */
    1466     2129334 :   if (lP == l) { set_avma(av); return matid(l-1); }
    1467     2129334 :   y = scalarmat_shallow(p, l-1);
    1468     5800779 :   for (i = 1; i < lP; i++) gel(y,P[i]) = gel(x,i);
    1469     2129359 :   return gerepilecopy(av, FpM_hnfend(av,y,p));
    1470             : }
    1471             : 
    1472             : static GEN
    1473      234936 : allhnfmod(GEN x, GEN dm, int flag)
    1474             : {
    1475      234936 :   if (typ(x)!=t_MAT) pari_err_TYPE("allhnfmod",x);
    1476      234936 :   RgM_check_ZM(x, "allhnfmod");
    1477      234929 :   if (isintzero(dm)) return ZM_hnf(x);
    1478      234929 :   return ZM_hnfmodall(x, dm, flag);
    1479             : }
    1480             : GEN
    1481          14 : hnfmod(GEN x, GEN d)
    1482             : {
    1483          14 :   if (typ(d) != t_INT) pari_err_TYPE("mathnfmod",d);
    1484          14 :   return allhnfmod(x, d, 0);
    1485             : }
    1486             : GEN
    1487      234922 : hnfmodid(GEN x, GEN d)
    1488             : {
    1489      234922 :   switch(typ(d))
    1490             :   {
    1491      233984 :     case t_INT: break;
    1492         938 :     case t_VEC: case t_COL:
    1493         938 :       if (RgV_is_ZV(d)) break;
    1494           0 :     default: pari_err_TYPE("mathnfmodid",d);
    1495             :   }
    1496      234922 :   return allhnfmod(x, d, hnf_MODID);
    1497             : }
    1498             : 
    1499             : /* M a ZM in HNF. Normalize with *centered* residues */
    1500             : GEN
    1501        2246 : ZM_hnfcenter(GEN M)
    1502             : {
    1503        2246 :   long i, j, k, N = lg(M)-1;
    1504        2246 :   pari_sp av = avma;
    1505             : 
    1506        8689 :   for (j=N-1; j>0; j--) /* skip last line */
    1507             :   {
    1508        6443 :     GEN Mj = gel(M,j), a = gel(Mj,j);
    1509       33691 :     for (k = j+1; k <= N; k++)
    1510             :     {
    1511       27248 :       GEN Mk = gel(M,k), q = diviiround(gel(Mk,j), a);
    1512       27248 :       long s = signe(q);
    1513       27248 :       if (!s) continue;
    1514       17536 :       if (is_pm1(q))
    1515             :       {
    1516        3809 :         if (s < 0)
    1517        1965 :           for (i = 1; i <= j; i++) gel(Mk,i) = addii(gel(Mk,i), gel(Mj,i));
    1518             :         else
    1519       12101 :           for (i = 1; i <= j; i++) gel(Mk,i) = subii(gel(Mk,i), gel(Mj,i));
    1520             :       }
    1521             :       else
    1522      123277 :         for (i = 1; i <= j; i++) gel(Mk,i) = subii(gel(Mk,i), mulii(q,gel(Mj,i)));
    1523       17536 :       if (gc_needed(av,1))
    1524             :       {
    1525           0 :         if (DEBUGMEM>1) pari_warn(warnmem,"ZM_hnfcenter, j = %ld",j);
    1526           0 :         M = gerepilecopy(av, M);
    1527             :       }
    1528             :     }
    1529             :   }
    1530        2246 :   return M;
    1531             : }
    1532             : 
    1533             : /***********************************************************************/
    1534             : /*                                                                     */
    1535             : /*                 HNFLLL (Havas, Majewski, Mathews)                   */
    1536             : /*                                                                     */
    1537             : /***********************************************************************/
    1538             : 
    1539             : static void
    1540      606807 : Minus(long j, GEN lambda)
    1541             : {
    1542      606807 :   long k, n = lg(lambda);
    1543             : 
    1544     4096015 :   for (k=1  ; k<j; k++) togglesign_safe(&gcoeff(lambda,k,j));
    1545     6298051 :   for (k=j+1; k<n; k++) togglesign_safe(&gcoeff(lambda,j,k));
    1546      606807 : }
    1547             : 
    1548             : /* index of first non-zero entry */
    1549             : static long
    1550    36365169 : findi(GEN M)
    1551             : {
    1552    36365169 :   long i, n = lg(M);
    1553   660367577 :   for (i=1; i<n; i++)
    1554   651695078 :     if (signe(gel(M,i))) return i;
    1555     8672499 :   return 0;
    1556             : }
    1557             : 
    1558             : static long
    1559    36365169 : findi_normalize(GEN Aj, GEN B, long j, GEN lambda)
    1560             : {
    1561    36365169 :   long r = findi(Aj);
    1562    36365169 :   if (r && signe(gel(Aj,r)) < 0)
    1563             :   {
    1564      606807 :     ZV_togglesign(Aj); if (B) ZV_togglesign(gel(B,j));
    1565      606807 :     Minus(j,lambda);
    1566             :   }
    1567    36365169 :   return r;
    1568             : }
    1569             : 
    1570             : static void
    1571    18182294 : reduce2(GEN A, GEN B, long k, long j, long *row0, long *row1, GEN lambda, GEN D)
    1572             : {
    1573             :   GEN q;
    1574             :   long i;
    1575             : 
    1576    18182294 :   *row0 = findi_normalize(gel(A,j), B,j,lambda);
    1577    18182294 :   *row1 = findi_normalize(gel(A,k), B,k,lambda);
    1578    18182294 :   if (*row0)
    1579    11883994 :     q = truedivii(gcoeff(A,*row0,k), gcoeff(A,*row0,j));
    1580     6298300 :   else if (abscmpii(shifti(gcoeff(lambda,j,k), 1), gel(D,j)) > 0)
    1581     2797205 :     q = diviiround(gcoeff(lambda,j,k), gel(D,j));
    1582             :   else
    1583     3501095 :     return;
    1584             : 
    1585    14681199 :   if (signe(q))
    1586             :   {
    1587     6257000 :     GEN Lk = gel(lambda,k), Lj = gel(lambda,j);
    1588     6257000 :     togglesign_safe(&q);
    1589     6257000 :     if (*row0) ZC_lincomb1_inplace(gel(A,k),gel(A,j),q);
    1590     6257000 :     if (B) ZC_lincomb1_inplace(gel(B,k),gel(B,j),q);
    1591     6257000 :     gel(Lk,j) = addmulii(gel(Lk,j), q, gel(D,j));
    1592     6257000 :     if (is_pm1(q))
    1593             :     {
    1594     1904783 :       if (signe(q) > 0)
    1595             :       {
    1596     2384597 :         for (i=1; i<j; i++)
    1597     1718561 :           if (signe(gel(Lj,i))) gel(Lk,i) = addii(gel(Lk,i), gel(Lj,i));
    1598             :       }
    1599             :       else
    1600             :       {
    1601     4226942 :         for (i=1; i<j; i++)
    1602     2988195 :           if (signe(gel(Lj,i))) gel(Lk,i) = subii(gel(Lk,i), gel(Lj,i));
    1603             :       }
    1604             :     }
    1605             :     else
    1606             :     {
    1607    18005721 :       for (i=1; i<j; i++)
    1608    13653504 :         if (signe(gel(Lj,i))) gel(Lk,i) = addmulii(gel(Lk,i), q, gel(Lj,i));
    1609             :     }
    1610             :   }
    1611             : }
    1612             : 
    1613             : static void
    1614     3236908 : hnfswap(GEN A, GEN B, long k, GEN lambda, GEN D)
    1615             : {
    1616     3236908 :   GEN t, p1, p2, Lk = gel(lambda,k);
    1617     3236908 :   long i,j,n = lg(A);
    1618             : 
    1619     3236908 :   swap(gel(A,k), gel(A,k-1));
    1620     3236908 :   if (B) swap(gel(B,k), gel(B,k-1));
    1621    16218263 :   for (j=k-2; j; j--) swap(gcoeff(lambda,j,k-1), gel(Lk,j));
    1622    41172260 :   for (i=k+1; i<n; i++)
    1623             :   {
    1624    37935352 :     GEN Li = gel(lambda,i);
    1625    37935352 :     p1 = mulii(gel(Li,k-1), gel(D,k));
    1626    37935352 :     p2 = mulii(gel(Li,k), gel(Lk,k-1));
    1627    37935352 :     t = subii(p1,p2);
    1628             : 
    1629    37935352 :     p1 = mulii(gel(Li,k), gel(D,k-2));
    1630    37935352 :     p2 = mulii(gel(Li,k-1), gel(Lk,k-1));
    1631    37935352 :     gel(Li,k-1) = diviiexact(addii(p1,p2), gel(D,k-1));
    1632    37935352 :     gel(Li,k) = diviiexact(t, gel(D,k-1));
    1633             :   }
    1634     3236908 :   p1 = mulii(gel(D,k-2), gel(D,k));
    1635     3236908 :   p2 = sqri(gel(Lk,k-1));
    1636     3236908 :   gel(D,k-1) = diviiexact(addii(p1,p2), gel(D,k-1));
    1637     3236908 : }
    1638             : 
    1639             : /* reverse row order in matrix A, IN PLACE */
    1640             : static GEN
    1641      113418 : reverse_rows(GEN A)
    1642             : {
    1643      113418 :   long i, j, h, n = lg(A);
    1644      113418 :   if (n == 1) return A;
    1645      113418 :   h = lgcols(A);
    1646      965614 :   for (j=1; j<n; j++)
    1647             :   {
    1648      852196 :     GEN c = gel(A,j);
    1649             :     /* start at (h-1) >>1 : if h = 2i even, no need to swap c[i] and itself */
    1650     3320602 :     for (i=(h-1)>>1; i; i--) swap(gel(c,i), gel(c,h-i));
    1651             :   }
    1652      113418 :   return A;
    1653             : }
    1654             : /* decide whether to swap */
    1655             : static int
    1656     1246817 : must_swap(long k, GEN lambda, GEN D)
    1657             : {
    1658     1246817 :   pari_sp av = avma;
    1659     1246817 :   GEN z = addii(mulii(gel(D,k-2),gel(D,k)), sqri(gcoeff(lambda,k-1,k)));
    1660     1246817 :   return gc_bool(av, cmpii(z, sqri(gel(D,k-1))) < 0);
    1661             : }
    1662             : 
    1663             : GEN
    1664       56709 : ZM_hnflll(GEN A, GEN *ptB, int remove)
    1665             : {
    1666       56709 :   pari_sp av = avma;
    1667             :   long n, k, kmax;
    1668             :   GEN B, lambda, D;
    1669             : 
    1670       56709 :   n = lg(A);
    1671       56709 :   A = reverse_rows(ZM_copy(A)); /* ZM_copy for in place findi_normalize() */
    1672       56709 :   B = ptB? matid(n-1): NULL;
    1673       56709 :   D = const_vec(n, gen_1) + 1;
    1674       56709 :   lambda = zeromatcopy(n-1,n-1);
    1675       56709 :   k = kmax = 2;
    1676     5413187 :   while (k < n)
    1677             :   {
    1678             :     long row0, row1;
    1679             :     int do_swap;
    1680     5356478 :     reduce2(A,B,k,k-1,&row0,&row1,lambda,D);
    1681     5356478 :     if      (row0) do_swap = (!row1 || row0 <= row1);
    1682     1443957 :     else if (row1) do_swap = 0;
    1683     1175259 :     else do_swap = must_swap(k,lambda,D);
    1684     5356478 :     if (do_swap)
    1685             :     {
    1686     3135273 :       hnfswap(A,B,k,lambda,D);
    1687     3135273 :       if (k > 2) k--;
    1688             :     }
    1689             :     else
    1690             :     {
    1691             :       long i;
    1692    15047021 :       for (i=k-2; i; i--)
    1693             :       {
    1694             :         long row0, row1;
    1695    12825816 :         reduce2(A,B,k,i,&row0,&row1,lambda,D);
    1696             :       }
    1697     2221205 :       if (++k > kmax) kmax = k;
    1698             :     }
    1699     5356478 :     if (gc_needed(av,3))
    1700             :     {
    1701         159 :       GEN b = D-1;
    1702         159 :       if (DEBUGMEM>1) pari_warn(warnmem,"hnflll, kmax = %ld / %ld",kmax,n-1);
    1703         159 :       gerepileall(av, B? 4: 3, &A, &lambda, &b, &B);
    1704         159 :       if (gc_needed(av,1)) paristack_resize(0); /* avoid desperation GC */
    1705         159 :       D = b+1;
    1706             :     }
    1707             :   }
    1708             :   /* handle trivial case: return negative diag coefficient otherwise */
    1709       56709 :   if (n == 2) (void)findi_normalize(gel(A,1), B,1,lambda);
    1710       56709 :   A = reverse_rows(A);
    1711       56709 :   if (remove)
    1712             :   {
    1713             :     long i;
    1714        3892 :     for (i = 1; i < n; i++)
    1715        3780 :       if (!ZV_equal0(gel(A,i))) break;
    1716        1407 :     remove_0cols(i-1, &A, &B, remove);
    1717             :   }
    1718       56709 :   gerepileall(av, B? 2: 1, &A, &B);
    1719       56709 :   if (B) *ptB = B;
    1720       56709 :   return A;
    1721             : }
    1722             : 
    1723             : GEN
    1724           7 : hnflll(GEN x)
    1725             : {
    1726           7 :   GEN z = cgetg(3, t_VEC);
    1727           7 :   gel(z,1) = ZM_hnflll(x, &gel(z,2), 1);
    1728           7 :   return z;
    1729             : }
    1730             : 
    1731             : /* Variation on HNFLLL: Extended GCD */
    1732             : 
    1733             : static void
    1734      472518 : reduce1(GEN A, GEN B, long k, long j, GEN lambda, GEN D)
    1735             : {
    1736             :   GEN q;
    1737             :   long i;
    1738             : 
    1739      472518 :   if (signe(gel(A,j)))
    1740       94270 :     q = diviiround(gel(A,k),gel(A,j));
    1741      378248 :   else if (abscmpii(shifti(gcoeff(lambda,j,k), 1), gel(D,j)) > 0)
    1742       14242 :     q = diviiround(gcoeff(lambda,j,k), gel(D,j));
    1743             :   else
    1744      364006 :     return;
    1745             : 
    1746      108512 :   if (signe(q))
    1747             :   {
    1748       91004 :     GEN Lk = gel(lambda,k), Lj = gel(lambda,j);
    1749       91004 :     togglesign_safe(&q);
    1750       91004 :     gel(A,k) = addmulii(gel(A,k), q, gel(A,j));
    1751       91004 :     ZC_lincomb1_inplace(gel(B,k),gel(B,j),q);
    1752       91004 :     gel(Lk,j) = addmulii(gel(Lk,j), q, gel(D,j));
    1753      123973 :     for (i=1; i<j; i++)
    1754       32969 :       if (signe(gel(Lj,i))) gel(Lk,i) = addmulii(gel(Lk,i), q, gel(Lj,i));
    1755             :   }
    1756             : }
    1757             : 
    1758             : static GEN
    1759       88316 : ZV_gcdext_i(GEN A)
    1760             : {
    1761       88316 :   long k, n = lg(A);
    1762             :   GEN B, lambda, D;
    1763             : 
    1764       88316 :   if (n == 1) retmkvec2(gen_1, cgetg(1,t_MAT));
    1765       88316 :   A = leafcopy(A);
    1766       88316 :   B = matid(n-1);
    1767       88316 :   lambda = zeromatcopy(n-1,n-1);
    1768       88316 :   D = const_vec(n, gen_1) + 1;
    1769       88316 :   k = 2;
    1770      372953 :   while (k < n)
    1771             :   {
    1772             :     int do_swap;
    1773             : 
    1774      284637 :     reduce1(A,B,k,k-1,lambda,D);
    1775      284637 :     if    (signe(gel(A,k-1))) do_swap = 1;
    1776      190367 :     else if (signe(gel(A,k))) do_swap = 0;
    1777       71558 :     else do_swap = must_swap(k,lambda,D);
    1778      284637 :     if (do_swap)
    1779             :     {
    1780      101635 :       hnfswap(A,B,k,lambda,D);
    1781      101635 :       if (k > 2) k--;
    1782             :     }
    1783             :     else
    1784             :     {
    1785             :       long i;
    1786      370883 :       for (i=k-2; i; i--) reduce1(A,B,k,i,lambda,D);
    1787      183002 :       k++;
    1788             :     }
    1789             :   }
    1790       88316 :   if (signe(gel(A,n-1)) < 0)
    1791             :   {
    1792        8048 :     gel(A,n-1) = negi(gel(A,n-1));
    1793        8048 :     ZV_togglesign(gel(B,n-1));
    1794             :   }
    1795       88316 :   return mkvec2(gel(A,n-1), B);
    1796             : }
    1797             : GEN
    1798       88302 : ZV_extgcd(GEN A)
    1799             : {
    1800       88302 :   pari_sp av = avma;
    1801       88302 :   return gerepilecopy(av, ZV_gcdext_i(A));
    1802             : }
    1803             : /* as ZV_extgcd, transforming the gcd into a t_MAT, for mathnf0 */
    1804             : static GEN
    1805          21 : ZV_hnfgcdext(GEN A)
    1806             : {
    1807          21 :   pari_sp av = avma;
    1808             :   GEN z;
    1809          21 :   if (lg(A) == 1) retmkvec2(cgetg(1,t_MAT),cgetg(1,t_MAT));
    1810          14 :   z = ZV_gcdext_i(A);
    1811          14 :   gel(z,1) = mkmat(mkcol(gel(z,1)));
    1812          14 :   return gerepilecopy(av, z);
    1813             : }
    1814             : 
    1815             : /* HNF with permutation. */
    1816             : GEN
    1817         357 : ZM_hnfperm(GEN A, GEN *ptU, GEN *ptperm)
    1818             : {
    1819             :   GEN U, c, l, perm, d, p, q, b;
    1820         357 :   pari_sp av = avma, av1;
    1821             :   long r, t, i, j, j1, k, m, n;
    1822             : 
    1823         357 :   n = lg(A)-1;
    1824         357 :   if (!n)
    1825             :   {
    1826           7 :     if (ptU) *ptU = cgetg(1,t_MAT);
    1827           7 :     if (ptperm) *ptperm = cgetg(1,t_VEC);
    1828           7 :     return cgetg(1, t_MAT);
    1829             :   }
    1830         350 :   m = nbrows(A);
    1831         350 :   c = zero_zv(m);
    1832         350 :   l = zero_zv(n);
    1833         350 :   perm = cgetg(m+1, t_VECSMALL);
    1834         350 :   av1 = avma;
    1835         350 :   A = RgM_shallowcopy(A);
    1836         350 :   U = ptU? matid(n): NULL;
    1837             :   /* U base change matrix : A0*U = A all along */
    1838        1575 :   for (r=0, k=1; k <= n; k++)
    1839             :   {
    1840        3633 :     for (j=1; j<k; j++)
    1841             :     {
    1842        2408 :       if (!l[j]) continue;
    1843        2275 :       t=l[j]; b=gcoeff(A,t,k);
    1844        2275 :       if (!signe(b)) continue;
    1845             : 
    1846         553 :       ZC_elem(b,gcoeff(A,t,j), A,U,k,j);
    1847         553 :       d = gcoeff(A,t,j);
    1848         553 :       if (signe(d) < 0)
    1849             :       {
    1850           0 :         ZV_neg_inplace(gel(A,j));
    1851           0 :         if (U) ZV_togglesign(gel(U,j));
    1852           0 :         d = gcoeff(A,t,j);
    1853             :       }
    1854        1428 :       for (j1=1; j1<j; j1++)
    1855             :       {
    1856         875 :         if (!l[j1]) continue;
    1857         861 :         q = truedivii(gcoeff(A,t,j1),d);
    1858         861 :         if (!signe(q)) continue;
    1859             : 
    1860         329 :         togglesign(q);
    1861         329 :         ZC_lincomb1_inplace(gel(A,j1), gel(A,j), q);
    1862         329 :         if (U) ZC_lincomb1_inplace(gel(U,j1), gel(U,j), q);
    1863             :       }
    1864             :     }
    1865        4081 :     t = m; while (t && (c[t] || !signe(gcoeff(A,t,k)))) t--;
    1866        1225 :     if (t)
    1867             :     {
    1868        1120 :       p = gcoeff(A,t,k);
    1869       18704 :       for (i=t-1; i; i--)
    1870             :       {
    1871       17584 :         q = gcoeff(A,i,k);
    1872       17584 :         if (signe(q) && abscmpii(p,q) > 0) { p = q; t = i; }
    1873             :       }
    1874        1120 :       perm[++r] = l[k] = t; c[t] = k;
    1875        1120 :       if (signe(p) < 0)
    1876             :       {
    1877         154 :         ZV_neg_inplace(gel(A,k));
    1878         154 :         if (U) ZV_togglesign(gel(U,k));
    1879         154 :         p = gcoeff(A,t,k);
    1880             :       }
    1881             :       /* p > 0 */
    1882        3185 :       for (j=1; j<k; j++)
    1883             :       {
    1884        2065 :         if (!l[j]) continue;
    1885        2044 :         q = truedivii(gcoeff(A,t,j),p);
    1886        2044 :         if (!signe(q)) continue;
    1887             : 
    1888         231 :         togglesign(q);
    1889         231 :         ZC_lincomb1_inplace(gel(A,j), gel(A,k), q);
    1890         231 :         if (U) ZC_lincomb1_inplace(gel(U,j), gel(U,k), q);
    1891             :       }
    1892             :     }
    1893        1225 :     if (gc_needed(av1,1))
    1894             :     {
    1895           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"hnfperm, k=%ld",k);
    1896           0 :       gerepileall(av1, U? 2: 1, &A, &U);
    1897             :     }
    1898             :   }
    1899         350 :   if (r < m)
    1900             :   {
    1901        4802 :     for (i=1,k=r; i<=m; i++)
    1902        4508 :       if (!c[i]) perm[++k] = i;
    1903             :   }
    1904             : 
    1905             : /* We have A0*U=A, U in Gl(n,Z)
    1906             :  * basis for Im(A):  columns of A s.t l[j]>0 (r   cols)
    1907             :  * basis for Ker(A): columns of U s.t l[j]=0 (n-r cols) */
    1908         350 :   p = cgetg(r+1,t_MAT);
    1909        2625 :   for (i=1; i<=m/2; i++) lswap(perm[i], perm[m+1-i]);
    1910         350 :   if (U)
    1911             :   {
    1912          70 :     GEN u = cgetg(n+1,t_MAT);
    1913         322 :     for (t=1,k=r,j=1; j<=n; j++)
    1914         252 :       if (l[j])
    1915             :       {
    1916         154 :         u[k + n-r] = U[j];
    1917         154 :         gel(p,k--) = vecpermute(gel(A,j), perm);
    1918             :       }
    1919             :       else
    1920          98 :         u[t++] = U[j];
    1921          70 :     *ptU = u;
    1922          70 :     if (ptperm) *ptperm = perm;
    1923          70 :     gerepileall(av, ptperm? 3: 2, &p, ptU, ptperm);
    1924             :   }
    1925             :   else
    1926             :   {
    1927        1253 :     for (k=r,j=1; j<=n; j++)
    1928         973 :       if (l[j]) gel(p,k--) = vecpermute(gel(A,j), perm);
    1929         280 :     if (ptperm) *ptperm = perm;
    1930         280 :     gerepileall(av, ptperm? 2: 1, &p, ptperm);
    1931             :   }
    1932         350 :   return p;
    1933             : }
    1934             : 
    1935             : GEN
    1936         224 : ZM_hnf_knapsack(GEN x)
    1937             : {
    1938         224 :   GEN t, perm, H = ZM_hnfperm(x,NULL,&perm);
    1939         224 :   long i,j, l = lg(H), h = lgcols(H);
    1940        3059 :   for (i=1; i<h; i++)
    1941             :   {
    1942        2919 :     int fl = 0;
    1943       15841 :     for (j=1; j<l; j++)
    1944             :     {
    1945       13006 :       t = gcoeff(H,i,j);
    1946       13006 :       if (signe(t))
    1947             :       {
    1948        2954 :         if (!is_pm1(t) || fl) return NULL;
    1949        2870 :         fl = 1;
    1950             :       }
    1951             :     }
    1952             :   }
    1953         140 :   return rowpermute(H, perm_inv(perm));
    1954             : }
    1955             : 
    1956             : GEN
    1957          14 : hnfperm(GEN A)
    1958             : {
    1959          14 :   GEN y = cgetg(4, t_VEC);
    1960          14 :   gel(y,1) = ZM_hnfperm(A, &gel(y,2), &gel(y,3));
    1961          14 :   return y;
    1962             : }
    1963             : 
    1964             : /* Hermite Normal Form, with base change matrix if ptB != NULL.
    1965             :  * If 'remove' = 1, remove 0 columns (do NOT update *ptB accordingly)
    1966             :  * If 'remove' = 2, remove 0 columns and update *ptB accordingly */
    1967             : GEN
    1968      133645 : ZM_hnfall_i(GEN A, GEN *ptB, long remove)
    1969             : {
    1970             :   pari_sp av;
    1971             :   long m, n, r, i, j, k, li;
    1972             :   GEN B, c, h, a;
    1973             : 
    1974      133645 :   RgM_dimensions(A, &m,&n);
    1975      133645 :   if (!n)
    1976             :   {
    1977           7 :     if (ptB) *ptB = cgetg(1,t_MAT);
    1978           7 :     return cgetg(1,t_MAT);
    1979             :   }
    1980      133638 :   c = zero_zv(m);
    1981      133638 :   h = const_vecsmall(n, m);
    1982      133638 :   av = avma;
    1983      133638 :   A = RgM_shallowcopy(A);
    1984      133638 :   B = ptB? matid(n): NULL;
    1985      133638 :   r = n+1;
    1986      422364 :   for (li=m; li; li--)
    1987             :   {
    1988      912868 :     for (j=1; j<r; j++)
    1989             :     {
    1990     1063927 :       for (i=h[j]; i>li; i--)
    1991             :       {
    1992      151325 :         a = gcoeff(A,i,j);
    1993      151325 :         k = c[i];
    1994             :         /* zero a = Aij  using  Aik */
    1995      151325 :         if (signe(a)) ZC_elem(a,gcoeff(A,i,k), A,B,j,k);
    1996      151325 :         ZM_reduce(A,B, i,k); /* ensure reduced entries even if a = 0 */
    1997             :       }
    1998      912602 :       if (gc_needed(av,1) && (j & 0x7f) == 0)
    1999             :       {
    2000           0 :         if (DEBUGMEM>1)
    2001           0 :           pari_warn(warnmem,"ZM_hnfall[1], li = %ld, j = %ld", li, j);
    2002           0 :         gerepileall(av, B? 2: 1, &A, &B);
    2003             :       }
    2004      912602 :       if (signe( gcoeff(A,li,j) )) break;
    2005      624142 :       h[j] = li-1;
    2006             :     }
    2007      288726 :     if (j == r) continue;
    2008      288460 :     r--;
    2009      288460 :     if (j < r) /* A[j] != 0 */
    2010             :     {
    2011      213958 :       swap(gel(A,j), gel(A,r));
    2012      213958 :       if (B) swap(gel(B,j), gel(B,r));
    2013      213958 :       h[j] = h[r]; h[r] = li; c[li] = r;
    2014             :     }
    2015      288460 :     if (signe(gcoeff(A,li,r)) < 0)
    2016             :     {
    2017       58005 :       ZV_neg_inplace(gel(A,r));
    2018       58005 :       if (B) ZV_togglesign(gel(B,r));
    2019             :     }
    2020      288460 :     ZM_reduce(A,B, li,r);
    2021      288460 :     if (gc_needed(av,1))
    2022             :     {
    2023           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"ZM_hnfall[2], li = %ld", li);
    2024           0 :       gerepileall(av, B? 2: 1, &A, &B);
    2025             :     }
    2026             :   }
    2027             : 
    2028      133638 :   if (DEBUGLEVEL>5) err_printf("\nhnfall, final phase: ");
    2029      133638 :   r--; /* first r cols are in the image the n-r (independent) last ones */
    2030      785689 :   for (j=1; j<=r; j++)
    2031             :   {
    2032     1587173 :     for (i=h[j]; i; i--)
    2033             :     {
    2034      935122 :       a = gcoeff(A,i,j);
    2035      935122 :       k = c[i];
    2036      935122 :       if (signe(a)) ZC_elem(a,gcoeff(A,i,k), A,B, j,k);
    2037      935122 :       ZM_reduce(A,B, i,k); /* ensure reduced entries, even if a = 0 */
    2038             :     }
    2039      652051 :     if (gc_needed(av,1) && (j & 0x7f) == 0)
    2040             :     {
    2041           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"ZM_hnfall[3], j = %ld", j);
    2042           0 :       gerepileall(av, B? 2: 1, &A, &B);
    2043             :     }
    2044             :   }
    2045      133638 :   if (DEBUGLEVEL>5) err_printf("\n");
    2046      133638 :   if (remove) remove_0cols(r, &A, &B, remove);
    2047      133638 :   if (ptB) *ptB = B;
    2048      133638 :   return A;
    2049             : }
    2050             : GEN
    2051        4824 : ZM_hnfall(GEN A, GEN *ptB, long remove)
    2052             : {
    2053        4824 :   pari_sp av = avma;
    2054        4824 :   A = ZM_hnfall_i(A, ptB, remove);
    2055        4824 :   gerepileall(av, ptB? 2: 1, &A, ptB);
    2056        4824 :   return A;
    2057             : }
    2058             : 
    2059             : GEN
    2060          14 : hnfall(GEN x)
    2061             : {
    2062          14 :   GEN z = cgetg(3, t_VEC);
    2063          14 :   gel(z,1) = ZM_hnfall(x, (GEN*)(z+2), 1);
    2064          14 :   return z;
    2065             : }
    2066             : GEN
    2067          14 : hnf(GEN x) { return mathnf0(x,0); }
    2068             : 
    2069             : /* C = A^(-1)t where A and C are integral, A is upper triangular, t t_INT */
    2070             : GEN
    2071       22820 : hnf_invscale(GEN A, GEN t)
    2072             : {
    2073       22820 :   long n = lg(A)-1, i,j,k;
    2074       22820 :   GEN m, c = cgetg(n+1,t_MAT);
    2075             : 
    2076       22820 :   if (!n) return c;
    2077       85568 :   for (k=1; k<=n; k++)
    2078             :   { /* cf hnf_divscale with B = id, thus b = e_k */
    2079       62748 :     GEN u = cgetg(n+1, t_COL);
    2080       62748 :     pari_sp av = avma;
    2081       62748 :     gel(c,k) = u;
    2082       62748 :     gel(u,n) = k == n? gerepileuptoint(av, diviiexact(t, gcoeff(A,n,n))): gen_0;
    2083      273924 :     for (i=n-1; i>0; i--)
    2084             :     {
    2085      211176 :       av = avma; m = i == k? t: gen_0;
    2086     1080135 :       for (j=i+1; j<=n; j++) m = subii(m, mulii(gcoeff(A,i,j),gel(u,j)));
    2087      211176 :       gel(u,i) = gerepileuptoint(av, diviiexact(m, gcoeff(A,i,i)));
    2088             :     }
    2089             :   }
    2090       22820 :   return c;
    2091             : }
    2092             : 
    2093             : /* C = A^(-1)(tB) where A, B, C are integral, A is upper triangular, t t_INT */
    2094             : GEN
    2095      111242 : hnf_divscale(GEN A, GEN B, GEN t)
    2096             : {
    2097      111242 :   long n = lg(A)-1, i,j,k;
    2098      111242 :   GEN m, c = cgetg(n+1,t_MAT);
    2099             : 
    2100      111242 :   if (!n) return c;
    2101      450366 :   for (k=1; k<=n; k++)
    2102             :   {
    2103      339124 :     GEN u = cgetg(n+1, t_COL), b = gel(B,k);
    2104      339124 :     pari_sp av = avma;
    2105      339124 :     gel(c,k) = u; m = mulii(gel(b,n),t);
    2106      339124 :     gel(u,n) = gerepileuptoint(av, diviiexact(m, gcoeff(A,n,n)));
    2107     1644410 :     for (i=n-1; i>0; i--)
    2108             :     {
    2109     1305286 :       av = avma; m = mulii(gel(b,i),t);
    2110     9064057 :       for (j=i+1; j<=n; j++) m = subii(m, mulii(gcoeff(A,i,j),gel(u,j)));
    2111     1305285 :       gel(u,i) = gerepileuptoint(av, diviiexact(m, gcoeff(A,i,i)));
    2112             :     }
    2113             :   }
    2114      111242 :   return c;
    2115             : }
    2116             : 
    2117             : /* A, B integral upper HNF. A^(-1) B integral ? */
    2118             : int
    2119          91 : hnfdivide(GEN A, GEN B)
    2120             : {
    2121          91 :   pari_sp av = avma;
    2122          91 :   long n = lg(A)-1, i,j,k;
    2123             :   GEN u, b, m, r;
    2124             : 
    2125          91 :   if (!n) return 1;
    2126          91 :   if (lg(B)-1 != n) pari_err_DIM("hnfdivide");
    2127          91 :   u = cgetg(n+1, t_COL);
    2128         413 :   for (k=1; k<=n; k++)
    2129             :   {
    2130         322 :     b = gel(B,k);
    2131         322 :     m = gel(b,k);
    2132         322 :     gel(u,k) = dvmdii(m, gcoeff(A,k,k), &r);
    2133         322 :     if (r != gen_0) return gc_long(av, 0);
    2134         798 :     for (i=k-1; i>0; i--)
    2135             :     {
    2136         476 :       m = gel(b,i);
    2137        1337 :       for (j=i+1; j<=k; j++) m = subii(m, mulii(gcoeff(A,i,j),gel(u,j)));
    2138         476 :       m = dvmdii(m, gcoeff(A,i,i), &r);
    2139         476 :       if (r != gen_0) return gc_long(av, 0);
    2140         476 :       gel(u,i) = m;
    2141             :     }
    2142             :   }
    2143          91 :   return gc_long(av, 1);
    2144             : }
    2145             : 
    2146             : /* A upper HNF, b integral vector. Return A^(-1) b if integral,
    2147             :  * NULL otherwise. Assume #A[,1] = #b. */
    2148             : GEN
    2149      114689 : hnf_invimage(GEN A, GEN b)
    2150             : {
    2151      114689 :   pari_sp av = avma;
    2152      114689 :   long n = lg(A)-1, m, i, k;
    2153             :   GEN u, r;
    2154             : 
    2155      114689 :   if (!n) return lg(b)==1? cgetg(1,t_COL):NULL;
    2156      114647 :   m = nbrows(A); /* m >= n */
    2157      114647 :   u = cgetg(n+1, t_COL);
    2158      324124 :   for (i = n, k = m; k > 0; k--)
    2159             :   {
    2160      324124 :     pari_sp av2 = avma;
    2161             :     long j;
    2162      324124 :     GEN t = gel(b,k), Aki = gcoeff(A,k,i);
    2163      324124 :     if (typ(t) != t_INT) pari_err_TYPE("hnf_invimage",t);
    2164     1114280 :     for (j=i+1; j<=n; j++) t = subii(t, mulii(gcoeff(A,k,j),gel(u,j)));
    2165      324124 :     if (!signe(Aki))
    2166             :     {
    2167           0 :       if (signe(t)) return gc_NULL(av);
    2168           0 :       set_avma(av2); gel(u,i) = gen_0; continue;
    2169             :     }
    2170      324124 :     t = dvmdii(t, Aki, &r);
    2171      324124 :     if (r != gen_0) return gc_NULL(av);
    2172      270426 :     gel(u,i) = gerepileuptoint(av2, t);
    2173      270426 :     if (--i == 0) break;
    2174             :   }
    2175             :   /* If there is a solution, it must be u. Check remaining equations */
    2176      121898 :   for (; k > 0; k--)
    2177             :   {
    2178       60949 :     pari_sp av2 = avma;
    2179             :     long j;
    2180       60949 :     GEN t = gel(b,k);
    2181       60949 :     if (typ(t) != t_INT) pari_err_TYPE("hnf_invimage",t);
    2182      249865 :     for (j=1; j<=n; j++) t = subii(t, mulii(gcoeff(A,k,j),gel(u,j)));
    2183       60949 :     if (signe(t)) return gc_NULL(av);
    2184       60949 :     set_avma(av2);
    2185             :   }
    2186       60949 :   return u;
    2187             : }
    2188             : 
    2189             : /* A upper HNF, B integral matrix or column. Return A^(-1) B if integral,
    2190             :  * NULL otherwise */
    2191             : GEN
    2192      105651 : hnf_solve(GEN A, GEN B)
    2193             : {
    2194             :   pari_sp av;
    2195             :   long i, l;
    2196             :   GEN C;
    2197             : 
    2198      105651 :   if (typ(B) == t_COL) return hnf_invimage(A, B);
    2199       73976 :   av = avma; C = cgetg_copy(B, &l);
    2200      115304 :   for (i = 1; i < l; i++)
    2201             :   {
    2202       76342 :     GEN c = hnf_invimage(A, gel(B,i));
    2203       76342 :     if (!c) return gc_NULL(av);
    2204       41328 :     gel(C,i) = c;
    2205             :   }
    2206       38962 :   return C;
    2207             : }
    2208             : 
    2209             : /***************************************************************/
    2210             : /**                                                           **/
    2211             : /**               SMITH NORMAL FORM REDUCTION                 **/
    2212             : /**                                                           **/
    2213             : /***************************************************************/
    2214             : 
    2215             : static GEN
    2216           0 : trivsmith(long all)
    2217             : {
    2218             :   GEN z;
    2219           0 :   if (!all) return cgetg(1,t_VEC);
    2220           0 :   z=cgetg(4,t_VEC);
    2221           0 :   gel(z,1) = cgetg(1,t_MAT);
    2222           0 :   gel(z,2) = cgetg(1,t_MAT);
    2223           0 :   gel(z,3) = cgetg(1,t_MAT); return z;
    2224             : }
    2225             : 
    2226             : static void
    2227          38 : snf_pile1(pari_sp av, GEN *x, GEN *U)
    2228             : {
    2229             :   GEN *gptr[2];
    2230          38 :   int c = 1; gptr[0]=x;
    2231          38 :   if (*U) gptr[c++] = U;
    2232          38 :   gerepilemany(av,gptr,c);
    2233          38 : }
    2234             : static void
    2235      194122 : snf_pile(pari_sp av, GEN *x, GEN *U, GEN *V)
    2236             : {
    2237             :   GEN *gptr[3];
    2238      194122 :   int c = 1; gptr[0]=x;
    2239      194122 :   if (*U) gptr[c++] = U;
    2240      194122 :   if (*V) gptr[c++] = V;
    2241      194122 :   gerepilemany(av,gptr,c);
    2242      194122 : }
    2243             : 
    2244             : static GEN
    2245      138951 : bezout_step(GEN *pa, GEN *pb, GEN *pu, GEN *pv)
    2246             : {
    2247      138951 :   GEN a = *pa, b = *pb, d;
    2248      138951 :   if (absequalii(a,b))
    2249             :   {
    2250       51934 :     long sa = signe(a), sb = signe(b);
    2251       51934 :     *pv = gen_0;
    2252       51934 :     if (sb == sa) {
    2253       49381 :       *pa = *pb = gen_1;
    2254       49381 :       if (sa > 0) {*pu=gen_1; return a;} else {*pu=gen_m1; return absi(a);}
    2255             :     }
    2256        2553 :     if (sa > 0) { *pa = *pu = gen_1; *pb = gen_m1; return a; }
    2257          63 :     *pa = *pu = gen_m1; *pb = gen_1; return b;
    2258             :   }
    2259       87017 :   d = bezout(a,b, pu,pv);
    2260       87017 :   *pa = diviiexact(a, d);
    2261       87017 :   *pb = diviiexact(b, d); return d;
    2262             : }
    2263             : 
    2264             : static int
    2265      186515 : negcmpii(void *E, GEN x, GEN y) { (void)E; return -cmpii(x,y); }
    2266             : 
    2267             : /* x square of maximal rank; does b = x[i,i] divide all entries in
    2268             :  * x[1..i-1, 1..i-1] ? If so, return 0; else the index of a problematic row */
    2269             : static long
    2270      141108 : ZM_snf_no_divide(GEN x, long i)
    2271             : {
    2272      141108 :   GEN b = gcoeff(x,i,i);
    2273             :   long j, k;
    2274             : 
    2275      141108 :   if (is_pm1(b)) return 0;
    2276      148342 :   for (k = 1; k < i; k++)
    2277      760237 :     for (j = 1; j < i; j++)
    2278      668385 :       if (!dvdii(gcoeff(x,k,j),b)) return k;
    2279       40432 :   return 0;
    2280             : }
    2281             : 
    2282             : static void
    2283      183444 : ZM_redpart(GEN x, GEN p, long I)
    2284             : {
    2285      183444 :   long l = lgefint(p), i, j;
    2286      928804 :   for (i = 1; i <= I; i++)
    2287     7673060 :     for (j = 1; j <= I; j++)
    2288             :     {
    2289     6927700 :       GEN c = gcoeff(x,i,j);
    2290     6927700 :       if (lgefint(c) > l) gcoeff(x,i,j) = remii(c, p);
    2291             :     }
    2292      183444 : }
    2293             : static void
    2294      122603 : ZMrow_divexact_inplace(GEN M, long i, GEN c)
    2295             : {
    2296      122603 :   long j, l = lg(M);
    2297      533812 :   for (j = 1; j < l; j++) gcoeff(M,i,j) = diviiexact(gcoeff(M,i,j), c);
    2298      122603 : }
    2299             : 
    2300             : /* Return the SNF D of matrix X. If ptU/ptV non-NULL set them to U/V
    2301             :  * to that D = UXV */
    2302             : GEN
    2303       99721 : ZM_snfall_i(GEN x, GEN *ptU, GEN *ptV, long flag)
    2304             : {
    2305       99721 :   pari_sp av0 = avma, av;
    2306       99721 :   const long return_vec = flag & 1;
    2307             :   long i, j, k, m0, m, n0, n;
    2308       99721 :   GEN u, v, U, V, V0, mdet, A = NULL, perm = NULL;
    2309             : 
    2310       99721 :   n0 = n = lg(x)-1;
    2311       99721 :   if (!n) {
    2312        4726 :     if (ptU) *ptU = cgetg(1,t_MAT);
    2313        4726 :     if (ptV) *ptV = cgetg(1,t_MAT);
    2314        4726 :     return cgetg(1, return_vec? t_VEC: t_MAT);
    2315             :   }
    2316       94995 :   m0 = m = nbrows(x);
    2317             : 
    2318       94995 :   U = V = V0 = NULL; /* U = TRANSPOSE of row transform matrix [act on columns]*/
    2319       94995 :   if (m == n && ZM_ishnf(x))
    2320             :   {
    2321       90060 :     mdet = ZM_det_triangular(x); /* != 0 */
    2322       90060 :     if (ptV) *ptV = matid(n);
    2323             :   }
    2324             :   else
    2325             :   {
    2326        4935 :     mdet = ZM_detmult(x);
    2327        4935 :     if (!signe(mdet))
    2328          70 :       x = ZM_hnfperm(x, ptV, ptU? &perm: NULL);
    2329             :     else
    2330             :     { /* m <= n */
    2331        4865 :       if (!ptV)
    2332         336 :         x = ZM_hnfmod(x,mdet);
    2333        4529 :       else if (m == n)
    2334             :       {
    2335        4508 :         GEN H = ZM_hnfmod(x,mdet);
    2336        4508 :         *ptV = ZM_gauss(x,H);
    2337        4508 :         x = H;
    2338             :       }
    2339             :       else
    2340          21 :         x = ZM_hnfperm(x, ptV, ptU? &perm: NULL);
    2341        4865 :       mdet = ZM_det_triangular(x); /* > 0 */
    2342             :     }
    2343        4935 :     n = lg(x)-1; /* n independent columns */
    2344        4935 :     if (ptV)
    2345             :     {
    2346        4543 :       V = *ptV;
    2347        4543 :       if (n != n0)
    2348             :       {
    2349          21 :         V0 = vecslice(V, 1, n0 - n); /* kernel */
    2350          21 :         V  = vecslice(V, n0-n+1, n0);
    2351             :       }
    2352             :     }
    2353        4935 :     if (!signe(mdet))
    2354             :     {
    2355          70 :       if (n)
    2356             :       {
    2357          63 :         x = ZM_snfall_i(shallowtrans(x), ptV, ptU, return_vec); /* swap V,U */
    2358          63 :         if (!return_vec && n != m) x = shallowtrans(x);
    2359          63 :         if (ptV) V = ZM_mul(V, shallowtrans(*ptV));
    2360          63 :         if (ptU) U = *ptU; /* TRANSPOSE */
    2361             :       }
    2362             :       else /* 0 matrix */
    2363             :       {
    2364           7 :         x = cgetg(1,t_MAT);
    2365           7 :         if (ptV) V = cgetg(1, t_MAT);
    2366           7 :         if (ptU) U = matid(m);
    2367             :       }
    2368          70 :       goto THEEND;
    2369             :     }
    2370             :   }
    2371       94925 :   if (ptV || ptU) U = matid(n); /* we will compute V in terms of U */
    2372       94925 :   if (DEBUGLEVEL>7) err_printf("starting SNF loop");
    2373             : 
    2374             :   /* square, maximal rank n */
    2375       94925 :   A = x; x = shallowcopy(x); av = avma;
    2376      219975 :   for (i = n; i > 1; i--)
    2377             :   {
    2378      125050 :     if (DEBUGLEVEL>7) err_printf("\ni = %ld: ",i);
    2379             :     for(;;)
    2380       95300 :     {
    2381      220350 :       int c = 0;
    2382             :       GEN a, b;
    2383      795392 :       for (j = i-1; j >= 1; j--)
    2384             :       {
    2385      575042 :         b = gcoeff(x,i,j); if (!signe(b)) continue;
    2386       60341 :         a = gcoeff(x,i,i);
    2387       60341 :         ZC_elem(b, a, x,NULL, j,i);
    2388       60341 :         if (gc_needed(av,1))
    2389             :         {
    2390           1 :           if (DEBUGMEM>1) pari_warn(warnmem,"[1]: ZM_snfall i = %ld", i);
    2391           1 :           snf_pile1(av, &x,&U);
    2392             :         }
    2393             :       }
    2394      220350 :       if (DEBUGLEVEL>7) err_printf("; ");
    2395      795392 :       for (j=i-1; j>=1; j--)
    2396             :       {
    2397             :         GEN d;
    2398      575042 :         b = gcoeff(x,j,i); if (!signe(b)) continue;
    2399      138951 :         a = gcoeff(x,i,i);
    2400      138951 :         d = bezout_step(&a, &b, &u, &v);
    2401     2411060 :         for (k = 1; k < i; k++)
    2402             :         {
    2403     2272109 :           GEN t = addii(mulii(u,gcoeff(x,i,k)),mulii(v,gcoeff(x,j,k)));
    2404     2272109 :           gcoeff(x,j,k) = subii(mulii(a,gcoeff(x,j,k)),
    2405     2272109 :                                 mulii(b,gcoeff(x,i,k)));
    2406     2272109 :           gcoeff(x,i,k) = t;
    2407             :         }
    2408      138951 :         gcoeff(x,j,i) = gen_0;
    2409      138951 :         gcoeff(x,i,i) = d;
    2410      138951 :         if (U) update(u,v,a,b,(GEN*)(U+i),(GEN*)(U+j));
    2411      138951 :         if (gc_needed(av,1))
    2412             :         {
    2413          37 :           if (DEBUGMEM>1) pari_warn(warnmem,"[2]: ZM_snfall, i = %ld", i);
    2414          37 :           snf_pile1(av, &x,&U);
    2415             :         }
    2416      138951 :         c = 1;
    2417             :       }
    2418      220350 :       if (!c)
    2419             :       {
    2420      141108 :         k = ZM_snf_no_divide(x, i);
    2421      141108 :         if (!k) break;
    2422             : 
    2423             :         /* x[k,j] != 0 mod b */
    2424       69692 :         for (j = 1; j <= i; j++)
    2425       53634 :           gcoeff(x,i,j) = addii(gcoeff(x,i,j),gcoeff(x,k,j));
    2426       16058 :         if (U) gel(U,i) = gadd(gel(U,i),gel(U,k));
    2427             :       }
    2428       95300 :       ZM_redpart(x, mdet, i);
    2429       95300 :       if (U && (flag & 2)) ZM_redpart(U, mdet, n);
    2430       95300 :       if (gc_needed(av,1))
    2431             :       {
    2432           0 :         if (DEBUGMEM>1) pari_warn(warnmem,"[3]: ZM_snfall");
    2433           0 :         snf_pile1(av, &x,&U);
    2434             :       }
    2435             :     }
    2436             :   }
    2437       94925 :   if (DEBUGLEVEL>7) err_printf("\n");
    2438      314900 :   for (k = n; k; k--)
    2439             :   {
    2440      219975 :     GEN d = gcdii(gcoeff(x,k,k), mdet);
    2441      219975 :     gcoeff(x,k,k) = d;
    2442      219975 :     if (!is_pm1(d)) mdet = diviiexact(mdet,d);
    2443             :   }
    2444       94925 : THEEND:
    2445       94995 :   if (U) U = shallowtrans(U);
    2446       94995 :   if (ptV && A)
    2447             :   { /* U A V = D => D^(-1) U A = V^(-1) */
    2448       90508 :     long l = lg(x);
    2449       90508 :     GEN W = ZM_mul(U, A);
    2450      213111 :     for (i = 1; i < l; i++)
    2451             :     {
    2452      180116 :       GEN c = gcoeff(x,i,i);
    2453      180116 :       if (is_pm1(c)) break; /* only 1 from now on */
    2454      122603 :       ZMrow_divexact_inplace(W, i, c);
    2455             :     }
    2456       90508 :     if (flag & 2)
    2457             :     {
    2458       82549 :       W = FpM_red(W, gcoeff(x,1,1));
    2459       82549 :       W = matinvmod(W, gcoeff(x,1,1));
    2460             :     }
    2461             :     else
    2462        7959 :       W = ZM_inv(W, NULL);
    2463       90508 :     V = V? ZM_mul(V, W): W;
    2464             :   }
    2465       94995 :   if (return_vec)
    2466             :   {
    2467       86854 :     long l = lg(x)-1;
    2468       86854 :     if (typ(x) == t_MAT) x = RgM_diagonal_shallow(x);
    2469       86854 :     if (m0 > l) x = shallowconcat(zerovec(m0-l), x);
    2470             :   }
    2471             : 
    2472       94995 :   if (V0)
    2473             :   { /* add kernel */
    2474          21 :     if (!return_vec) x = shallowconcat(zeromat(m,n0-n), x);
    2475          21 :     if (ptV) V = shallowconcat(V0, V);
    2476             :   }
    2477       94995 :   if (perm && U) U = vecpermute(U, perm_inv(perm));
    2478       94995 :   snf_pile(av0, &x,&U,&V);
    2479       94995 :   if (ptU) *ptU = U;
    2480       94995 :   if (ptV) *ptV = V;
    2481       94995 :   return x;
    2482             : }
    2483             : GEN
    2484       12215 : ZM_snfall(GEN x, GEN *U, GEN *V) { return ZM_snfall_i(x, U, V, 0); }
    2485             : GEN
    2486        1330 : ZM_snf(GEN x) { return ZM_snfall_i(x, NULL,NULL, 1); }
    2487             : GEN
    2488         385 : smith(GEN x) { return ZM_snfall_i(x, NULL,NULL, 1); }
    2489             : GEN
    2490          28 : smithall(GEN x)
    2491             : {
    2492          28 :   GEN z = cgetg(4, t_VEC);
    2493          28 :   gel(z,3) = ZM_snfall_i(x, (GEN*)(z+1),(GEN*)(z+2), 0);
    2494          28 :   return z;
    2495             : }
    2496             : 
    2497             : void
    2498      182307 : ZM_snfclean(GEN d, GEN u, GEN v)
    2499             : {
    2500      182307 :   long i, c, l = lg(d);
    2501             : 
    2502      182307 :   if (typ(d) == t_VEC)
    2503      519628 :     for (c=1; c<l; c++) { GEN t = gel(d,c); if (is_pm1(t)) break; }
    2504             :   else
    2505             :   {
    2506           0 :     for (c=1; c<l; c++) { GEN t = gcoeff(d,c,c); if (is_pm1(t)) break; }
    2507           0 :     if (c < l) for (i = 1; i < c; i++) setlg(gel(d,i), c);
    2508             :   }
    2509      182307 :   setlg(d, c);
    2510      575119 :   if (u) for (i=1; i<l; i++) setlg(gel(u,i), c);
    2511      182307 :   if (v) setlg(v, c);
    2512      182307 : }
    2513             : 
    2514             : /* Assume z was computed by [g]smithall(). Remove the 1s on the diagonal */
    2515             : GEN
    2516         539 : smithclean(GEN z)
    2517             : {
    2518             :   long i, j, h, l, c, d;
    2519             :   GEN U, V, y, D, t;
    2520             : 
    2521         539 :   if (typ(z) != t_VEC) pari_err_TYPE("smithclean",z);
    2522         539 :   l = lg(z); if (l == 1) return cgetg(1,t_VEC);
    2523         532 :   U = gel(z,1);
    2524         532 :   if (l != 4 || typ(U) != t_MAT)
    2525             :   { /* assume z = vector of elementary divisors */
    2526        1470 :     for (c=1; c<l; c++)
    2527        1309 :       if (gequal1(gel(z,c))) break;
    2528         511 :     return gcopy_lg(z, c);
    2529             :   }
    2530          21 :   V = gel(z,2);
    2531          21 :   D = gel(z,3);
    2532          21 :   l = lg(D);
    2533          21 :   if (l == 1) return gcopy(z);
    2534          21 :   h = lgcols(D);
    2535          21 :   if (h > l)
    2536             :   { /* D = vconcat(zero matrix, diagonal matrix) */
    2537          21 :     for (c=1+h-l, d=1; c<h; c++,d++)
    2538          21 :       if (gequal1(gcoeff(D,c,d))) break;
    2539             :   }
    2540           7 :   else if (h < l)
    2541             :   { /* D = concat(zero matrix, diagonal matrix) */
    2542           7 :     for (c=1, d=1+l-h; d<l; c++,d++)
    2543           7 :       if (gequal1(gcoeff(D,c,d))) break;
    2544             :   }
    2545             :   else
    2546             :   { /* D diagonal */
    2547           0 :     for (c=1; c<l; c++)
    2548           0 :       if (gequal1(gcoeff(D,c,c))) break;
    2549           0 :     d = c;
    2550             :   }
    2551             :   /* U was (h-1)x(h-1), V was (l-1)x(l-1), D was (h-1)x(l-1) */
    2552          21 :   y = cgetg(4,t_VEC);
    2553             :   /* truncate U to (c-1) x (h-1) */
    2554          21 :   gel(y,1) = t = cgetg(h,t_MAT);
    2555          77 :   for (j=1; j<h; j++) gel(t,j) = gcopy_lg(gel(U,j), c);
    2556             :   /* truncate V to (l-1) x (d-1) */
    2557          21 :   gel(y,2) = gcopy_lg(V, d);
    2558          21 :   gel(y,3) = t = zeromatcopy(c-1, d-1);
    2559             :   /* truncate D to a (c-1) x (d-1) matrix */
    2560          21 :   if (d > 1)
    2561             :   {
    2562          14 :     if (h > l)
    2563             :     {
    2564          14 :       for (i=1+h-l, j=1; i<c; i++,j++)
    2565           7 :         gcoeff(t,i,j) = gcopy(gcoeff(D,i,j));
    2566             :     }
    2567           7 :     else if (h < l)
    2568             :     {
    2569           7 :       for (i=1, j=1+l-h; j<d; i++,j++)
    2570           0 :         gcoeff(t,i,j) = gcopy(gcoeff(D,i,j));
    2571             :     }
    2572             :     else
    2573             :     {
    2574           0 :       for (j=1; j<d; j++)
    2575           0 :         gcoeff(t,j,j) = gcopy(gcoeff(D,j,j));
    2576             :     }
    2577             :   }
    2578          21 :   return y;
    2579             : }
    2580             : 
    2581             : /* does b = x[i,i] divide all entries in x[1..i-1,1..i-1] ? If so, return 0;
    2582             :  * else return the index of a problematic row */
    2583             : static long
    2584         112 : gsnf_no_divide(GEN x, long i, long vx)
    2585             : {
    2586         112 :   GEN b = gcoeff(x,i,i);
    2587             :   long j, k;
    2588             : 
    2589         112 :   if (gequal0(b))
    2590             :   {
    2591          14 :     for (k = 1; k < i; k++)
    2592          14 :       for (j = 1; j < i; j++)
    2593          14 :         if (!gequal0(gcoeff(x,k,j))) return k;
    2594           0 :     return 0;
    2595             :   }
    2596             : 
    2597          98 :   if (!is_RgX(b,vx) || degpol(b)<=0) return 0;
    2598         133 :   for (k = 1; k < i; k++)
    2599         231 :     for (j = 1; j < i; j++)
    2600             :     {
    2601         154 :       GEN z = gcoeff(x,k,j), r;
    2602         154 :       if (!is_RgX(z,vx)) z = scalarpol(z, vx);
    2603         154 :       r = RgX_rem(z, b);
    2604         154 :       if (signe(r) && (! isinexactreal(r) ||
    2605           0 :              gexpo(r) > 16 + gexpo(b) - prec2nbits(gprecision(r)))
    2606           7 :          ) return k;
    2607             :     }
    2608          49 :   return 0;
    2609             : }
    2610             : 
    2611             : /* Hermite Normal Form, with base change matrix if ptB != NULL.
    2612             :  * If 'remove' = 1, remove 0 columns (do NOT update *ptB accordingly)
    2613             :  * If 'remove' = 2, remove 0 columns and update *ptB accordingly */
    2614             : GEN
    2615          35 : RgM_hnfall(GEN A, GEN *pB, long remove)
    2616             : {
    2617             :   pari_sp av;
    2618             :   long li, j, k, m, n, def, ldef;
    2619             :   GEN B;
    2620          35 :   long vx = gvar(A);
    2621             : 
    2622          35 :   n = lg(A)-1;
    2623          35 :   if (vx==NO_VARIABLE || !n)
    2624             :   {
    2625           0 :     RgM_check_ZM(A, "mathnf0");
    2626           0 :     return ZM_hnfall(A, pB, remove);
    2627             :   }
    2628          35 :   m = nbrows(A);
    2629          35 :   av = avma;
    2630          35 :   A = RgM_shallowcopy(A);
    2631          35 :   B = pB? matid(n): NULL;
    2632          35 :   def = n; ldef = (m>n)? m-n: 0;
    2633          77 :   for (li=m; li>ldef; li--)
    2634             :   {
    2635             :     GEN d, T;
    2636          63 :     for (j=def-1; j; j--)
    2637             :     {
    2638          21 :       GEN a = gcoeff(A,li,j);
    2639          21 :       if (gequal0(a)) continue;
    2640             : 
    2641           7 :       k = (j==1)? def: j-1;
    2642           7 :       RgC_elem(a,gcoeff(A,li,k), A,B, j,k, li, vx);
    2643             :     }
    2644          42 :     T = normalize_as_RgX(gcoeff(A,li,def), vx, &d);
    2645          42 :     if (gequal0(T))
    2646           0 :     { if (ldef) ldef--; }
    2647             :     else
    2648             :     {
    2649          42 :       gcoeff(A,li,def) = T;
    2650          42 :       if (B && !gequal1(d)) gel(B, def) = RgC_Rg_div(gel(B, def), d);
    2651          42 :       RgM_reduce(A, B, li, def, vx);
    2652          42 :       def--;
    2653             :     }
    2654          42 :     if (gc_needed(av,1))
    2655             :     {
    2656           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"ghnfall");
    2657           0 :       gerepileall(av, B? 2: 1, &A, &B);
    2658             :     }
    2659             :   }
    2660             :   /* rank A = n - def */
    2661          35 :   if (remove) remove_0cols(def, &A, &B, remove);
    2662          35 :   gerepileall(av, B? 2: 1, &A, &B);
    2663          35 :   if (B) *pB = B;
    2664          35 :   return A;
    2665             : }
    2666             : 
    2667             : static GEN
    2668          35 : RgXM_snf(GEN x,long all)
    2669             : {
    2670             :   pari_sp av;
    2671             :   long i, j, k, n;
    2672             :   GEN z, u, v, U, V;
    2673          35 :   long vx = gvar(x);
    2674          35 :   n = lg(x)-1; if (!n) return trivsmith(all);
    2675          35 :   if (vx==NO_VARIABLE) pari_err_TYPE("RgXM_snf",x);
    2676          35 :   if (lgcols(x) != n+1) pari_err_DIM("gsmithall");
    2677          35 :   av = avma;
    2678          35 :   x = RgM_shallowcopy(x);
    2679          35 :   if (all) { U = matid(n); V = matid(n); }
    2680         126 :   for (i=n; i>=2; i--)
    2681             :   {
    2682             :     for(;;)
    2683          98 :     {
    2684             :       GEN a, b, d;
    2685         189 :       int c = 0;
    2686         546 :       for (j=i-1; j>=1; j--)
    2687             :       {
    2688         357 :         b = gcoeff(x,i,j); if (gequal0(b)) continue;
    2689         140 :         a = gcoeff(x,i,i);
    2690         140 :         d = gbezout_step(&b, &a, &v, &u, vx);
    2691         469 :         for (k = 1; k < i; k++)
    2692             :         {
    2693         329 :           GEN t = gadd(gmul(u,gcoeff(x,k,i)),gmul(v,gcoeff(x,k,j)));
    2694         329 :           gcoeff(x,k,j) = gsub(gmul(a,gcoeff(x,k,j)),gmul(b,gcoeff(x,k,i)));
    2695         329 :           gcoeff(x,k,i) = t;
    2696             :         }
    2697         140 :         gcoeff(x,i,j) = gen_0;
    2698         140 :         gcoeff(x,i,i) = d;
    2699         140 :         if (all) update(u,v,a,b,(GEN*)(V+i),(GEN*)(V+j));
    2700             :       }
    2701         546 :       for (j=i-1; j>=1; j--)
    2702             :       {
    2703         357 :         b = gcoeff(x,j,i); if (gequal0(b)) continue;
    2704         119 :         a = gcoeff(x,i,i);
    2705         119 :         d = gbezout_step(&b, &a, &v, &u, vx);
    2706         406 :         for (k = 1; k < i; k++)
    2707             :         {
    2708         287 :           GEN t = gadd(gmul(u,gcoeff(x,i,k)),gmul(v,gcoeff(x,j,k)));
    2709         287 :           gcoeff(x,j,k) = gsub(gmul(a,gcoeff(x,j,k)),gmul(b,gcoeff(x,i,k)));
    2710         287 :           gcoeff(x,i,k) = t;
    2711             :         }
    2712         119 :         gcoeff(x,j,i) = gen_0;
    2713         119 :         gcoeff(x,i,i) = d;
    2714         119 :         if (all) update(u,v,a,b,(GEN*)(U+i),(GEN*)(U+j));
    2715         119 :         c = 1;
    2716             :       }
    2717         189 :       if (!c)
    2718             :       {
    2719         112 :         k = gsnf_no_divide(x, i, vx);
    2720         112 :         if (!k) break;
    2721             : 
    2722          70 :         for (j=1; j<=i; j++)
    2723          49 :           gcoeff(x,i,j) = gadd(gcoeff(x,i,j),gcoeff(x,k,j));
    2724          21 :         if (all) gel(U,i) = gadd(gel(U,i),gel(U,k));
    2725             :       }
    2726          98 :       if (gc_needed(av,1))
    2727             :       {
    2728           0 :         if (DEBUGMEM>1) pari_warn(warnmem,"gsmithall");
    2729           0 :         gerepileall(av, all? 3: 1, &x, &U, &V);
    2730             :       }
    2731             :     }
    2732             :   }
    2733         161 :   for (k=1; k<=n; k++)
    2734             :   {
    2735         126 :     GEN d, T = normalize_as_RgX(gcoeff(x,k,k), vx, &d);
    2736         126 :     if (gequal0(T)) continue;
    2737         112 :     if (all && !gequal1(d)) gel(V,k) = RgC_Rg_div(gel(V,k), d);
    2738         112 :     gcoeff(x,k,k) = T;
    2739             :   }
    2740          35 :   z = all? mkvec3(shallowtrans(U), V, x): RgM_diagonal_shallow(x);
    2741          35 :   return gerepilecopy(av, z);
    2742             : }
    2743             : 
    2744             : GEN
    2745         448 : matsnf0(GEN x,long flag)
    2746             : {
    2747         448 :   pari_sp av = avma;
    2748         448 :   if (flag > 7) pari_err_FLAG("matsnf");
    2749         448 :   if (typ(x) == t_VEC && flag & 4) return smithclean(x);
    2750         448 :   if (typ(x)!=t_MAT) pari_err_TYPE("matsnf",x);
    2751         448 :   if (RgM_is_ZM(x)) x = flag&1 ? smithall(x): smith(x);
    2752          35 :   else              x = RgXM_snf(x, flag&1);
    2753         448 :   if (flag & 4) x = gerepileupto(av, smithclean(x));
    2754         448 :   return x;
    2755             : }
    2756             : GEN
    2757           0 : gsmith(GEN x) { return RgXM_snf(x,0); }
    2758             : GEN
    2759           0 : gsmithall(GEN x) { return RgXM_snf(x,1); }
    2760             : 
    2761             : /* H is a relation matrix, either in HNF or a t_VEC (diagonal HNF) */
    2762             : static GEN
    2763      182307 : snf_group(GEN H, GEN D, GEN *newU, GEN *newUi)
    2764             : {
    2765             :   long i, j, l;
    2766             : 
    2767      182307 :   ZM_snfclean(D, newU? *newU: NULL, newUi? *newUi: NULL);
    2768      182307 :   l = lg(D);
    2769      182307 :   if (newU) {
    2770      165081 :     GEN U = *newU;
    2771      473011 :     for (i = 1; i < l; i++)
    2772             :     {
    2773      307930 :       GEN d = gel(D,i), d2 = shifti(d, 1);
    2774     1212582 :       for (j = 1; j < lg(U); j++)
    2775      904652 :         gcoeff(U,i,j) = centermodii(gcoeff(U,i,j), d, d2);
    2776             :     }
    2777      165081 :     *newU = U;
    2778             :   }
    2779      182307 :   if (newUi && l > 1)
    2780             :   { /* UHV=D -> U^-1 = (HV)D^-1 -> U^-1 = H(VD^-1 mod 1) mod H */
    2781             :     /* Ui = ZM_inv(U, NULL); setlg(Ui, l); */
    2782      174620 :     GEN V = *newUi, Ui;
    2783      174620 :     int Hvec = (typ(H) == t_VEC);
    2784      511528 :     for (i = 1; i < l; i++) gel(V,i) = FpC_red(gel(V,i), gel(D,i));
    2785      174620 :     if (!Hvec)
    2786             :     {
    2787       75689 :       if (ZM_isdiagonal(H)) { H = RgM_diagonal_shallow(H); Hvec = 1; }
    2788             :     }
    2789      174620 :     Ui = Hvec? ZM_diag_mul(H, V): ZM_mul(H, V);
    2790      511528 :     for (i = 1; i < l; i++) gel(Ui,i) = ZC_Z_divexact(gel(Ui,i), gel(D,i));
    2791      174620 :     if (Hvec)
    2792      401990 :     { for (i = 1; i < l; i++) gel(Ui,i) = vecmodii(gel(Ui,i), H); }
    2793             :     else
    2794       42166 :       Ui = ZM_hnfrem(Ui, H);
    2795      174620 :     *newUi = Ui;
    2796             :   }
    2797      182307 :   return D;
    2798             : }
    2799             : /* H relation matrix among row of generators g in HNF.  Let URV = D its SNF,
    2800             :  * newU R newV = newD its clean SNF (no 1 in Dnew). Return the diagonal of
    2801             :  * newD, newU and newUi such that  1/U = (newUi, ?).
    2802             :  * Rationale: let (G,0) = g Ui be the new generators then
    2803             :  * 0 = G U R --> G D = 0,  g = G newU  and  G = g newUi */
    2804             : GEN
    2805       83180 : ZM_snf_group(GEN H, GEN *newU, GEN *newUi)
    2806             : {
    2807       83180 :   GEN D = ZM_snfall_i(H, newU, newUi, 1 + 2);
    2808       83180 :   return snf_group(H, D, newU, newUi);
    2809             : }
    2810             : 
    2811             : /* D a ZV: SNF for matdiagonal(D). Faster because we only ensure elementary
    2812             :  * divisors condition: d[n] | ... | d[1] and need not convert D to matrix form*/
    2813             : GEN
    2814       99127 : ZV_snfall(GEN D, GEN *pU, GEN *pV)
    2815             : {
    2816       99127 :   pari_sp av = avma;
    2817       99127 :   long j, n = lg(D)-1;
    2818       99127 :   GEN U = pU? matid(n): NULL;
    2819       99127 :   GEN V = pV? matid(n): NULL;
    2820             :   GEN p;
    2821             : 
    2822       99127 :   D = leafcopy(D);
    2823      325192 :   for (j = n; j > 0; j--)
    2824             :   {
    2825      226065 :     GEN b = gel(D,j);
    2826      226065 :     if (signe(b) < 0)
    2827             :     {
    2828           0 :       gel(D,j) = negi(b);
    2829           0 :       if (V) ZV_togglesign(gel(V,j));
    2830             :     }
    2831             :   }
    2832             :   /* entries are non-negative integers */
    2833       99127 :   p = gen_indexsort(D, NULL, &negcmpii);
    2834       99127 :   D = vecpermute(D, p);
    2835       99127 :   if (U) U = vecpermute(U, p);
    2836       99127 :   if (V) V = vecpermute(V, p);
    2837             :   /* entries are sorted by decreasing value */
    2838      325192 :   for (j = n; j > 0; j--)
    2839             :   {
    2840      226065 :     GEN b = gel(D,j);
    2841             :     long i;
    2842      416892 :     for (i = j-1; i > 0; i--)
    2843             :     { /* invariant: a >= b. If au+bv = d is a Bezout relation, A=a/d and B=b/d
    2844             :        * we have [B,-A;u,v]*diag(a,b)*[1-u*A,1; -u*A,1]] = diag(Ab, d) */
    2845      192031 :       GEN a = gel(D,i), u,v, d = bezout(a,b, &u,&v), A, Wi, Wj;
    2846      192031 :       if (equalii(d,b)) continue;
    2847       25403 :       A = diviiexact(a,d);
    2848       25403 :       if (V)
    2849             :       {
    2850       25361 :         GEN t = mulii(u,A);
    2851       25361 :         Wi = ZC_lincomb(subui(1,t), negi(t), gel(V,i), gel(V,j));
    2852       25361 :         Wj = ZC_add(gel(V,i), gel(V,j));
    2853       25361 :         gel(V,i) = Wi;
    2854       25361 :         gel(V,j) = Wj;
    2855             :       }
    2856       25403 :       if (U)
    2857             :       {
    2858       25403 :         GEN B = diviiexact(b,d);
    2859       25403 :         Wi = ZC_lincomb(B, negi(A), gel(U,i), gel(U,j));
    2860       25403 :         Wj = ZC_lincomb(u, v, gel(U,i), gel(U,j));
    2861       25403 :         gel(U,i) = Wi;
    2862       25403 :         gel(U,j) = Wj;
    2863             :       }
    2864       25403 :       gel(D,i) = mulii(A,b); /* lcm(a,b) */
    2865       25403 :       gel(D,j) = d; /* gcd(a,b) */
    2866       25403 :       b = gel(D,j); if (equali1(b)) break;
    2867             :     }
    2868             :   }
    2869       99127 :   snf_pile(av, &D,&U,&V);
    2870       99127 :   if (U) *pU = shallowtrans(U);
    2871       99127 :   if (V) *pV = V;
    2872       99127 :   return D;
    2873             : }
    2874             : GEN
    2875       99127 : ZV_snf_group(GEN d, GEN *newU, GEN *newUi)
    2876             : {
    2877       99127 :   GEN D = ZV_snfall(d, newU, newUi);
    2878       99127 :   return snf_group(d, D, newU, newUi);
    2879             : }
    2880             : 
    2881             : /* D a vector of elementary divisors. Truncate (setlg) to leave out trivial
    2882             :  * entries (= 1) */
    2883             : void
    2884           0 : ZV_snf_trunc(GEN D)
    2885             : {
    2886           0 :   long i, l = lg(D);
    2887           0 :   for (i = 1; i < l; i++)
    2888           0 :     if (is_pm1(gel(D,i))) { setlg(D,i); break; }
    2889           0 : }

Generated by: LCOV version 1.13