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 - ZV.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.12.1 lcov report (development 25406-bf255ab81b) Lines: 801 901 88.9 %
Date: 2020-06-04 05:59:24 Functions: 120 130 92.3 %
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             : static int
      18     1275187 : check_ZV(GEN x, long l)
      19             : {
      20             :   long i;
      21     8090537 :   for (i=1; i<l; i++)
      22     6815441 :     if (typ(gel(x,i)) != t_INT) return 0;
      23     1275096 :   return 1;
      24             : }
      25             : void
      26      512414 : RgV_check_ZV(GEN A, const char *s)
      27             : {
      28      512414 :   if (!RgV_is_ZV(A)) pari_err_TYPE(stack_strcat(s," [integer vector]"), A);
      29      512407 : }
      30             : void
      31      368665 : RgM_check_ZM(GEN A, const char *s)
      32             : {
      33      368665 :   long n = lg(A);
      34      368665 :   if (n != 1)
      35             :   {
      36      368609 :     long j, m = lgcols(A);
      37     1643705 :     for (j=1; j<n; j++)
      38     1275187 :       if (!check_ZV(gel(A,j), m))
      39          91 :         pari_err_TYPE(stack_strcat(s," [integer matrix]"), A);
      40             :   }
      41      368574 : }
      42             : 
      43             : static long
      44    62677751 : ZV_max_lg_i(GEN x, long m)
      45             : {
      46    62677751 :   long i, prec = 2;
      47   476771873 :   for (i=1; i<m; i++) { long l = lgefint(gel(x,i)); if (l > prec) prec = l; }
      48    62677751 :   return prec;
      49             : }
      50             : 
      51             : long
      52       10612 : ZV_max_lg(GEN x)
      53       10612 : { return ZV_max_lg_i(x, lg(x)); }
      54             : 
      55             : static long
      56    21681811 : ZM_max_lg_i(GEN x, long n, long m)
      57             : {
      58    21681811 :   long prec = 2;
      59    21681811 :   if (n != 1)
      60             :   {
      61             :     long j;
      62    84348950 :     for (j=1; j<n; j++)
      63             :     {
      64    62667139 :       long l = ZV_max_lg_i(gel(x,j), m);
      65    62667139 :       if (l > prec) prec = l;
      66             :     }
      67             :   }
      68    21681811 :   return prec;
      69             : }
      70             : 
      71             : long
      72      442191 : ZM_max_lg(GEN x)
      73             : {
      74      442191 :   long n = lg(x);
      75      442191 :   if (n==1) return 2;
      76      442191 :   return ZM_max_lg_i(x, n, lgcols(x));
      77             : }
      78             : 
      79             : GEN
      80        3394 : ZM_supnorm(GEN x)
      81             : {
      82        3394 :   long i, j, h, lx = lg(x);
      83        3394 :   GEN s = gen_0;
      84        3394 :   if (lx == 1) return gen_1;
      85        3394 :   h = lgcols(x);
      86       20533 :   for (j=1; j<lx; j++)
      87             :   {
      88       17139 :     GEN xj = gel(x,j);
      89      240192 :     for (i=1; i<h; i++)
      90             :     {
      91      223053 :       GEN c = gel(xj,i);
      92      223053 :       if (abscmpii(c, s) > 0) s = c;
      93             :     }
      94             :   }
      95        3394 :   return absi(s);
      96             : }
      97             : 
      98             : /********************************************************************/
      99             : /**                                                                **/
     100             : /**                           MULTIPLICATION                       **/
     101             : /**                                                                **/
     102             : /********************************************************************/
     103             : /* x non-empty ZM, y a compatible nc (dimension > 0). */
     104             : static GEN
     105           0 : ZM_nc_mul_i(GEN x, GEN y, long c, long l)
     106             : {
     107             :   long i, j;
     108             :   pari_sp av;
     109           0 :   GEN z = cgetg(l,t_COL), s;
     110             : 
     111           0 :   for (i=1; i<l; i++)
     112             :   {
     113           0 :     av = avma; s = muliu(gcoeff(x,i,1),y[1]);
     114           0 :     for (j=2; j<c; j++)
     115           0 :       if (y[j]) s = addii(s, muliu(gcoeff(x,i,j),y[j]));
     116           0 :     gel(z,i) = gerepileuptoint(av,s);
     117             :   }
     118           0 :   return z;
     119             : }
     120             : 
     121             : /* x ZV, y a compatible zc. */
     122             : GEN
     123     2214660 : ZV_zc_mul(GEN x, GEN y)
     124             : {
     125     2214660 :   long j, l = lg(x);
     126     2214660 :   pari_sp av = avma;
     127     2214660 :   GEN s = mulis(gel(x,1),y[1]);
     128    50028650 :   for (j=2; j<l; j++)
     129    47813990 :     if (y[j]) s = addii(s, mulis(gel(x,j),y[j]));
     130     2214660 :   return gerepileuptoint(av,s);
     131             : }
     132             : 
     133             : /* x non-empty ZM, y a compatible zc (dimension > 0). */
     134             : static GEN
     135     5058706 : ZM_zc_mul_i(GEN x, GEN y, long c, long l)
     136             : {
     137             :   long i, j;
     138     5058706 :   GEN z = cgetg(l,t_COL);
     139             : 
     140    28716711 :   for (i=1; i<l; i++)
     141             :   {
     142    23658005 :     pari_sp av = avma;
     143    23658005 :     GEN s = mulis(gcoeff(x,i,1),y[1]);
     144   476033906 :     for (j=2; j<c; j++)
     145   452375901 :       if (y[j]) s = addii(s, mulis(gcoeff(x,i,j),y[j]));
     146    23658005 :     gel(z,i) = gerepileuptoint(av,s);
     147             :   }
     148     5058706 :   return z;
     149             : }
     150             : GEN
     151     3508584 : ZM_zc_mul(GEN x, GEN y) {
     152     3508584 :   long l = lg(x);
     153     3508584 :   if (l == 1) return cgetg(1, t_COL);
     154     3508584 :   return ZM_zc_mul_i(x,y, l, lgcols(x));
     155             : }
     156             : 
     157             : /* y non-empty ZM, x a compatible zv (dimension > 0). */
     158             : GEN
     159        1260 : zv_ZM_mul(GEN x, GEN y) {
     160        1260 :   long i,j, lx = lg(x), ly = lg(y);
     161             :   GEN z;
     162        1260 :   if (lx == 1) return zerovec(ly-1);
     163        1260 :   z = cgetg(ly,t_VEC);
     164        3080 :   for (j=1; j<ly; j++)
     165             :   {
     166        1820 :     pari_sp av = avma;
     167        1820 :     GEN s = mulsi(x[1], gcoeff(y,1,j));
     168        3472 :     for (i=2; i<lx; i++)
     169        1652 :       if (x[i]) s = addii(s, mulsi(x[i], gcoeff(y,i,j)));
     170        1820 :     gel(z,j) = gerepileuptoint(av,s);
     171             :   }
     172        1260 :   return z;
     173             : }
     174             : 
     175             : /* x ZM, y a compatible zm (dimension > 0). */
     176             : GEN
     177      755578 : ZM_zm_mul(GEN x, GEN y)
     178             : {
     179      755578 :   long j, c, l = lg(x), ly = lg(y);
     180      755578 :   GEN z = cgetg(ly, t_MAT);
     181      755578 :   if (l == 1) return z;
     182      755578 :   c = lgcols(x);
     183     2305700 :   for (j = 1; j < ly; j++) gel(z,j) = ZM_zc_mul_i(x, gel(y,j), l,c);
     184      755578 :   return z;
     185             : }
     186             : /* x ZM, y a compatible zn (dimension > 0). */
     187             : GEN
     188           0 : ZM_nm_mul(GEN x, GEN y)
     189             : {
     190           0 :   long j, c, l = lg(x), ly = lg(y);
     191           0 :   GEN z = cgetg(ly, t_MAT);
     192           0 :   if (l == 1) return z;
     193           0 :   c = lgcols(x);
     194           0 :   for (j = 1; j < ly; j++) gel(z,j) = ZM_nc_mul_i(x, gel(y,j), l,c);
     195           0 :   return z;
     196             : }
     197             : 
     198             : /* Strassen-Winograd algorithm */
     199             : 
     200             : /*
     201             :   Return A[ma+1..ma+da, na+1..na+ea] - B[mb+1..mb+db, nb+1..nb+eb]
     202             :   as an (m x n)-matrix, padding the input with zeroes as necessary.
     203             : */
     204             : static GEN
     205      345328 : add_slices(long m, long n,
     206             :            GEN A, long ma, long da, long na, long ea,
     207             :            GEN B, long mb, long db, long nb, long eb)
     208             : {
     209      345328 :   long min_d = minss(da, db), min_e = minss(ea, eb), i, j;
     210      345328 :   GEN M = cgetg(n + 1, t_MAT), C;
     211             : 
     212     2671391 :   for (j = 1; j <= min_e; j++) {
     213     2326063 :     gel(M, j) = C = cgetg(m + 1, t_COL);
     214    36247456 :     for (i = 1; i <= min_d; i++)
     215    33921393 :       gel(C, i) = addii(gcoeff(A, ma + i, na + j),
     216    33921393 :                         gcoeff(B, mb + i, nb + j));
     217     2354979 :     for (; i <= da; i++)
     218       28916 :       gel(C, i) = gcoeff(A, ma + i, na + j);
     219     2326063 :     for (; i <= db; i++)
     220           0 :       gel(C, i) = gcoeff(B, mb + i, nb + j);
     221     2326063 :     for (; i <= m; i++)
     222           0 :       gel(C, i) = gen_0;
     223             :   }
     224      365979 :   for (; j <= ea; j++) {
     225       20651 :     gel(M, j) = C = cgetg(m + 1, t_COL);
     226       76826 :     for (i = 1; i <= da; i++)
     227       56175 :       gel(C, i) = gcoeff(A, ma + i, na + j);
     228       20651 :     for (; i <= m; i++)
     229           0 :       gel(C, i) = gen_0;
     230             :   }
     231      345328 :   for (; j <= eb; j++) {
     232           0 :     gel(M, j) = C = cgetg(m + 1, t_COL);
     233           0 :     for (i = 1; i <= db; i++)
     234           0 :       gel(C, i) = gcoeff(B, mb + i, nb + j);
     235           0 :     for (; i <= m; i++)
     236           0 :       gel(C, i) = gen_0;
     237             :   }
     238      345328 :   for (; j <= n; j++)
     239           0 :     gel(M, j) = zerocol(m);
     240      345328 :   return M;
     241             : }
     242             : 
     243             : /*
     244             :   Return A[ma+1..ma+da, na+1..na+ea] - B[mb+1..mb+db, nb+1..nb+eb]
     245             :   as an (m x n)-matrix, padding the input with zeroes as necessary.
     246             : */
     247             : static GEN
     248      302162 : subtract_slices(long m, long n,
     249             :                 GEN A, long ma, long da, long na, long ea,
     250             :                 GEN B, long mb, long db, long nb, long eb)
     251             : {
     252      302162 :   long min_d = minss(da, db), min_e = minss(ea, eb), i, j;
     253      302162 :   GEN M = cgetg(n + 1, t_MAT), C;
     254             : 
     255     2325370 :   for (j = 1; j <= min_e; j++) {
     256     2023208 :     gel(M, j) = C = cgetg(m + 1, t_COL);
     257    32769770 :     for (i = 1; i <= min_d; i++)
     258    30746562 :       gel(C, i) = subii(gcoeff(A, ma + i, na + j),
     259    30746562 :                         gcoeff(B, mb + i, nb + j));
     260     2043388 :     for (; i <= da; i++)
     261       20180 :       gel(C, i) = gcoeff(A, ma + i, na + j);
     262     2055539 :     for (; i <= db; i++)
     263       32331 :       gel(C, i) = negi(gcoeff(B, mb + i, nb + j));
     264     2023208 :     for (; i <= m; i++)
     265           0 :       gel(C, i) = gen_0;
     266             :   }
     267      302162 :   for (; j <= ea; j++) {
     268           0 :     gel(M, j) = C = cgetg(m + 1, t_COL);
     269           0 :     for (i = 1; i <= da; i++)
     270           0 :       gel(C, i) = gcoeff(A, ma + i, na + j);
     271           0 :     for (; i <= m; i++)
     272           0 :       gel(C, i) = gen_0;
     273             :   }
     274      317870 :   for (; j <= eb; j++) {
     275       15708 :     gel(M, j) = C = cgetg(m + 1, t_COL);
     276       49918 :     for (i = 1; i <= db; i++)
     277       34210 :       gel(C, i) = negi(gcoeff(B, mb + i, nb + j));
     278       15708 :     for (; i <= m; i++)
     279           0 :       gel(C, i) = gen_0;
     280             :   }
     281      317870 :   for (; j <= n; j++)
     282       15708 :     gel(M, j) = zerocol(m);
     283      302162 :   return M;
     284             : }
     285             : 
     286             : static GEN ZM_mul_i(GEN x, GEN y, long l, long lx, long ly);
     287             : 
     288             : /* Strassen-Winograd matrix product A (m x n) * B (n x p) */
     289             : static GEN
     290       43166 : ZM_mul_sw(GEN A, GEN B, long m, long n, long p)
     291             : {
     292       43166 :   pari_sp av = avma;
     293       43166 :   long m1 = (m + 1)/2, m2 = m/2,
     294       43166 :     n1 = (n + 1)/2, n2 = n/2,
     295       43166 :     p1 = (p + 1)/2, p2 = p/2;
     296             :   GEN A11, A12, A22, B11, B21, B22,
     297             :     S1, S2, S3, S4, T1, T2, T3, T4,
     298             :     M1, M2, M3, M4, M5, M6, M7,
     299             :     V1, V2, V3, C11, C12, C21, C22, C;
     300             : 
     301       43166 :   T2 = subtract_slices(n1, p2, B, 0, n1, p1, p2, B, n1, n2, p1, p2);
     302       43166 :   S1 = subtract_slices(m2, n1, A, m1, m2, 0, n1, A, 0, m2, 0, n1);
     303       43166 :   M2 = ZM_mul_i(S1, T2, m2 + 1, n1 + 1, p2 + 1);
     304       43166 :   if (gc_needed(av, 1))
     305           0 :     gerepileall(av, 2, &T2, &M2);  /* destroy S1 */
     306       43166 :   T3 = subtract_slices(n1, p1, T2, 0, n1, 0, p2, B, 0, n1, 0, p1);
     307       43166 :   if (gc_needed(av, 1))
     308           0 :     gerepileall(av, 2, &M2, &T3);  /* destroy T2 */
     309       43166 :   S2 = add_slices(m2, n1, A, m1, m2, 0, n1, A, m1, m2, n1, n2);
     310       43166 :   T1 = subtract_slices(n1, p1, B, 0, n1, p1, p2, B, 0, n1, 0, p2);
     311       43166 :   M3 = ZM_mul_i(S2, T1, m2 + 1, n1 + 1, p2 + 1);
     312       43166 :   if (gc_needed(av, 1))
     313           0 :     gerepileall(av, 4, &M2, &T3, &S2, &M3);  /* destroy T1 */
     314       43166 :   S3 = subtract_slices(m1, n1, S2, 0, m2, 0, n1, A, 0, m1, 0, n1);
     315       43166 :   if (gc_needed(av, 1))
     316           0 :     gerepileall(av, 4, &M2, &T3, &M3, &S3);  /* destroy S2 */
     317       43166 :   A11 = matslice(A, 1, m1, 1, n1);
     318       43166 :   B11 = matslice(B, 1, n1, 1, p1);
     319       43166 :   M1 = ZM_mul_i(A11, B11, m1 + 1, n1 + 1, p1 + 1);
     320       43166 :   if (gc_needed(av, 1))
     321           0 :     gerepileall(av, 5, &M2, &T3, &M3, &S3, &M1);  /* destroy A11, B11 */
     322       43166 :   A12 = matslice(A, 1, m1, n1 + 1, n);
     323       43166 :   B21 = matslice(B, n1 + 1, n, 1, p1);
     324       43166 :   M4 = ZM_mul_i(A12, B21, m1 + 1, n2 + 1, p1 + 1);
     325       43166 :   if (gc_needed(av, 1))
     326           0 :     gerepileall(av, 6, &M2, &T3, &M3, &S3, &M1, &M4);  /* destroy A12, B21 */
     327       43166 :   C11 = add_slices(m1, p1, M1, 0, m1, 0, p1, M4, 0, m1, 0, p1);
     328       43166 :   if (gc_needed(av, 1))
     329           0 :     gerepileall(av, 6, &M2, &T3, &M3, &S3, &M1, &C11);  /* destroy M4 */
     330       43166 :   M5 = ZM_mul_i(S3, T3, m1 + 1, n1 + 1, p1 + 1);
     331       43166 :   S4 = subtract_slices(m1, n2, A, 0, m1, n1, n2, S3, 0, m1, 0, n2);
     332       43166 :   if (gc_needed(av, 1))
     333           5 :     gerepileall(av, 7, &M2, &T3, &M3, &M1, &C11, &M5, &S4);  /* destroy S3 */
     334       43166 :   T4 = add_slices(n2, p1, B, n1, n2, 0, p1, T3, 0, n2, 0, p1);
     335       43166 :   if (gc_needed(av, 1))
     336           0 :     gerepileall(av, 7, &M2, &M3, &M1, &C11, &M5, &S4, &T4);  /* destroy T3 */
     337       43166 :   V1 = subtract_slices(m1, p1, M1, 0, m1, 0, p1, M5, 0, m1, 0, p1);
     338       43166 :   if (gc_needed(av, 1))
     339           1 :     gerepileall(av, 6, &M2, &M3, &S4, &T4, &C11, &V1);  /* destroy M1, M5 */
     340       43166 :   B22 = matslice(B, n1 + 1, n, p1 + 1, p);
     341       43166 :   M6 = ZM_mul_i(S4, B22, m1 + 1, n2 + 1, p2 + 1);
     342       43166 :   if (gc_needed(av, 1))
     343           6 :     gerepileall(av, 6, &M2, &M3, &T4, &C11, &V1, &M6);  /* destroy S4, B22 */
     344       43166 :   A22 = matslice(A, m1 + 1, m, n1 + 1, n);
     345       43166 :   M7 = ZM_mul_i(A22, T4, m2 + 1, n2 + 1, p1 + 1);
     346       43166 :   if (gc_needed(av, 1))
     347           0 :     gerepileall(av, 6, &M2, &M3, &C11, &V1, &M6, &M7);  /* destroy A22, T4 */
     348       43166 :   V3 = add_slices(m1, p2, V1, 0, m1, 0, p2, M3, 0, m2, 0, p2);
     349       43166 :   C12 = add_slices(m1, p2, V3, 0, m1, 0, p2, M6, 0, m1, 0, p2);
     350       43166 :   if (gc_needed(av, 1))
     351           6 :     gerepileall(av, 6, &M2, &M3, &C11, &V1, &M7, &C12);  /* destroy V3, M6 */
     352       43166 :   V2 = add_slices(m2, p1, V1, 0, m2, 0, p1, M2, 0, m2, 0, p2);
     353       43166 :   if (gc_needed(av, 1))
     354           0 :     gerepileall(av, 5, &M3, &C11, &M7, &C12, &V2);  /* destroy V1, M2 */
     355       43166 :   C21 = add_slices(m2, p1, V2, 0, m2, 0, p1, M7, 0, m2, 0, p1);
     356       43166 :   if (gc_needed(av, 1))
     357           7 :     gerepileall(av, 5, &M3, &C11, &C12, &V2, &C21);  /* destroy M7 */
     358       43166 :   C22 = add_slices(m2, p2, V2, 0, m2, 0, p2, M3, 0, m2, 0, p2);
     359       43166 :   if (gc_needed(av, 1))
     360           0 :     gerepileall(av, 4, &C11, &C12, &C21, &C22);  /* destroy V2, M3 */
     361       43166 :   C = shallowconcat(vconcat(C11, C21), vconcat(C12, C22));
     362       43166 :   return gerepilecopy(av, C);
     363             : }
     364             : 
     365             : /* x[i,]*y. Assume lg(x) > 1 and 0 < i < lgcols(x) */
     366             : static GEN
     367   283389479 : ZMrow_ZC_mul_i(GEN x, GEN y, long i, long lx)
     368             : {
     369   283389479 :   pari_sp av = avma;
     370   283389479 :   GEN c = mulii(gcoeff(x,i,1), gel(y,1)), ZERO = gen_0;
     371             :   long k;
     372  3070183095 :   for (k = 2; k < lx; k++)
     373             :   {
     374  2786794028 :     GEN t = mulii(gcoeff(x,i,k), gel(y,k));
     375  2786793704 :     if (t != ZERO) c = addii(c, t);
     376             :   }
     377   283389067 :   return gerepileuptoint(av, c);
     378             : }
     379             : GEN
     380    45657886 : ZMrow_ZC_mul(GEN x, GEN y, long i)
     381    45657886 : { return ZMrow_ZC_mul_i(x, y, i, lg(x)); }
     382             : 
     383             : /* return x * y, 1 < lx = lg(x), l = lgcols(x) */
     384             : static GEN
     385    45557583 : ZM_ZC_mul_i(GEN x, GEN y, long lx, long l)
     386             : {
     387    45557583 :   GEN z = cgetg(l,t_COL);
     388             :   long i;
     389   283289185 :   for (i=1; i<l; i++) gel(z,i) = ZMrow_ZC_mul_i(x,y,i,lx);
     390    45557580 :   return z;
     391             : }
     392             : 
     393             : static GEN
     394    10579227 : ZM_mul_classical(GEN x, GEN y, long l, long lx, long ly)
     395             : {
     396             :   long j;
     397    10579227 :   GEN z = cgetg(ly, t_MAT);
     398    39591187 :   for (j = 1; j < ly; j++)
     399    29011960 :     gel(z, j) = ZM_ZC_mul_i(x, gel(y, j), lx, l);
     400    10579227 :   return z;
     401             : }
     402             : 
     403             : static GEN
     404        1477 : ZM_mul_slice(GEN A, GEN B, GEN P, GEN *mod)
     405             : {
     406        1477 :   pari_sp av = avma;
     407        1477 :   long i, n = lg(P)-1;
     408             :   GEN H, T;
     409        1477 :   if (n == 1)
     410             :   {
     411         847 :     ulong p = uel(P,1);
     412         847 :     GEN a = ZM_to_Flm(A, p);
     413         847 :     GEN b = ZM_to_Flm(B, p);
     414         847 :     GEN Hp = gerepileupto(av, Flm_to_ZM(Flm_mul(a, b, p)));
     415         847 :     *mod = utoi(p); return Hp;
     416             :   }
     417         630 :   T = ZV_producttree(P);
     418         630 :   A = ZM_nv_mod_tree(A, P, T);
     419         630 :   B = ZM_nv_mod_tree(B, P, T);
     420         630 :   H = cgetg(n+1, t_VEC);
     421       10426 :   for(i=1; i <= n; i++)
     422        9796 :     gel(H,i) = Flm_mul(gel(A,i),gel(B,i),P[i]);
     423         630 :   H = nmV_chinese_center_tree_seq(H, P, T, ZV_chinesetree(P, T));
     424         630 :   *mod = gmael(T, lg(T)-1, 1);
     425         630 :   gerepileall(av, 2, &H, mod);
     426         630 :   return H;
     427             : }
     428             : 
     429             : GEN
     430        1477 : ZM_mul_worker(GEN P, GEN A, GEN B)
     431             : {
     432        1477 :   GEN V = cgetg(3, t_VEC);
     433        1477 :   gel(V,1) = ZM_mul_slice(A, B, P, &gel(V,2));
     434        1477 :   return V;
     435             : }
     436             : 
     437             : static long
     438      189519 : ZC_expi(GEN x)
     439             : {
     440      189519 :   long i, l = lg(x), m=-1;
     441    19643480 :   for(i = 1; i < l; i++)
     442             :   {
     443    19453961 :     GEN xi = gel(x,i);
     444    19453961 :     long e = signe(xi) ? expi(xi): -1;
     445    19453961 :     if (e > m) m = e;
     446             :   }
     447      189519 :   return m;
     448             : }
     449             : 
     450             : static long
     451        1950 : ZM_expi(GEN x)
     452             : {
     453        1950 :   long i, l = lg(x), m=-1;
     454      191469 :   for(i = 1; i < l; i++)
     455             :   {
     456      189519 :     long e = ZC_expi(gel(x,i));
     457      189519 :     if (e > m) m = e;
     458             :   }
     459        1950 :   return m;
     460             : }
     461             : 
     462             : static GEN
     463         975 : ZM_mul_fast(GEN A, GEN B, long lA, long lB)
     464             : {
     465         975 :   pari_sp av = avma;
     466             :   forprime_t S;
     467             :   GEN H, worker;
     468         975 :   long h, mA = ZM_expi(A), mB = ZM_expi(B);
     469         975 :   if (mA==-1 || mB==-1) return zeromat(nbrows(A),lB-1);
     470         961 :   h = 3 + mA + mB + expu(lA-1);
     471         961 :   init_modular_big(&S);
     472         961 :   worker = snm_closure(is_entry("_ZM_mul_worker"), mkvec2(A,B));
     473         961 :   H = gen_crt("ZM_mul", worker, &S, NULL, h, 0, NULL,
     474             :               nmV_chinese_center, FpM_center);
     475         961 :   return gerepileupto(av, H);
     476             : }
     477             : 
     478             : /* Strassen-Winograd used for dim >= ZM_sw_bound */
     479             : static GEN
     480    10618202 : ZM_mul_i(GEN x, GEN y, long l, long lx, long ly)
     481             : {
     482    10618202 :   if (lx > 70 && ly > 70 && l > 70)
     483         975 :     return ZM_mul_fast(x, y, lx, ly);
     484             :   else
     485             :   {
     486    10617227 :     long s = maxss(ZM_max_lg_i(x,lx,l), ZM_max_lg_i(y,ly,lx));
     487    10617227 :     long ZM_sw_bound = s > 60 ? 2: s > 25 ? 4: s>15 ? 8 : s > 8 ? 16 : 32;
     488    10617227 :     if (l <= ZM_sw_bound || lx <= ZM_sw_bound || ly <= ZM_sw_bound)
     489    10574070 :       return ZM_mul_classical(x, y, l, lx, ly);
     490             :     else
     491       43157 :       return ZM_mul_sw(x, y, l - 1, lx - 1, ly - 1);
     492             :   }
     493             : }
     494             : 
     495             : GEN
     496    10337146 : ZM_mul(GEN x, GEN y)
     497             : {
     498    10337146 :   long lx=lg(x), ly=lg(y);
     499    10337146 :   if (ly==1) return cgetg(1,t_MAT);
     500    10317076 :   if (lx==1) return zeromat(0, ly-1);
     501    10316040 :   return ZM_mul_i(x, y, lgcols(x), lx, ly);
     502             : }
     503             : 
     504             : GEN
     505      139165 : QM_mul(GEN x, GEN y)
     506             : {
     507      139165 :   GEN dx, nx = Q_primitive_part(x, &dx);
     508      139165 :   GEN dy, ny = Q_primitive_part(y, &dy);
     509      139165 :   GEN z = ZM_mul(nx, ny);
     510      139165 :   if (dx || dy)
     511             :   {
     512      139165 :     GEN d = dx ? dy ? gmul(dx, dy): dx : dy;
     513      139165 :     if (!gequal1(d)) z = ZM_Q_mul(z, d);
     514             :   }
     515      139165 :   return z;
     516             : }
     517             : 
     518             : GEN
     519         700 : QM_sqr(GEN x)
     520             : {
     521         700 :   GEN dx, nx = Q_primitive_part(x, &dx);
     522         700 :   GEN z = ZM_sqr(nx);
     523         700 :   if (dx)
     524         700 :     z = ZM_Q_mul(z, gsqr(dx));
     525         700 :   return z;
     526             : }
     527             : 
     528             : GEN
     529      116899 : QM_QC_mul(GEN x, GEN y)
     530             : {
     531      116899 :   GEN dx, nx = Q_primitive_part(x, &dx);
     532      116899 :   GEN dy, ny = Q_primitive_part(y, &dy);
     533      116899 :   GEN z = ZM_ZC_mul(nx, ny);
     534      116899 :   if (dx || dy)
     535             :   {
     536      116871 :     GEN d = dx ? dy ? gmul(dx, dy): dx : dy;
     537      116871 :     if (!gequal1(d)) z = ZC_Q_mul(z, d);
     538             :   }
     539      116899 :   return z;
     540             : }
     541             : 
     542             : /* assume result is symmetric */
     543             : GEN
     544           0 : ZM_multosym(GEN x, GEN y)
     545             : {
     546           0 :   long j, lx, ly = lg(y);
     547             :   GEN M;
     548           0 :   if (ly == 1) return cgetg(1,t_MAT);
     549           0 :   lx = lg(x); /* = lgcols(y) */
     550           0 :   if (lx == 1) return cgetg(1,t_MAT);
     551             :   /* ly = lgcols(x) */
     552           0 :   M = cgetg(ly, t_MAT);
     553           0 :   for (j=1; j<ly; j++)
     554             :   {
     555           0 :     GEN z = cgetg(ly,t_COL), yj = gel(y,j);
     556             :     long i;
     557           0 :     for (i=1; i<j; i++) gel(z,i) = gcoeff(M,j,i);
     558           0 :     for (i=j; i<ly; i++)gel(z,i) = ZMrow_ZC_mul_i(x,yj,i,lx);
     559           0 :     gel(M,j) = z;
     560             :   }
     561           0 :   return M;
     562             : }
     563             : 
     564             : /* compute m*diagonal(d), assume lg(d) = lg(m). Shallow */
     565             : GEN
     566          14 : ZM_mul_diag(GEN m, GEN d)
     567             : {
     568             :   long j, l;
     569          14 :   GEN y = cgetg_copy(m, &l);
     570          56 :   for (j=1; j<l; j++)
     571             :   {
     572          42 :     GEN c = gel(d,j);
     573          42 :     gel(y,j) = equali1(c)? gel(m,j): ZC_Z_mul(gel(m,j), c);
     574             :   }
     575          14 :   return y;
     576             : }
     577             : /* compute diagonal(d)*m, assume lg(d) = lg(m~). Shallow */
     578             : GEN
     579      132404 : ZM_diag_mul(GEN d, GEN m)
     580             : {
     581      132404 :   long i, j, l = lg(d), lm = lg(m);
     582      132404 :   GEN y = cgetg(lm, t_MAT);
     583      402261 :   for (j=1; j<lm; j++) gel(y,j) = cgetg(l, t_COL);
     584      420160 :   for (i=1; i<l; i++)
     585             :   {
     586      287756 :     GEN c = gel(d,i);
     587      287756 :     if (equali1(c))
     588       33698 :       for (j=1; j<lm; j++) gcoeff(y,i,j) = gcoeff(m,i,j);
     589             :     else
     590      971182 :       for (j=1; j<lm; j++) gcoeff(y,i,j) = mulii(gcoeff(m,i,j), c);
     591             :   }
     592      132404 :   return y;
     593             : }
     594             : 
     595             : /* assume lx > 1 is lg(x) = lg(y) */
     596             : static GEN
     597    14714090 : ZV_dotproduct_i(GEN x, GEN y, long lx)
     598             : {
     599    14714090 :   pari_sp av = avma;
     600    14714090 :   GEN c = mulii(gel(x,1), gel(y,1));
     601             :   long i;
     602   114186383 :   for (i = 2; i < lx; i++)
     603             :   {
     604    99472293 :     GEN t = mulii(gel(x,i), gel(y,i));
     605    99472293 :     if (t != gen_0) c = addii(c, t);
     606             :   }
     607    14714090 :   return gerepileuptoint(av, c);
     608             : }
     609             : 
     610             : /* x~ * y, assuming result is symmetric */
     611             : GEN
     612      161726 : ZM_transmultosym(GEN x, GEN y)
     613             : {
     614      161726 :   long i, j, l, ly = lg(y);
     615             :   GEN M;
     616      161726 :   if (ly == 1) return cgetg(1,t_MAT);
     617             :   /* lg(x) = ly */
     618      161726 :   l = lgcols(y); /* = lgcols(x) */
     619      161726 :   M = cgetg(ly, t_MAT);
     620     1266212 :   for (i=1; i<ly; i++)
     621             :   {
     622     1104486 :     GEN xi = gel(x,i), c = cgetg(ly,t_COL);
     623     1104486 :     gel(M,i) = c;
     624     5066630 :     for (j=1; j<i; j++)
     625     3962144 :       gcoeff(M,i,j) = gel(c,j) = ZV_dotproduct_i(xi,gel(y,j),l);
     626     1104486 :     gel(c,i) = ZV_dotproduct_i(xi,gel(y,i),l);
     627             :   }
     628      161726 :   return M;
     629             : }
     630             : /* x~ * y */
     631             : GEN
     632         497 : ZM_transmul(GEN x, GEN y)
     633             : {
     634         497 :   long i, j, l, lx, ly = lg(y);
     635             :   GEN M;
     636         497 :   if (ly == 1) return cgetg(1,t_MAT);
     637         497 :   lx = lg(x);
     638         497 :   l = lgcols(y);
     639         497 :   if (lgcols(x) != l) pari_err_OP("operation 'ZM_transmul'", x,y);
     640         497 :   M = cgetg(ly, t_MAT);
     641        1617 :   for (i=1; i<ly; i++)
     642             :   {
     643        1120 :     GEN yi = gel(y,i), c = cgetg(lx,t_COL);
     644        1120 :     gel(M,i) = c;
     645        5061 :     for (j=1; j<lx; j++) gel(c,j) = ZV_dotproduct_i(yi,gel(x,j),l);
     646             :   }
     647         497 :   return M;
     648             : }
     649             : 
     650             : static GEN
     651        5166 : ZM_sqr_i(GEN x, long l)
     652             : {
     653        5166 :   if (l > 70)
     654           0 :     return ZM_mul_fast(x, x, l, l);
     655             :   else
     656             :   {
     657        5166 :     long s = ZM_max_lg_i(x,l,l);
     658        5166 :     long ZM_sw_bound = s > 60 ? 2: s > 25 ? 4: s>15 ? 8 : s > 8 ? 16 : 32;
     659        5166 :     if (l <= ZM_sw_bound)
     660        5157 :       return ZM_mul_classical(x, x, l, l, l);
     661             :     else
     662           9 :       return ZM_mul_sw(x, x, l - 1, l - 1, l - 1);
     663             :   }
     664             : }
     665             : 
     666             : GEN
     667        5166 : ZM_sqr(GEN x)
     668             : {
     669        5166 :   long lx=lg(x);
     670        5166 :   if (lx==1) return cgetg(1,t_MAT);
     671        5166 :   return ZM_sqr_i(x, lx);
     672             : }
     673             : GEN
     674    16545623 : ZM_ZC_mul(GEN x, GEN y)
     675             : {
     676    16545623 :   long lx = lg(x);
     677    16545623 :   return lx==1? cgetg(1,t_COL): ZM_ZC_mul_i(x, y, lx, lgcols(x));
     678             : }
     679             : 
     680             : GEN
     681     1294162 : ZC_Z_div(GEN x, GEN c)
     682     6304667 : { pari_APPLY_type(t_COL, Qdivii(gel(x,i), c)) }
     683             : 
     684             : GEN
     685       14546 : ZM_Z_div(GEN x, GEN c)
     686      164262 : { pari_APPLY_same(ZC_Z_div(gel(x, i), c)) }
     687             : 
     688             : GEN
     689     1315201 : ZC_Q_mul(GEN A, GEN z)
     690             : {
     691     1315201 :   pari_sp av = avma;
     692     1315201 :   long i, l = lg(A);
     693             :   GEN d, n, Ad, B, u;
     694     1315201 :   if (typ(z)==t_INT) return ZC_Z_mul(A,z);
     695     1312492 :   n = gel(z, 1); d = gel(z, 2);
     696     1312492 :   Ad = FpC_red(A, d);
     697     1312492 :   u = gcdii(d, FpV_factorback(Ad, NULL, d));
     698     1312492 :   B = cgetg(l, t_COL);
     699     1312492 :   if (equali1(u))
     700             :   {
     701      249355 :     for(i=1; i<l; i++)
     702      205945 :       gel(B, i) = mkfrac(mulii(n, gel(A,i)), d);
     703             :   } else
     704             :   {
     705    12234250 :     for(i=1; i<l; i++)
     706             :     {
     707    10965168 :       GEN di = gcdii(gel(Ad, i), u), ni = mulii(n, diviiexact(gel(A,i), di));
     708    10965168 :       if (equalii(d, di))
     709     7861490 :         gel(B, i) = ni;
     710             :       else
     711     3103678 :         gel(B, i) = mkfrac(ni, diviiexact(d, di));
     712             :     }
     713             :   }
     714     1312492 :   return gerepilecopy(av, B);
     715             : }
     716             : 
     717             : GEN
     718      280765 : ZM_Q_mul(GEN x, GEN z)
     719             : {
     720      280765 :   if (typ(z)==t_INT) return ZM_Z_mul(x,z);
     721     1394349 :   pari_APPLY_same(ZC_Q_mul(gel(x, i), z));
     722             : }
     723             : 
     724             : long
     725   182226205 : zv_dotproduct(GEN x, GEN y)
     726             : {
     727   182226205 :   long i, lx = lg(x);
     728             :   ulong c;
     729   182226205 :   if (lx == 1) return 0;
     730   182226205 :   c = uel(x,1)*uel(y,1);
     731  2825755660 :   for (i = 2; i < lx; i++)
     732  2643529455 :     c += uel(x,i)*uel(y,i);
     733   182226205 :   return c;
     734             : }
     735             : 
     736             : GEN
     737      166031 : ZV_ZM_mul(GEN x, GEN y)
     738             : {
     739      166031 :   long i, lx = lg(x), ly = lg(y);
     740             :   GEN z;
     741      166031 :   if (lx == 1) return zerovec(ly-1);
     742      166017 :   z = cgetg(ly, t_VEC);
     743      624961 :   for (i = 1; i < ly; i++) gel(z,i) = ZV_dotproduct_i(x, gel(y,i), lx);
     744      166017 :   return z;
     745             : }
     746             : 
     747             : GEN
     748           0 : ZC_ZV_mul(GEN x, GEN y)
     749             : {
     750           0 :   long i,j, lx=lg(x), ly=lg(y);
     751             :   GEN z;
     752           0 :   if (ly==1) return cgetg(1,t_MAT);
     753           0 :   z = cgetg(ly,t_MAT);
     754           0 :   for (j=1; j < ly; j++)
     755             :   {
     756           0 :     gel(z,j) = cgetg(lx,t_COL);
     757           0 :     for (i=1; i<lx; i++) gcoeff(z,i,j) = mulii(gel(x,i),gel(y,j));
     758             :   }
     759           0 :   return z;
     760             : }
     761             : 
     762             : GEN
     763     6930718 : ZV_dotsquare(GEN x)
     764             : {
     765             :   long i, lx;
     766             :   pari_sp av;
     767             :   GEN z;
     768     6930718 :   lx = lg(x);
     769     6930718 :   if (lx == 1) return gen_0;
     770     6930718 :   av = avma; z = sqri(gel(x,1));
     771    23511638 :   for (i=2; i<lx; i++) z = addii(z, sqri(gel(x,i)));
     772     6930718 :   return gerepileuptoint(av,z);
     773             : }
     774             : 
     775             : GEN
     776    13463937 : ZV_dotproduct(GEN x,GEN y)
     777             : {
     778             :   long lx;
     779    13463937 :   if (x == y) return ZV_dotsquare(x);
     780     9184575 :   lx = lg(x);
     781     9184575 :   if (lx == 1) return gen_0;
     782     9184575 :   return ZV_dotproduct_i(x, y, lx);
     783             : }
     784             : 
     785             : static GEN
     786         294 : _ZM_mul(void *data /*ignored*/, GEN x, GEN y)
     787         294 : { (void)data; return ZM_mul(x,y); }
     788             : static GEN
     789        3227 : _ZM_sqr(void *data /*ignored*/, GEN x)
     790        3227 : { (void)data; return ZM_sqr(x); }
     791             : GEN
     792           0 : ZM_pow(GEN x, GEN n)
     793             : {
     794           0 :   pari_sp av = avma;
     795           0 :   if (!signe(n)) return matid(lg(x)-1);
     796           0 :   return gerepilecopy(av, gen_pow_i(x, n, NULL, &_ZM_sqr, &_ZM_mul));
     797             : }
     798             : GEN
     799        2779 : ZM_powu(GEN x, ulong n)
     800             : {
     801        2779 :   pari_sp av = avma;
     802        2779 :   if (!n) return matid(lg(x)-1);
     803        2779 :   return gerepilecopy(av, gen_powu_i(x, n, NULL, &_ZM_sqr, &_ZM_mul));
     804             : }
     805             : /********************************************************************/
     806             : /**                                                                **/
     807             : /**                           ADD, SUB                             **/
     808             : /**                                                                **/
     809             : /********************************************************************/
     810             : static GEN
     811    17627503 : ZC_add_i(GEN x, GEN y, long lx)
     812             : {
     813    17627503 :   GEN A = cgetg(lx, t_COL);
     814             :   long i;
     815   341784726 :   for (i=1; i<lx; i++) gel(A,i) = addii(gel(x,i), gel(y,i));
     816    17626962 :   return A;
     817             : }
     818             : GEN
     819     9510650 : ZC_add(GEN x, GEN y) { return ZC_add_i(x, y, lg(x)); }
     820             : GEN
     821      146937 : ZC_Z_add(GEN x, GEN y)
     822             : {
     823      146937 :   long k, lx = lg(x);
     824      146937 :   GEN z = cgetg(lx, t_COL);
     825      146937 :   if (lx == 1) pari_err_TYPE2("+",x,y);
     826      146937 :   gel(z,1) = addii(y,gel(x,1));
     827      833210 :   for (k = 2; k < lx; k++) gel(z,k) = icopy(gel(x,k));
     828      146937 :   return z;
     829             : }
     830             : 
     831             : static GEN
     832     3055762 : ZC_sub_i(GEN x, GEN y, long lx)
     833             : {
     834             :   long i;
     835     3055762 :   GEN A = cgetg(lx, t_COL);
     836    33024092 :   for (i=1; i<lx; i++) gel(A,i) = subii(gel(x,i), gel(y,i));
     837     3055760 :   return A;
     838             : }
     839             : GEN
     840     2927248 : ZC_sub(GEN x, GEN y) { return ZC_sub_i(x, y, lg(x)); }
     841             : GEN
     842           0 : ZC_Z_sub(GEN x, GEN y)
     843             : {
     844           0 :   long k, lx = lg(x);
     845           0 :   GEN z = cgetg(lx, t_COL);
     846           0 :   if (lx == 1) pari_err_TYPE2("+",x,y);
     847           0 :   gel(z,1) = subii(gel(x,1), y);
     848           0 :   for (k = 2; k < lx; k++) gel(z,k) = icopy(gel(x,k));
     849           0 :   return z;
     850             : }
     851             : GEN
     852      111895 : Z_ZC_sub(GEN a, GEN x)
     853             : {
     854      111895 :   long k, lx = lg(x);
     855      111895 :   GEN z = cgetg(lx, t_COL);
     856      111905 :   if (lx == 1) pari_err_TYPE2("-",a,x);
     857      111905 :   gel(z,1) = subii(a, gel(x,1));
     858      308425 :   for (k = 2; k < lx; k++) gel(z,k) = negi(gel(x,k));
     859      111924 :   return z;
     860             : }
     861             : 
     862             : GEN
     863      634873 : ZM_add(GEN x, GEN y)
     864             : {
     865      634873 :   long lx = lg(x), l, j;
     866             :   GEN z;
     867      634873 :   if (lx == 1) return cgetg(1, t_MAT);
     868      626683 :   z = cgetg(lx, t_MAT); l = lgcols(x);
     869     8743207 :   for (j = 1; j < lx; j++) gel(z,j) = ZC_add_i(gel(x,j), gel(y,j), l);
     870      626683 :   return z;
     871             : }
     872             : GEN
     873       58695 : ZM_sub(GEN x, GEN y)
     874             : {
     875       58695 :   long lx = lg(x), l, j;
     876             :   GEN z;
     877       58695 :   if (lx == 1) return cgetg(1, t_MAT);
     878       58695 :   z = cgetg(lx, t_MAT); l = lgcols(x);
     879      187209 :   for (j = 1; j < lx; j++) gel(z,j) = ZC_sub_i(gel(x,j), gel(y,j), l);
     880       58695 :   return z;
     881             : }
     882             : /********************************************************************/
     883             : /**                                                                **/
     884             : /**                         LINEAR COMBINATION                     **/
     885             : /**                                                                **/
     886             : /********************************************************************/
     887             : /* return X/c assuming division is exact */
     888             : GEN
     889     3709829 : ZC_Z_divexact(GEN x, GEN c)
     890    48908862 : { pari_APPLY_type(t_COL, diviiexact(gel(x,i), c)) }
     891             : GEN
     892           0 : ZC_u_divexact(GEN x, ulong c)
     893           0 : { pari_APPLY_type(t_COL, diviuexact(gel(x,i), c)) }
     894             : 
     895             : GEN
     896      656231 : ZM_Z_divexact(GEN x, GEN c)
     897     3734493 : { pari_APPLY_same(ZC_Z_divexact(gel(x,i), c)) }
     898             : 
     899             : GEN
     900    16567515 : ZC_Z_mul(GEN x, GEN c)
     901             : {
     902    16567515 :   if (!signe(c)) return zerocol(lg(x)-1);
     903    16099704 :   if (is_pm1(c)) return (signe(c) > 0)? ZC_copy(x): ZC_neg(x);
     904   177009222 :   pari_APPLY_type(t_COL, mulii(gel(x,i), c))
     905             : }
     906             : 
     907             : GEN
     908       23453 : ZC_z_mul(GEN x, long c)
     909             : {
     910       23453 :   if (!c) return zerocol(lg(x)-1);
     911       16229 :   if (c == 1) return ZC_copy(x);
     912       11707 :   if (c ==-1) return ZC_neg(x);
     913      131534 :   pari_APPLY_type(t_COL, mulsi(c, gel(x,i)))
     914             : }
     915             : 
     916             : GEN
     917       18040 : zv_z_mul(GEN x, long n)
     918      414785 : { pari_APPLY_long(x[i]*n) }
     919             : 
     920             : /* return a ZM */
     921             : GEN
     922           0 : nm_Z_mul(GEN X, GEN c)
     923             : {
     924           0 :   long i, j, h, l = lg(X), s = signe(c);
     925             :   GEN A;
     926           0 :   if (l == 1) return cgetg(1, t_MAT);
     927           0 :   h = lgcols(X);
     928           0 :   if (!s) return zeromat(h-1, l-1);
     929           0 :   if (is_pm1(c)) {
     930           0 :     if (s > 0) return Flm_to_ZM(X);
     931           0 :     X = Flm_to_ZM(X); ZM_togglesign(X); return X;
     932             :   }
     933           0 :   A = cgetg(l, t_MAT);
     934           0 :   for (j = 1; j < l; j++)
     935             :   {
     936           0 :     GEN a = cgetg(h, t_COL), x = gel(X, j);
     937           0 :     for (i = 1; i < h; i++) gel(a,i) = muliu(c, x[i]);
     938           0 :     gel(A,j) = a;
     939             :   }
     940           0 :   return A;
     941             : }
     942             : GEN
     943     1506102 : ZM_Z_mul(GEN X, GEN c)
     944             : {
     945     1506102 :   long i, j, h, l = lg(X);
     946             :   GEN A;
     947     1506102 :   if (l == 1) return cgetg(1, t_MAT);
     948     1506102 :   h = lgcols(X);
     949     1506102 :   if (!signe(c)) return zeromat(h-1, l-1);
     950     1504597 :   if (is_pm1(c)) return (signe(c) > 0)? ZM_copy(X): ZM_neg(X);
     951     1270105 :   A = cgetg(l, t_MAT);
     952    11134260 :   for (j = 1; j < l; j++)
     953             :   {
     954     9864155 :     GEN a = cgetg(h, t_COL), x = gel(X, j);
     955   181107956 :     for (i = 1; i < h; i++) gel(a,i) = mulii(c, gel(x,i));
     956     9864155 :     gel(A,j) = a;
     957             :   }
     958     1270105 :   return A;
     959             : }
     960             : void
     961    54723005 : ZC_lincomb1_inplace_i(GEN X, GEN Y, GEN v, long n)
     962             : {
     963             :   long i;
     964  1543586777 :   for (i = n; i; i--) gel(X,i) = addmulii_inplace(gel(X,i), gel(Y,i), v);
     965    54714963 : }
     966             : /* X <- X + v Y (elementary col operation) */
     967             : void
     968    48892551 : ZC_lincomb1_inplace(GEN X, GEN Y, GEN v)
     969             : {
     970    48892551 :   if (lgefint(v) != 2) return ZC_lincomb1_inplace_i(X, Y, v, lg(X)-1);
     971             : }
     972             : void
     973    11450002 : Flc_lincomb1_inplace(GEN X, GEN Y, ulong v, ulong q)
     974             : {
     975             :   long i;
     976    11450002 :   if (!v) return; /* v = 0 */
     977   327026516 :   for (i = lg(X)-1; i; i--) X[i] = Fl_add(X[i], Fl_mul(Y[i], v, q), q);
     978             : }
     979             : 
     980             : /* X + v Y, wasteful if (v = 0) */
     981             : static GEN
     982     3704736 : ZC_lincomb1(GEN v, GEN x, GEN y)
     983    50439645 : { pari_APPLY_type(t_COL, addmulii(gel(x,i), gel(y,i), v)) }
     984             : 
     985             : /* -X + vY */
     986             : static GEN
     987      236502 : ZC_lincomb_1(GEN v, GEN x, GEN y)
     988     1595256 : { pari_APPLY_type(t_COL, mulsubii(gel(y,i), v, gel(x,i))) }
     989             : 
     990             : /* X,Y compatible ZV; u,v in Z. Returns A = u*X + v*Y */
     991             : GEN
     992    13617666 : ZC_lincomb(GEN u, GEN v, GEN X, GEN Y)
     993             : {
     994             :   long su, sv;
     995             :   GEN A;
     996             : 
     997    13617666 :   su = signe(u); if (!su) return ZC_Z_mul(Y, v);
     998    13617344 :   sv = signe(v); if (!sv) return ZC_Z_mul(X, u);
     999    13492982 :   if (is_pm1(v))
    1000             :   {
    1001     3592095 :     if (is_pm1(u))
    1002             :     {
    1003     3147634 :       if (su != sv) A = ZC_sub(X, Y);
    1004      814563 :       else          A = ZC_add(X, Y);
    1005     3147634 :       if (su < 0) ZV_togglesign(A); /* in place but was created above */
    1006             :     }
    1007             :     else
    1008             :     {
    1009      444461 :       if (sv > 0) A = ZC_lincomb1 (u, Y, X);
    1010      184899 :       else        A = ZC_lincomb_1(u, Y, X);
    1011             :     }
    1012             :   }
    1013     9900946 :   else if (is_pm1(u))
    1014             :   {
    1015     3496777 :     if (su > 0) A = ZC_lincomb1 (v, X, Y);
    1016       51603 :     else        A = ZC_lincomb_1(v, X, Y);
    1017             :   }
    1018             :   else
    1019             :   { /* not cgetg_copy: x may be a t_VEC */
    1020     6404210 :     long i, lx = lg(X);
    1021     6404210 :     A = cgetg(lx,t_COL);
    1022    30460984 :     for (i=1; i<lx; i++) gel(A,i) = lincombii(u,v,gel(X,i),gel(Y,i));
    1023             :   }
    1024    13492505 :   return A;
    1025             : }
    1026             : 
    1027             : /********************************************************************/
    1028             : /**                                                                **/
    1029             : /**                           CONVERSIONS                          **/
    1030             : /**                                                                **/
    1031             : /********************************************************************/
    1032             : GEN
    1033       74984 : ZV_to_nv(GEN x)
    1034      228781 : { pari_APPLY_ulong(itou(gel(x,i))) }
    1035             : 
    1036             : GEN
    1037       33535 : zm_to_ZM(GEN x)
    1038      389443 : { pari_APPLY_type(t_MAT, zc_to_ZC(gel(x,i))) }
    1039             : 
    1040             : GEN
    1041         112 : zmV_to_ZMV(GEN x)
    1042         812 : { pari_APPLY_type(t_VEC, zm_to_ZM(gel(x,i))) }
    1043             : 
    1044             : /* same as Flm_to_ZM but do not assume positivity */
    1045             : GEN
    1046        1218 : ZM_to_zm(GEN x)
    1047      257922 : { pari_APPLY_same(ZV_to_zv(gel(x,i))) }
    1048             : 
    1049             : GEN
    1050      284480 : zv_to_Flv(GEN x, ulong p)
    1051     4486328 : { pari_APPLY_ulong(umodsu(x[i], p)) }
    1052             : 
    1053             : GEN
    1054       19894 : zm_to_Flm(GEN x, ulong p)
    1055      304374 : { pari_APPLY_same(zv_to_Flv(gel(x,i),p)) }
    1056             : 
    1057             : GEN
    1058          35 : ZMV_to_zmV(GEN x)
    1059         357 : { pari_APPLY_type(t_VEC, ZM_to_zm(gel(x,i))) }
    1060             : 
    1061             : /********************************************************************/
    1062             : /**                                                                **/
    1063             : /**                         COPY, NEGATION                         **/
    1064             : /**                                                                **/
    1065             : /********************************************************************/
    1066             : GEN
    1067     5246333 : ZC_copy(GEN x)
    1068             : {
    1069     5246333 :   long i, lx = lg(x);
    1070     5246333 :   GEN y = cgetg(lx, t_COL);
    1071    45589792 :   for (i=1; i<lx; i++)
    1072             :   {
    1073    40343427 :     GEN c = gel(x,i);
    1074    40343427 :     gel(y,i) = lgefint(c) == 2? gen_0: icopy(c);
    1075             :   }
    1076     5246365 :   return y;
    1077             : }
    1078             : 
    1079             : GEN
    1080      287152 : ZM_copy(GEN x)
    1081     2799112 : { pari_APPLY_same(ZC_copy(gel(x,i))) }
    1082             : 
    1083             : void
    1084       67505 : ZV_neg_inplace(GEN M)
    1085             : {
    1086       67505 :   long l = lg(M);
    1087      283779 :   while (--l > 0) gel(M,l) = negi(gel(M,l));
    1088       67505 : }
    1089             : GEN
    1090     1813893 : ZC_neg(GEN x)
    1091    16510808 : { pari_APPLY_type(t_COL, negi(gel(x,i))) }
    1092             : GEN
    1093       31644 : zv_neg(GEN x)
    1094      369873 : { pari_APPLY_long(-x[i]) }
    1095             : GEN
    1096         261 : zv_neg_inplace(GEN M)
    1097             : {
    1098         261 :   long l = lg(M);
    1099        1079 :   while (--l > 0) M[l] = -M[l];
    1100         261 :   return M;
    1101             : }
    1102             : GEN
    1103      303859 : ZM_neg(GEN x)
    1104     1250250 : { pari_APPLY_same(ZC_neg(gel(x,i))) }
    1105             : 
    1106             : void
    1107     1687353 : ZV_togglesign(GEN M)
    1108             : {
    1109     1687353 :   long l = lg(M);
    1110    41111168 :   while (--l > 0) togglesign_safe(&gel(M,l));
    1111     1687353 : }
    1112             : void
    1113           0 : ZM_togglesign(GEN M)
    1114             : {
    1115           0 :   long l = lg(M);
    1116           0 :   while (--l > 0) ZV_togglesign(gel(M,l));
    1117           0 : }
    1118             : 
    1119             : /********************************************************************/
    1120             : /**                                                                **/
    1121             : /**                        "DIVISION" mod HNF                      **/
    1122             : /**                                                                **/
    1123             : /********************************************************************/
    1124             : /* Reduce ZC x modulo ZM y in HNF, may return x itself (not a copy) */
    1125             : GEN
    1126     1958817 : ZC_hnfremdiv(GEN x, GEN y, GEN *Q)
    1127             : {
    1128     1958817 :   long i, l = lg(x);
    1129             :   GEN q;
    1130             : 
    1131     1958817 :   if (Q) *Q = cgetg(l,t_COL);
    1132     1958817 :   if (l == 1) return cgetg(1,t_COL);
    1133    17488659 :   for (i = l-1; i>0; i--)
    1134             :   {
    1135    15529842 :     q = diviiround(gel(x,i), gcoeff(y,i,i));
    1136    15529842 :     if (signe(q)) {
    1137     6184289 :       togglesign(q);
    1138     6184289 :       x = ZC_lincomb(gen_1, q, x, gel(y,i));
    1139             :     }
    1140    15529842 :     if (Q) gel(*Q, i) = q;
    1141             :   }
    1142     1958817 :   return x;
    1143             : }
    1144             : 
    1145             : /* x = y Q + R, may return some columns of x (not copies) */
    1146             : GEN
    1147       65883 : ZM_hnfdivrem(GEN x, GEN y, GEN *Q)
    1148             : {
    1149       65883 :   long lx = lg(x), i;
    1150       65883 :   GEN R = cgetg(lx, t_MAT);
    1151       65883 :   if (Q)
    1152             :   {
    1153       23856 :     GEN q = cgetg(lx, t_MAT); *Q = q;
    1154       44958 :     for (i=1; i<lx; i++) gel(R,i) = ZC_hnfremdiv(gel(x,i),y,(GEN*)(q+i));
    1155             :   }
    1156             :   else
    1157      109708 :     for (i=1; i<lx; i++)
    1158             :     {
    1159       67681 :       pari_sp av = avma;
    1160       67681 :       GEN z = ZC_hnfrem(gel(x,i),y);
    1161       67681 :       gel(R,i) = (avma == av)? ZC_copy(z): gerepileupto(av, z);
    1162             :     }
    1163       65883 :   return R;
    1164             : }
    1165             : 
    1166             : /********************************************************************/
    1167             : /**                                                                **/
    1168             : /**                               TESTS                            **/
    1169             : /**                                                                **/
    1170             : /********************************************************************/
    1171             : int
    1172    21184744 : zv_equal0(GEN V)
    1173             : {
    1174    21184744 :   long l = lg(V);
    1175    34377022 :   while (--l > 0)
    1176    28051218 :     if (V[l]) return 0;
    1177     6325804 :   return 1;
    1178             : }
    1179             : 
    1180             : int
    1181     1278263 : ZV_equal0(GEN V)
    1182             : {
    1183     1278263 :   long l = lg(V);
    1184     4646251 :   while (--l > 0)
    1185     4245234 :     if (signe(gel(V,l))) return 0;
    1186      401017 :   return 1;
    1187             : }
    1188             : int
    1189       13806 : ZMrow_equal0(GEN V, long i)
    1190             : {
    1191       13806 :   long l = lg(V);
    1192       22979 :   while (--l > 0)
    1193       19247 :     if (signe(gcoeff(V,i,l))) return 0;
    1194        3732 :   return 1;
    1195             : }
    1196             : 
    1197             : static int
    1198    15864528 : ZV_equal_lg(GEN V, GEN W, long l)
    1199             : {
    1200    24136820 :   while (--l > 0)
    1201    22127920 :     if (!equalii(gel(V,l), gel(W,l))) return 0;
    1202     2008900 :   return 1;
    1203             : }
    1204             : int
    1205    14624336 : ZV_equal(GEN V, GEN W)
    1206             : {
    1207    14624336 :   long l = lg(V);
    1208    14624336 :   if (lg(W) != l) return 0;
    1209    14624336 :   return ZV_equal_lg(V, W, l);
    1210             : }
    1211             : int
    1212      478595 : ZM_equal(GEN A, GEN B)
    1213             : {
    1214      478595 :   long i, m, l = lg(A);
    1215      478595 :   if (lg(B) != l) return 0;
    1216      478595 :   if (l == 1) return 1;
    1217      478595 :   m = lgcols(A);
    1218      478595 :   if (lgcols(B) != m) return 0;
    1219     1651470 :   for (i = 1; i < l; i++)
    1220     1240188 :     if (!ZV_equal_lg(gel(A,i), gel(B,i), m)) return 0;
    1221      411282 :   return 1;
    1222             : }
    1223             : int
    1224       68590 : ZM_equal0(GEN A)
    1225             : {
    1226       68590 :   long i, j, m, l = lg(A);
    1227       68590 :   if (l == 1) return 1;
    1228       68590 :   m = lgcols(A);
    1229      114783 :   for (j = 1; j < l; j++)
    1230     2465180 :     for (i = 1; i < m; i++)
    1231     2418987 :       if (signe(gcoeff(A,i,j))) return 0;
    1232       13922 :   return 1;
    1233             : }
    1234             : int
    1235    30122377 : zv_equal(GEN V, GEN W)
    1236             : {
    1237    30122377 :   long l = lg(V);
    1238    30122377 :   if (lg(W) != l) return 0;
    1239   254087978 :   while (--l > 0)
    1240   224482297 :     if (V[l] != W[l]) return 0;
    1241    29605681 :   return 1;
    1242             : }
    1243             : 
    1244             : int
    1245        1337 : zvV_equal(GEN V, GEN W)
    1246             : {
    1247        1337 :   long l = lg(V);
    1248        1337 :   if (lg(W) != l) return 0;
    1249       77707 :   while (--l > 0)
    1250       77301 :     if (!zv_equal(gel(V,l),gel(W,l))) return 0;
    1251         406 :   return 1;
    1252             : }
    1253             : 
    1254             : int
    1255      155879 : ZM_ishnf(GEN x)
    1256             : {
    1257      155879 :   long i,j, lx = lg(x);
    1258      442125 :   for (i=1; i<lx; i++)
    1259             :   {
    1260      328141 :     GEN xii = gcoeff(x,i,i);
    1261      328141 :     if (signe(xii) <= 0) return 0;
    1262      627103 :     for (j=1; j<i; j++)
    1263      324358 :       if (signe(gcoeff(x,i,j))) return 0;
    1264      657784 :     for (j=i+1; j<lx; j++)
    1265             :     {
    1266      371538 :       GEN xij = gcoeff(x,i,j);
    1267      371538 :       if (signe(xij)<0 || cmpii(xij,xii)>=0) return 0;
    1268             :     }
    1269             :   }
    1270      113984 :   return 1;
    1271             : }
    1272             : int
    1273      122841 : ZM_isidentity(GEN x)
    1274             : {
    1275      122841 :   long i,j, lx = lg(x);
    1276             : 
    1277      122841 :   if (lx == 1) return 1;
    1278      122834 :   if (lx != lgcols(x)) return 0;
    1279      883250 :   for (j=1; j<lx; j++)
    1280             :   {
    1281      760465 :     GEN c = gel(x,j);
    1282     3533596 :     for (i=1; i<j; )
    1283     2773152 :       if (signe(gel(c,i++))) return 0;
    1284             :     /* i = j */
    1285      760444 :     if (!equali1(gel(c,i++))) return 0;
    1286     3533610 :     for (   ; i<lx; )
    1287     2773187 :       if (signe(gel(c,i++))) return 0;
    1288             :   }
    1289      122785 :   return 1;
    1290             : }
    1291             : int
    1292       75241 : ZM_isdiagonal(GEN x)
    1293             : {
    1294       75241 :   long i,j, lx = lg(x);
    1295       75241 :   if (lx == 1) return 1;
    1296       75241 :   if (lx != lgcols(x)) return 0;
    1297             : 
    1298      189907 :   for (j=1; j<lx; j++)
    1299             :   {
    1300      156721 :     GEN c = gel(x,j);
    1301      222570 :     for (i=1; i<j; i++)
    1302      107869 :       if (signe(gel(c,i))) return 0;
    1303      280911 :     for (i++; i<lx; i++)
    1304      166245 :       if (signe(gel(c,i))) return 0;
    1305             :   }
    1306       33186 :   return 1;
    1307             : }
    1308             : int
    1309       43816 : ZM_isscalar(GEN x, GEN s)
    1310             : {
    1311       43816 :   long i, j, lx = lg(x);
    1312             : 
    1313       43816 :   if (lx == 1) return 1;
    1314       43816 :   if (!s) s = gcoeff(x,1,1);
    1315       43816 :   if (equali1(s)) return ZM_isidentity(x);
    1316       43816 :   if (lx != lgcols(x)) return 0;
    1317      235642 :   for (j=1; j<lx; j++)
    1318             :   {
    1319      215516 :     GEN c = gel(x,j);
    1320     1521162 :     for (i=1; i<j; )
    1321     1328542 :       if (signe(gel(c,i++))) return 0;
    1322             :     /* i = j */
    1323      192620 :     if (!equalii(gel(c,i++), s)) return 0;
    1324     1534844 :     for (   ; i<lx; )
    1325     1343018 :       if (signe(gel(c,i++))) return 0;
    1326             :   }
    1327       20126 :   return 1;
    1328             : }
    1329             : 
    1330             : long
    1331      102002 : ZC_is_ei(GEN x)
    1332             : {
    1333      102002 :   long i, j = 0, l = lg(x);
    1334      960317 :   for (i = 1; i < l; i++)
    1335             :   {
    1336      858316 :     GEN c = gel(x,i);
    1337      858316 :     long s = signe(c);
    1338      858316 :     if (!s) continue;
    1339      101996 :     if (s < 0 || !is_pm1(c) || j) return 0;
    1340      101995 :     j = i;
    1341             :   }
    1342      102001 :   return j;
    1343             : }
    1344             : 
    1345             : /********************************************************************/
    1346             : /**                                                                **/
    1347             : /**                       MISCELLANEOUS                            **/
    1348             : /**                                                                **/
    1349             : /********************************************************************/
    1350             : /* assume lg(x) = lg(y), x,y in Z^n */
    1351             : int
    1352      822682 : ZV_cmp(GEN x, GEN y)
    1353             : {
    1354      822682 :   long fl,i, lx = lg(x);
    1355     1505491 :   for (i=1; i<lx; i++)
    1356     1218618 :     if (( fl = cmpii(gel(x,i), gel(y,i)) )) return fl;
    1357      286873 :   return 0;
    1358             : }
    1359             : /* assume lg(x) = lg(y), x,y in Z^n */
    1360             : int
    1361        4523 : ZV_abscmp(GEN x, GEN y)
    1362             : {
    1363        4523 :   long fl,i, lx = lg(x);
    1364       21566 :   for (i=1; i<lx; i++)
    1365       21486 :     if (( fl = abscmpii(gel(x,i), gel(y,i)) )) return fl;
    1366          80 :   return 0;
    1367             : }
    1368             : 
    1369             : long
    1370     3631047 : zv_content(GEN x)
    1371             : {
    1372     3631047 :   long i, s, l = lg(x);
    1373     3631047 :   if (l == 1) return 0;
    1374     3631040 :   s = labs(x[1]);
    1375     6371724 :   for (i = 2; i < l && s != 1; i++) s = ugcd(s, labs(x[i]));
    1376     3631040 :   return s;
    1377             : }
    1378             : GEN
    1379      303050 : ZV_content(GEN x)
    1380             : {
    1381      303050 :   long i, l = lg(x);
    1382      303050 :   pari_sp av = avma;
    1383             :   GEN c;
    1384      303050 :   if (l == 1) return gen_0;
    1385      303050 :   if (l == 2) return absi(gel(x,1));
    1386      234982 :   c = gel(x,1);
    1387      553760 :   for (i = 2; i < l; i++)
    1388             :   {
    1389      405610 :     c = gcdii(c, gel(x,i));
    1390      405610 :     if (is_pm1(c)) { set_avma(av); return gen_1; }
    1391             :   }
    1392      148150 :   return gerepileuptoint(av, c);
    1393             : }
    1394             : 
    1395             : GEN
    1396     1088131 : ZM_det_triangular(GEN mat)
    1397             : {
    1398             :   pari_sp av;
    1399     1088131 :   long i,l = lg(mat);
    1400             :   GEN s;
    1401             : 
    1402     1088131 :   if (l<3) return l<2? gen_1: icopy(gcoeff(mat,1,1));
    1403      980222 :   av = avma; s = gcoeff(mat,1,1);
    1404     3150561 :   for (i=2; i<l; i++) s = mulii(s,gcoeff(mat,i,i));
    1405      980222 :   return gerepileuptoint(av,s);
    1406             : }
    1407             : 
    1408             : /* assumes no overflow */
    1409             : long
    1410      394611 : zv_prod(GEN v)
    1411             : {
    1412      394611 :   long n, i, l = lg(v);
    1413      394611 :   if (l == 1) return 1;
    1414      434070 :   n = v[1]; for (i = 2; i < l; i++) n *= v[i];
    1415      295652 :   return n;
    1416             : }
    1417             : 
    1418             : static GEN
    1419   278325396 : _mulii(void *E, GEN a, GEN b)
    1420   278325396 : { (void) E; return mulii(a, b); }
    1421             : 
    1422             : /* product of ulongs */
    1423             : GEN
    1424      375270 : zv_prod_Z(GEN v)
    1425             : {
    1426      375270 :   pari_sp av = avma;
    1427      375270 :   long k, m, n = lg(v)-1;
    1428             :   GEN V;
    1429      375270 :   switch(n) {
    1430       19607 :     case 0: return gen_1;
    1431       59787 :     case 1: return utoi(v[1]);
    1432      130683 :     case 2: return muluu(v[1], v[2]);
    1433             :   }
    1434      165193 :   m = n >> 1;
    1435      165193 :   V = cgetg(m + (odd(n)? 2: 1), t_VEC);
    1436     6654769 :   for (k = 1; k <= m; k++) gel(V,k) = muluu(v[k<<1], v[(k<<1)-1]);
    1437      165193 :   if (odd(n)) gel(V,k) = utoipos(v[n]);
    1438      165193 :   return gerepileuptoint(av, gen_product(V, NULL, &_mulii));
    1439             : }
    1440             : GEN
    1441    14694393 : vecsmall_prod(GEN v)
    1442             : {
    1443    14694393 :   pari_sp av = avma;
    1444    14694393 :   long k, m, n = lg(v)-1;
    1445             :   GEN V;
    1446    14694393 :   switch (n) {
    1447           0 :     case 0: return gen_1;
    1448           0 :     case 1: return stoi(v[1]);
    1449          21 :     case 2: return mulss(v[1], v[2]);
    1450             :   }
    1451    14694372 :   m = n >> 1;
    1452    14694372 :   V = cgetg(m + (odd(n)? 2: 1), t_VEC);
    1453   161556906 :   for (k = 1; k <= m; k++) gel(V,k) = mulss(v[k<<1], v[(k<<1)-1]);
    1454    14694372 :   if (odd(n)) gel(V,k) = stoi(v[n]);
    1455    14694372 :   return gerepileuptoint(av, gen_product(V, NULL, &_mulii));
    1456             : }
    1457             : 
    1458             : GEN
    1459     7606812 : ZV_prod(GEN v)
    1460             : {
    1461     7606812 :   pari_sp av = avma;
    1462     7606812 :   long i, l = lg(v);
    1463             :   GEN n;
    1464     7606812 :   if (l == 1) return gen_1;
    1465     7584265 :   if (l > 7) return gerepileuptoint(av, gen_product(v, NULL, _mulii));
    1466      242576 :   n = gel(v,1);
    1467      242576 :   if (l == 2) return icopy(n);
    1468      365437 :   for (i = 2; i < l; i++) n = mulii(n, gel(v,i));
    1469      144995 :   return gerepileuptoint(av, n);
    1470             : }
    1471             : /* assumes no overflow */
    1472             : long
    1473        2030 : zv_sum(GEN v)
    1474             : {
    1475        2030 :   long n, i, l = lg(v);
    1476        2030 :   if (l == 1) return 0;
    1477       23898 :   n = v[1]; for (i = 2; i < l; i++) n += v[i];
    1478        2009 :   return n;
    1479             : }
    1480             : GEN
    1481          63 : ZV_sum(GEN v)
    1482             : {
    1483          63 :   pari_sp av = avma;
    1484          63 :   long i, l = lg(v);
    1485             :   GEN n;
    1486          63 :   if (l == 1) return gen_0;
    1487          63 :   n = gel(v,1);
    1488          63 :   if (l == 2) return icopy(n);
    1489         532 :   for (i = 2; i < l; i++) n = addii(n, gel(v,i));
    1490          63 :   return gerepileuptoint(av, n);
    1491             : }
    1492             : 
    1493             : /********************************************************************/
    1494             : /**                                                                **/
    1495             : /**         GRAM SCHMIDT REDUCTION (integer matrices)              **/
    1496             : /**                                                                **/
    1497             : /********************************************************************/
    1498             : 
    1499             : /* L[k,] += q * L[l,], l < k. Inefficient if q = 0 */
    1500             : static void
    1501       74220 : Zupdate_row(long k, long l, GEN q, GEN L, GEN B)
    1502             : {
    1503       74220 :   long i, qq = itos_or_0(q);
    1504       74220 :   if (!qq)
    1505             :   {
    1506        6848 :     for(i=1;i<l;i++)  gcoeff(L,k,i) = addii(gcoeff(L,k,i),mulii(q,gcoeff(L,l,i)));
    1507        1293 :     gcoeff(L,k,l) = addii(gcoeff(L,k,l), mulii(q,B));
    1508        1293 :     return;
    1509             :   }
    1510       72927 :   if (qq == 1) {
    1511       31291 :     for (i=1;i<l; i++) gcoeff(L,k,i) = addii(gcoeff(L,k,i),gcoeff(L,l,i));
    1512       19747 :     gcoeff(L,k,l) = addii(gcoeff(L,k,l), B);
    1513       53180 :   } else if (qq == -1) {
    1514       39474 :     for (i=1;i<l; i++) gcoeff(L,k,i) = subii(gcoeff(L,k,i),gcoeff(L,l,i));
    1515       22430 :     gcoeff(L,k,l) = addii(gcoeff(L,k,l), negi(B));
    1516             :   } else {
    1517       65315 :     for(i=1;i<l;i++) gcoeff(L,k,i) = addii(gcoeff(L,k,i),mulsi(qq,gcoeff(L,l,i)));
    1518       30750 :     gcoeff(L,k,l) = addii(gcoeff(L,k,l), mulsi(qq,B));
    1519             :   }
    1520             : }
    1521             : 
    1522             : /* update L[k,] */
    1523             : static void
    1524      177387 : ZRED(long k, long l, GEN x, GEN L, GEN B)
    1525             : {
    1526      177387 :   GEN q = truedivii(addii(B,shifti(gcoeff(L,k,l),1)), shifti(B,1));
    1527      177387 :   if (!signe(q)) return;
    1528       74220 :   q = negi(q);
    1529       74220 :   Zupdate_row(k,l,q,L,B);
    1530       74220 :   gel(x,k) = ZC_lincomb(gen_1, q, gel(x,k), gel(x,l));
    1531             : }
    1532             : 
    1533             : /* Gram-Schmidt reduction, x a ZM */
    1534             : static void
    1535      201131 : ZincrementalGS(GEN x, GEN L, GEN B, long k)
    1536             : {
    1537             :   long i, j;
    1538      660877 :   for (j=1; j<=k; j++)
    1539             :   {
    1540      459746 :     pari_sp av = avma;
    1541      459746 :     GEN u = ZV_dotproduct(gel(x,k), gel(x,j));
    1542      962619 :     for (i=1; i<j; i++)
    1543             :     {
    1544      502873 :       u = subii(mulii(gel(B,i+1), u), mulii(gcoeff(L,k,i), gcoeff(L,j,i)));
    1545      502873 :       u = diviiexact(u, gel(B,i));
    1546             :     }
    1547      459746 :     gcoeff(L,k,j) = gerepileuptoint(av, u);
    1548             :   }
    1549      201131 :   gel(B,k+1) = gcoeff(L,k,k); gcoeff(L,k,k) = gen_1;
    1550      201131 : }
    1551             : 
    1552             : /* Variant reducemodinvertible(ZC v, ZM y), when y singular.
    1553             :  * Very inefficient if y is not LLL-reduced of maximal rank */
    1554             : static GEN
    1555       21231 : ZC_reducemodmatrix_i(GEN v, GEN y)
    1556             : {
    1557       21231 :   GEN B, L, x = shallowconcat(y, v);
    1558       21231 :   long k, lx = lg(x), nx = lx-1;
    1559             : 
    1560       21231 :   B = scalarcol_shallow(gen_1, lx);
    1561       21231 :   L = zeromatcopy(nx, nx);
    1562       90965 :   for (k=1; k <= nx; k++) ZincrementalGS(x, L, B, k);
    1563       69734 :   for (k = nx-1; k >= 1; k--) ZRED(nx,k, x,L,gel(B,k+1));
    1564       21231 :   return gel(x,nx);
    1565             : }
    1566             : GEN
    1567       21231 : ZC_reducemodmatrix(GEN v, GEN y) {
    1568       21231 :   pari_sp av = avma;
    1569       21231 :   return gerepilecopy(av, ZC_reducemodmatrix_i(v,y));
    1570             : }
    1571             : 
    1572             : /* Variant reducemodinvertible(ZM v, ZM y), when y singular.
    1573             :  * Very inefficient if y is not LLL-reduced of maximal rank */
    1574             : static GEN
    1575       37149 : ZM_reducemodmatrix_i(GEN v, GEN y)
    1576             : {
    1577             :   GEN B, L, V;
    1578       37149 :   long j, k, lv = lg(v), nx = lg(y), lx = nx+1;
    1579             : 
    1580       37149 :   V = cgetg(lv, t_MAT);
    1581       37149 :   B = scalarcol_shallow(gen_1, lx);
    1582       37149 :   L = zeromatcopy(nx, nx);
    1583       97713 :   for (k=1; k < nx; k++) ZincrementalGS(y, L, B, k);
    1584      107982 :   for (j = 1; j < lg(v); j++)
    1585             :   {
    1586       70833 :     GEN x = shallowconcat(y, gel(v,j));
    1587       70833 :     ZincrementalGS(x, L, B, nx); /* overwrite last */
    1588      199717 :     for (k = nx-1; k >= 1; k--) ZRED(nx,k, x,L,gel(B,k+1));
    1589       70833 :     gel(V,j) = gel(x,nx);
    1590             :   }
    1591       37149 :   return V;
    1592             : }
    1593             : GEN
    1594       37149 : ZM_reducemodmatrix(GEN v, GEN y) {
    1595       37149 :   pari_sp av = avma;
    1596       37149 :   return gerepilecopy(av, ZM_reducemodmatrix_i(v,y));
    1597             : }
    1598             : 
    1599             : GEN
    1600       11193 : ZC_reducemodlll(GEN x,GEN y)
    1601             : {
    1602       11193 :   pari_sp av = avma;
    1603       11193 :   GEN z = ZC_reducemodmatrix(x, ZM_lll(y, 0.75, LLL_INPLACE));
    1604       11193 :   return gerepilecopy(av, z);
    1605             : }
    1606             : GEN
    1607           0 : ZM_reducemodlll(GEN x,GEN y)
    1608             : {
    1609           0 :   pari_sp av = avma;
    1610           0 :   GEN z = ZM_reducemodmatrix(x, ZM_lll(y, 0.75, LLL_INPLACE));
    1611           0 :   return gerepilecopy(av, z);
    1612             : }

Generated by: LCOV version 1.13