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.0 lcov report (development 23222-06b1652be) Lines: 765 839 91.2 %
Date: 2018-11-15 05:40:53 Functions: 116 122 95.1 %
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     1363509 : check_ZV(GEN x, long l)
      19             : {
      20             :   long i;
      21     8693001 :   for (i=1; i<l; i++)
      22     7329576 :     if (typ(gel(x,i)) != t_INT) return 0;
      23     1363425 :   return 1;
      24             : }
      25             : void
      26      509105 : RgV_check_ZV(GEN A, const char *s)
      27             : {
      28      509105 :   if (!RgV_is_ZV(A)) pari_err_TYPE(stack_strcat(s," [integer vector]"), A);
      29      509098 : }
      30             : void
      31      395773 : RgM_check_ZM(GEN A, const char *s)
      32             : {
      33      395773 :   long n = lg(A);
      34      395773 :   if (n != 1)
      35             :   {
      36      395717 :     long j, m = lgcols(A);
      37     1759142 :     for (j=1; j<n; j++)
      38     1363509 :       if (!check_ZV(gel(A,j), m))
      39          84 :         pari_err_TYPE(stack_strcat(s," [integer matrix]"), A);
      40             :   }
      41      395689 : }
      42             : 
      43             : static long
      44    53000273 : ZV_max_lg_i(GEN x, long m)
      45             : {
      46    53000273 :   long i, prec = 2;
      47    53000273 :   for (i=1; i<m; i++) { long l = lgefint(gel(x,i)); if (l > prec) prec = l; }
      48    53000273 :   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    17092334 : ZM_max_lg_i(GEN x, long n, long m)
      57             : {
      58    17092334 :   long prec = 2;
      59    17092334 :   if (n != 1)
      60             :   {
      61             :     long j;
      62    70081995 :     for (j=1; j<n; j++)
      63             :     {
      64    52989661 :       long l = ZV_max_lg_i(gel(x,j), m);
      65    52989661 :       if (l > prec) prec = l;
      66             :     }
      67             :   }
      68    17092334 :   return prec;
      69             : }
      70             : 
      71             : long
      72        3846 : ZM_max_lg(GEN x)
      73             : {
      74        3846 :   long n = lg(x);
      75        3846 :   if (n==1) return 2;
      76        3846 :   return ZM_max_lg_i(x, n, lgcols(x));
      77             : }
      78             : 
      79             : GEN
      80        3269 : ZM_supnorm(GEN x)
      81             : {
      82        3269 :   long i, j, h, lx = lg(x);
      83        3269 :   GEN s = gen_0;
      84        3269 :   if (lx == 1) return gen_1;
      85        3269 :   h = lgcols(x);
      86       20265 :   for (j=1; j<lx; j++)
      87             :   {
      88       16996 :     GEN xj = gel(x,j);
      89      241038 :     for (i=1; i<h; i++)
      90             :     {
      91      224042 :       GEN c = gel(xj,i);
      92      224042 :       if (abscmpii(c, s) > 0) s = c;
      93             :     }
      94             :   }
      95        3269 :   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      688911 : ZM_nc_mul_i(GEN x, GEN y, long c, long l)
     106             : {
     107             :   long i, j;
     108             :   pari_sp av;
     109      688911 :   GEN z = cgetg(l,t_COL), s;
     110             : 
     111     3225512 :   for (i=1; i<l; i++)
     112             :   {
     113     2536601 :     av = avma; s = muliu(gcoeff(x,i,1),y[1]);
     114    20280099 :     for (j=2; j<c; j++)
     115    17743498 :       if (y[j]) s = addii(s, muliu(gcoeff(x,i,j),y[j]));
     116     2536601 :     gel(z,i) = gerepileuptoint(av,s);
     117             :   }
     118      688911 :   return z;
     119             : }
     120             : 
     121             : /* x ZV, y a compatible zc. */
     122             : GEN
     123        2065 : ZV_zc_mul(GEN x, GEN y)
     124             : {
     125        2065 :   long j, l = lg(x);
     126        2065 :   pari_sp av = avma;
     127        2065 :   GEN s = mulis(gel(x,1),y[1]);
     128      166775 :   for (j=2; j<l; j++)
     129      164710 :     if (y[j]) s = addii(s, mulis(gel(x,j),y[j]));
     130        2065 :   return gerepileuptoint(av,s);
     131             : }
     132             : 
     133             : /* x non-empty ZM, y a compatible zc (dimension > 0). */
     134             : static GEN
     135     4697940 : ZM_zc_mul_i(GEN x, GEN y, long c, long l)
     136             : {
     137             :   long i, j;
     138     4697940 :   GEN z = cgetg(l,t_COL);
     139             : 
     140    23757583 :   for (i=1; i<l; i++)
     141             :   {
     142    19059643 :     pari_sp av = avma;
     143    19059643 :     GEN s = mulis(gcoeff(x,i,1),y[1]);
     144   294140668 :     for (j=2; j<c; j++)
     145   275081025 :       if (y[j]) s = addii(s, mulis(gcoeff(x,i,j),y[j]));
     146    19059643 :     gel(z,i) = gerepileuptoint(av,s);
     147             :   }
     148     4697940 :   return z;
     149             : }
     150             : GEN
     151     3309268 : ZM_zc_mul(GEN x, GEN y) {
     152     3309268 :   long l = lg(x);
     153     3309268 :   if (l == 1) return cgetg(1, t_COL);
     154     3309268 :   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        1204 : zv_ZM_mul(GEN x, GEN y) {
     160        1204 :   long i,j, lx = lg(x), ly = lg(y);
     161             :   GEN z;
     162        1204 :   if (lx == 1) return zerovec(ly-1);
     163        1204 :   z = cgetg(ly,t_VEC);
     164        2968 :   for (j=1; j<ly; j++)
     165             :   {
     166        1764 :     pari_sp av = avma;
     167        1764 :     GEN s = mulsi(x[1], gcoeff(y,1,j));
     168        3416 :     for (i=2; i<lx; i++)
     169        1652 :       if (x[i]) s = addii(s, mulsi(x[i], gcoeff(y,i,j)));
     170        1764 :     gel(z,j) = gerepileuptoint(av,s);
     171             :   }
     172        1204 :   return z;
     173             : }
     174             : 
     175             : /* x ZM, y a compatible zm (dimension > 0). */
     176             : GEN
     177      678509 : ZM_zm_mul(GEN x, GEN y)
     178             : {
     179      678509 :   long j, c, l = lg(x), ly = lg(y);
     180      678509 :   GEN z = cgetg(ly, t_MAT);
     181      678509 :   if (l == 1) return z;
     182      678509 :   c = lgcols(x);
     183      678509 :   for (j = 1; j < ly; j++) gel(z,j) = ZM_zc_mul_i(x, gel(y,j), l,c);
     184      678509 :   return z;
     185             : }
     186             : /* x ZM, y a compatible zn (dimension > 0). */
     187             : GEN
     188      671839 : ZM_nm_mul(GEN x, GEN y)
     189             : {
     190      671839 :   long j, c, l = lg(x), ly = lg(y);
     191      671839 :   GEN z = cgetg(ly, t_MAT);
     192      671839 :   if (l == 1) return z;
     193      671839 :   c = lgcols(x);
     194      671839 :   for (j = 1; j < ly; j++) gel(z,j) = ZM_nc_mul_i(x, gel(y,j), l,c);
     195      671839 :   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      175216 : 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      175216 :   long min_d = minss(da, db), min_e = minss(ea, eb), i, j;
     210      175216 :   GEN M = cgetg(n + 1, t_MAT), C;
     211             : 
     212     3175904 :   for (j = 1; j <= min_e; j++) {
     213     3000688 :     gel(M, j) = C = cgetg(m + 1, t_COL);
     214    67250563 :     for (i = 1; i <= min_d; i++)
     215   128499750 :       gel(C, i) = addii(gcoeff(A, ma + i, na + j),
     216    64249875 :                         gcoeff(B, mb + i, nb + j));
     217     3132273 :     for (; i <= da; i++)
     218      131585 :       gel(C, i) = gcoeff(A, ma + i, na + j);
     219     3000688 :     for (; i <= db; i++)
     220           0 :       gel(C, i) = gcoeff(B, mb + i, nb + j);
     221     3000688 :     for (; i <= m; i++)
     222           0 :       gel(C, i) = gen_0;
     223             :   }
     224      189866 :   for (; j <= ea; j++) {
     225       14650 :     gel(M, j) = C = cgetg(m + 1, t_COL);
     226      282426 :     for (i = 1; i <= da; i++)
     227      267776 :       gel(C, i) = gcoeff(A, ma + i, na + j);
     228       14650 :     for (; i <= m; i++)
     229           0 :       gel(C, i) = gen_0;
     230             :   }
     231      175216 :   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      175216 :   for (; j <= n; j++)
     239           0 :     gel(M, j) = zerocol(m);
     240      175216 :   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      153314 : 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      153314 :   long min_d = minss(da, db), min_e = minss(ea, eb), i, j;
     253      153314 :   GEN M = cgetg(n + 1, t_MAT), C;
     254             : 
     255     2943858 :   for (j = 1; j <= min_e; j++) {
     256     2790544 :     gel(M, j) = C = cgetg(m + 1, t_COL);
     257    68297775 :     for (i = 1; i <= min_d; i++)
     258   131014462 :       gel(C, i) = subii(gcoeff(A, ma + i, na + j),
     259    65507231 :                         gcoeff(B, mb + i, nb + j));
     260     2918660 :     for (; i <= da; i++)
     261      128116 :       gel(C, i) = gcoeff(A, ma + i, na + j);
     262     2981685 :     for (; i <= db; i++)
     263      191141 :       gel(C, i) = negi(gcoeff(B, mb + i, nb + j));
     264     2790544 :     for (; i <= m; i++)
     265           0 :       gel(C, i) = gen_0;
     266             :   }
     267      153314 :   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      160942 :   for (; j <= eb; j++) {
     275        7628 :     gel(M, j) = C = cgetg(m + 1, t_COL);
     276      200236 :     for (i = 1; i <= db; i++)
     277      192608 :       gel(C, i) = negi(gcoeff(B, mb + i, nb + j));
     278        7628 :     for (; i <= m; i++)
     279           0 :       gel(C, i) = gen_0;
     280             :   }
     281      160942 :   for (; j <= n; j++)
     282        7628 :     gel(M, j) = zerocol(m);
     283      153314 :   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       21902 : ZM_mul_sw(GEN A, GEN B, long m, long n, long p)
     291             : {
     292       21902 :   pari_sp av = avma;
     293       21902 :   long m1 = (m + 1)/2, m2 = m/2,
     294       21902 :     n1 = (n + 1)/2, n2 = n/2,
     295       21902 :     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       21902 :   T2 = subtract_slices(n1, p2, B, 0, n1, p1, p2, B, n1, n2, p1, p2);
     302       21902 :   S1 = subtract_slices(m2, n1, A, m1, m2, 0, n1, A, 0, m2, 0, n1);
     303       21902 :   M2 = ZM_mul_i(S1, T2, m2 + 1, n1 + 1, p2 + 1);
     304       21902 :   if (gc_needed(av, 1))
     305           0 :     gerepileall(av, 2, &T2, &M2);  /* destroy S1 */
     306       21902 :   T3 = subtract_slices(n1, p1, T2, 0, n1, 0, p2, B, 0, n1, 0, p1);
     307       21902 :   if (gc_needed(av, 1))
     308           0 :     gerepileall(av, 2, &M2, &T3);  /* destroy T2 */
     309       21902 :   S2 = add_slices(m2, n1, A, m1, m2, 0, n1, A, m1, m2, n1, n2);
     310       21902 :   T1 = subtract_slices(n1, p1, B, 0, n1, p1, p2, B, 0, n1, 0, p2);
     311       21902 :   M3 = ZM_mul_i(S2, T1, m2 + 1, n1 + 1, p2 + 1);
     312       21902 :   if (gc_needed(av, 1))
     313           0 :     gerepileall(av, 4, &M2, &T3, &S2, &M3);  /* destroy T1 */
     314       21902 :   S3 = subtract_slices(m1, n1, S2, 0, m2, 0, n1, A, 0, m1, 0, n1);
     315       21902 :   if (gc_needed(av, 1))
     316           0 :     gerepileall(av, 4, &M2, &T3, &M3, &S3);  /* destroy S2 */
     317       21902 :   A11 = matslice(A, 1, m1, 1, n1);
     318       21902 :   B11 = matslice(B, 1, n1, 1, p1);
     319       21902 :   M1 = ZM_mul_i(A11, B11, m1 + 1, n1 + 1, p1 + 1);
     320       21902 :   if (gc_needed(av, 1))
     321           0 :     gerepileall(av, 5, &M2, &T3, &M3, &S3, &M1);  /* destroy A11, B11 */
     322       21902 :   A12 = matslice(A, 1, m1, n1 + 1, n);
     323       21902 :   B21 = matslice(B, n1 + 1, n, 1, p1);
     324       21902 :   M4 = ZM_mul_i(A12, B21, m1 + 1, n2 + 1, p1 + 1);
     325       21902 :   if (gc_needed(av, 1))
     326           0 :     gerepileall(av, 6, &M2, &T3, &M3, &S3, &M1, &M4);  /* destroy A12, B21 */
     327       21902 :   C11 = add_slices(m1, p1, M1, 0, m1, 0, p1, M4, 0, m1, 0, p1);
     328       21902 :   if (gc_needed(av, 1))
     329           0 :     gerepileall(av, 6, &M2, &T3, &M3, &S3, &M1, &C11);  /* destroy M4 */
     330       21902 :   M5 = ZM_mul_i(S3, T3, m1 + 1, n1 + 1, p1 + 1);
     331       21902 :   S4 = subtract_slices(m1, n2, A, 0, m1, n1, n2, S3, 0, m1, 0, n2);
     332       21902 :   if (gc_needed(av, 1))
     333           6 :     gerepileall(av, 7, &M2, &T3, &M3, &M1, &C11, &M5, &S4);  /* destroy S3 */
     334       21902 :   T4 = add_slices(n2, p1, B, n1, n2, 0, p1, T3, 0, n2, 0, p1);
     335       21902 :   if (gc_needed(av, 1))
     336           0 :     gerepileall(av, 7, &M2, &M3, &M1, &C11, &M5, &S4, &T4);  /* destroy T3 */
     337       21902 :   V1 = subtract_slices(m1, p1, M1, 0, m1, 0, p1, M5, 0, m1, 0, p1);
     338       21902 :   if (gc_needed(av, 1))
     339           0 :     gerepileall(av, 6, &M2, &M3, &S4, &T4, &C11, &V1);  /* destroy M1, M5 */
     340       21902 :   B22 = matslice(B, n1 + 1, n, p1 + 1, p);
     341       21902 :   M6 = ZM_mul_i(S4, B22, m1 + 1, n2 + 1, p2 + 1);
     342       21902 :   if (gc_needed(av, 1))
     343           0 :     gerepileall(av, 6, &M2, &M3, &T4, &C11, &V1, &M6);  /* destroy S4, B22 */
     344       21902 :   A22 = matslice(A, m1 + 1, m, n1 + 1, n);
     345       21902 :   M7 = ZM_mul_i(A22, T4, m2 + 1, n2 + 1, p1 + 1);
     346       21902 :   if (gc_needed(av, 1))
     347           0 :     gerepileall(av, 6, &M2, &M3, &C11, &V1, &M6, &M7);  /* destroy A22, T4 */
     348       21902 :   V3 = add_slices(m1, p2, V1, 0, m1, 0, p2, M3, 0, m2, 0, p2);
     349       21902 :   C12 = add_slices(m1, p2, V3, 0, m1, 0, p2, M6, 0, m1, 0, p2);
     350       21902 :   if (gc_needed(av, 1))
     351           6 :     gerepileall(av, 6, &M2, &M3, &C11, &V1, &M7, &C12);  /* destroy V3, M6 */
     352       21902 :   V2 = add_slices(m2, p1, V1, 0, m2, 0, p1, M2, 0, m2, 0, p2);
     353       21902 :   if (gc_needed(av, 1))
     354           0 :     gerepileall(av, 5, &M3, &C11, &M7, &C12, &V2);  /* destroy V1, M2 */
     355       21902 :   C21 = add_slices(m2, p1, V2, 0, m2, 0, p1, M7, 0, m2, 0, p1);
     356       21902 :   if (gc_needed(av, 1))
     357           0 :     gerepileall(av, 5, &M3, &C11, &C12, &V2, &C21);  /* destroy M7 */
     358       21902 :   C22 = add_slices(m2, p2, V2, 0, m2, 0, p2, M3, 0, m2, 0, p2);
     359       21902 :   if (gc_needed(av, 1))
     360           0 :     gerepileall(av, 4, &C11, &C12, &C21, &C22);  /* destroy V2, M3 */
     361       21902 :   C = mkmat2(mkcol2(C11, C21), mkcol2(C12, C22));
     362       21902 :   return gerepilecopy(av, shallowmatconcat(C));
     363             : }
     364             : 
     365             : /* x[i,]*y. Assume lg(x) > 1 and 0 < i < lgcols(x) */
     366             : static GEN
     367   187774920 : ZMrow_ZC_mul_i(GEN x, GEN y, long i, long lx)
     368             : {
     369   187774920 :   pari_sp av = avma;
     370   187774920 :   GEN c = mulii(gcoeff(x,i,1), gel(y,1)), ZERO = gen_0;
     371             :   long k;
     372  2378926189 :   for (k = 2; k < lx; k++)
     373             :   {
     374  2191151556 :     GEN t = mulii(gcoeff(x,i,k), gel(y,k));
     375  2191151310 :     if (t != ZERO) c = addii(c, t);
     376             :   }
     377   187774633 :   return gerepileuptoint(av, c);
     378             : }
     379             : GEN
     380    27175191 : ZMrow_ZC_mul(GEN x, GEN y, long i)
     381    27175191 : { 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    33525965 : ZM_ZC_mul_i(GEN x, GEN y, long lx, long l)
     386             : {
     387    33525965 :   GEN z = cgetg(l,t_COL);
     388             :   long i;
     389    33525986 :   for (i=1; i<l; i++) gel(z,i) = ZMrow_ZC_mul_i(x,y,i,lx);
     390    33526002 :   return z;
     391             : }
     392             : 
     393             : static GEN
     394     8524295 : ZM_mul_classical(GEN x, GEN y, long l, long lx, long ly)
     395             : {
     396             :   long j;
     397     8524295 :   GEN z = cgetg(ly, t_MAT);
     398    33523406 :   for (j = 1; j < ly; j++)
     399    24999111 :     gel(z, j) = ZM_ZC_mul_i(x, gel(y, j), lx, l);
     400     8524295 :   return z;
     401             : }
     402             : 
     403             : /* Strassen-Winograd used for dim >= ZM_sw_bound */
     404             : static GEN
     405     8542291 : ZM_mul_i(GEN x, GEN y, long l, long lx, long ly)
     406             : {
     407     8542291 :   long s = maxss(ZM_max_lg_i(x,lx,l), ZM_max_lg_i(y,ly,lx));
     408     8542291 :   long ZM_sw_bound = s > 60 ? 2: s > 25 ? 4: s>15 ? 8 : s > 8 ? 16 : 32;
     409     8542291 :   if (l <= ZM_sw_bound || lx <= ZM_sw_bound || ly <= ZM_sw_bound)
     410     8520398 :     return ZM_mul_classical(x, y, l, lx, ly);
     411             :   else
     412       21893 :     return ZM_mul_sw(x, y, l - 1, lx - 1, ly - 1);
     413             : }
     414             : 
     415             : GEN
     416     8453265 : ZM_mul(GEN x, GEN y)
     417             : {
     418     8453265 :   long lx=lg(x), ly=lg(y);
     419     8453265 :   if (ly==1) return cgetg(1,t_MAT);
     420     8389880 :   if (lx==1) return zeromat(0, ly-1);
     421     8388977 :   return ZM_mul_i(x, y, lgcols(x), lx, ly);
     422             : }
     423             : 
     424             : GEN
     425      248075 : QM_mul(GEN x, GEN y)
     426             : {
     427      248075 :   GEN dx, nx = Q_primitive_part(x, &dx);
     428      248075 :   GEN dy, ny = Q_primitive_part(y, &dy);
     429      248075 :   GEN z = ZM_mul(nx, ny);
     430      248075 :   if (dx || dy)
     431             :   {
     432      191203 :     GEN d = dx ? dy ? gmul(dx, dy): dx : dy;
     433      191203 :     if (!gequal1(d)) z = ZM_Q_mul(z, d);
     434             :   }
     435      248075 :   return z;
     436             : }
     437             : 
     438             : GEN
     439          91 : QM_QC_mul(GEN x, GEN y)
     440             : {
     441          91 :   GEN dx, nx = Q_primitive_part(x, &dx);
     442          91 :   GEN dy, ny = Q_primitive_part(y, &dy);
     443          91 :   GEN z = ZM_ZC_mul(nx, ny);
     444          91 :   if (dx || dy)
     445             :   {
     446          91 :     GEN d = dx ? dy ? gmul(dx, dy): dx : dy;
     447          91 :     if (!gequal1(d)) z = ZC_Q_mul(z, d);
     448             :   }
     449          91 :   return z;
     450             : }
     451             : 
     452             : /* assume result is symmetric */
     453             : GEN
     454           0 : ZM_multosym(GEN x, GEN y)
     455             : {
     456           0 :   long j, lx, ly = lg(y);
     457             :   GEN M;
     458           0 :   if (ly == 1) return cgetg(1,t_MAT);
     459           0 :   lx = lg(x); /* = lgcols(y) */
     460           0 :   if (lx == 1) return cgetg(1,t_MAT);
     461             :   /* ly = lgcols(x) */
     462           0 :   M = cgetg(ly, t_MAT);
     463           0 :   for (j=1; j<ly; j++)
     464             :   {
     465           0 :     GEN z = cgetg(ly,t_COL), yj = gel(y,j);
     466             :     long i;
     467           0 :     for (i=1; i<j; i++) gel(z,i) = gcoeff(M,j,i);
     468           0 :     for (i=j; i<ly; i++)gel(z,i) = ZMrow_ZC_mul_i(x,yj,i,lx);
     469           0 :     gel(M,j) = z;
     470             :   }
     471           0 :   return M;
     472             : }
     473             : 
     474             : /* compute m*diagonal(d), assume lg(d) = lg(m). Shallow */
     475             : GEN
     476          14 : ZM_mul_diag(GEN m, GEN d)
     477             : {
     478             :   long j, l;
     479          14 :   GEN y = cgetg_copy(m, &l);
     480          56 :   for (j=1; j<l; j++)
     481             :   {
     482          42 :     GEN c = gel(d,j);
     483          42 :     gel(y,j) = equali1(c)? gel(m,j): ZC_Z_mul(gel(m,j), c);
     484             :   }
     485          14 :   return y;
     486             : }
     487             : /* compute diagonal(d)*m, assume lg(d) = lg(m~). Shallow */
     488             : GEN
     489      132649 : ZM_diag_mul(GEN d, GEN m)
     490             : {
     491      132649 :   long i, j, l = lg(d), lm = lg(m);
     492      132649 :   GEN y = cgetg(lm, t_MAT);
     493      132649 :   for (j=1; j<lm; j++) gel(y,j) = cgetg(l, t_COL);
     494      432186 :   for (i=1; i<l; i++)
     495             :   {
     496      299537 :     GEN c = gel(d,i);
     497      299537 :     if (equali1(c))
     498       13944 :       for (j=1; j<lm; j++) gcoeff(y,i,j) = gcoeff(m,i,j);
     499             :     else
     500      285593 :       for (j=1; j<lm; j++) gcoeff(y,i,j) = mulii(gcoeff(m,i,j), c);
     501             :   }
     502      132649 :   return y;
     503             : }
     504             : 
     505             : /* assume lx > 1 is lg(x) = lg(y) */
     506             : static GEN
     507    13231216 : ZV_dotproduct_i(GEN x, GEN y, long lx)
     508             : {
     509    13231216 :   pari_sp av = avma;
     510    13231216 :   GEN c = mulii(gel(x,1), gel(y,1));
     511             :   long i;
     512   101866905 :   for (i = 2; i < lx; i++)
     513             :   {
     514    88635689 :     GEN t = mulii(gel(x,i), gel(y,i));
     515    88635689 :     if (t != gen_0) c = addii(c, t);
     516             :   }
     517    13231216 :   return gerepileuptoint(av, c);
     518             : }
     519             : 
     520             : /* x~ * y, assuming result is symmetric */
     521             : GEN
     522      161595 : ZM_transmultosym(GEN x, GEN y)
     523             : {
     524      161595 :   long i, j, l, ly = lg(y);
     525             :   GEN M;
     526      161595 :   if (ly == 1) return cgetg(1,t_MAT);
     527             :   /* lg(x) = ly */
     528      161595 :   l = lgcols(y); /* = lgcols(x) */
     529      161595 :   M = cgetg(ly, t_MAT);
     530     1265890 :   for (i=1; i<ly; i++)
     531             :   {
     532     1104295 :     GEN xi = gel(x,i), c = cgetg(ly,t_COL);
     533     1104295 :     gel(M,i) = c;
     534     5065656 :     for (j=1; j<i; j++)
     535     3961361 :       gcoeff(M,i,j) = gel(c,j) = ZV_dotproduct_i(xi,gel(y,j),l);
     536     1104295 :     gel(c,i) = ZV_dotproduct_i(xi,gel(y,i),l);
     537             :   }
     538      161595 :   return M;
     539             : }
     540             : /* x~ * y */
     541             : GEN
     542         497 : ZM_transmul(GEN x, GEN y)
     543             : {
     544         497 :   long i, j, l, lx, ly = lg(y);
     545             :   GEN M;
     546         497 :   if (ly == 1) return cgetg(1,t_MAT);
     547         497 :   lx = lg(x);
     548         497 :   l = lgcols(y);
     549         497 :   if (lgcols(x) != l) pari_err_OP("operation 'ZM_transmul'", x,y);
     550         497 :   M = cgetg(ly, t_MAT);
     551        1617 :   for (i=1; i<ly; i++)
     552             :   {
     553        1120 :     GEN yi = gel(y,i), c = cgetg(lx,t_COL);
     554        1120 :     gel(M,i) = c;
     555        1120 :     for (j=1; j<lx; j++) gel(c,j) = ZV_dotproduct_i(yi,gel(x,j),l);
     556             :   }
     557         497 :   return M;
     558             : }
     559             : 
     560             : static GEN
     561        3906 : ZM_sqr_i(GEN x, long l)
     562             : {
     563        3906 :   long s = ZM_max_lg_i(x,l,l);
     564        3906 :   long ZM_sw_bound = s > 60 ? 2: s > 25 ? 4: s>15 ? 8 : s > 8 ? 16 : 32;
     565        3906 :   if (l <= ZM_sw_bound)
     566        3897 :     return ZM_mul_classical(x, x, l, l, l);
     567             :   else
     568           9 :     return ZM_mul_sw(x, x, l - 1, l - 1, l - 1);
     569             : }
     570             : 
     571             : GEN
     572        3906 : ZM_sqr(GEN x)
     573             : {
     574        3906 :   long lx=lg(x);
     575        3906 :   if (lx==1) return cgetg(1,t_MAT);
     576        3906 :   return ZM_sqr_i(x, lx);
     577             : }
     578             : GEN
     579     8620232 : ZM_ZC_mul(GEN x, GEN y)
     580             : {
     581     8620232 :   long lx = lg(x);
     582     8620232 :   return lx==1? cgetg(1,t_COL): ZM_ZC_mul_i(x, y, lx, lgcols(x));
     583             : }
     584             : 
     585             : GEN
     586     1139753 : ZC_Z_div(GEN x, GEN c)
     587     1139753 : { pari_APPLY_type(t_COL, Qdivii(gel(x,i), c)) }
     588             : 
     589             : GEN
     590       14707 : ZM_Z_div(GEN x, GEN c)
     591       14707 : { pari_APPLY_same(ZC_Z_div(gel(x, i), c)) }
     592             : 
     593             : GEN
     594     1230630 : ZC_Q_mul(GEN A, GEN z)
     595             : {
     596     1230630 :   pari_sp av = avma;
     597     1230630 :   long i, l = lg(A);
     598             :   GEN d, n, Ad, B, u;
     599     1230630 :   if (typ(z)==t_INT) return ZC_Z_mul(A,z);
     600     1230252 :   n = gel(z, 1); d = gel(z, 2);
     601     1230252 :   Ad = FpC_red(A, d);
     602     1230252 :   u = gcdii(d, FpV_factorback(Ad, NULL, d));
     603     1230252 :   B = cgetg(l, t_COL);
     604     1230252 :   if (equali1(u))
     605             :   {
     606      191008 :     for(i=1; i<l; i++)
     607      152323 :       gel(B, i) = mkfrac(mulii(n, gel(A,i)), d);
     608             :   } else
     609             :   {
     610    11985289 :     for(i=1; i<l; i++)
     611             :     {
     612    10793722 :       GEN di = gcdii(gel(Ad, i), u), ni = mulii(n, diviiexact(gel(A,i), di));
     613    10793722 :       if (equalii(d, di))
     614     8052620 :         gel(B, i) = ni;
     615             :       else
     616     2741102 :         gel(B, i) = mkfrac(ni, diviiexact(d, di));
     617             :     }
     618             :   }
     619     1230252 :   return gerepilecopy(av, B);
     620             : }
     621             : 
     622             : GEN
     623      334154 : ZM_Q_mul(GEN x, GEN z)
     624             : {
     625      334154 :   if (typ(z)==t_INT) return ZM_Z_mul(x,z);
     626      212913 :   pari_APPLY_same(ZC_Q_mul(gel(x, i), z));
     627             : }
     628             : 
     629             : long
     630   182225757 : zv_dotproduct(GEN x, GEN y)
     631             : {
     632   182225757 :   long i, lx = lg(x);
     633             :   ulong c;
     634   182225757 :   if (lx == 1) return 0;
     635   182225757 :   c = uel(x,1)*uel(y,1);
     636  2825754764 :   for (i = 2; i < lx; i++)
     637  2643529007 :     c += uel(x,i)*uel(y,i);
     638   182225757 :   return c;
     639             : }
     640             : 
     641             : GEN
     642      163107 : ZV_ZM_mul(GEN x, GEN y)
     643             : {
     644      163107 :   long i, lx = lg(x), ly = lg(y);
     645             :   GEN z;
     646      163107 :   if (lx == 1) return zerovec(ly-1);
     647      163093 :   z = cgetg(ly, t_VEC);
     648      163093 :   for (i = 1; i < ly; i++) gel(z,i) = ZV_dotproduct_i(x, gel(y,i), lx);
     649      163093 :   return z;
     650             : }
     651             : 
     652             : GEN
     653           0 : ZC_ZV_mul(GEN x, GEN y)
     654             : {
     655           0 :   long i,j, lx=lg(x), ly=lg(y);
     656             :   GEN z;
     657           0 :   if (ly==1) return cgetg(1,t_MAT);
     658           0 :   z = cgetg(ly,t_MAT);
     659           0 :   for (j=1; j < ly; j++)
     660             :   {
     661           0 :     gel(z,j) = cgetg(lx,t_COL);
     662           0 :     for (i=1; i<lx; i++) gcoeff(z,i,j) = mulii(gel(x,i),gel(y,j));
     663             :   }
     664           0 :   return z;
     665             : }
     666             : 
     667             : GEN
     668     6115247 : ZV_dotsquare(GEN x)
     669             : {
     670             :   long i, lx;
     671             :   pari_sp av;
     672             :   GEN z;
     673     6115247 :   lx = lg(x);
     674     6115247 :   if (lx == 1) return gen_0;
     675     6115247 :   av = avma; z = sqri(gel(x,1));
     676     6115247 :   for (i=2; i<lx; i++) z = addii(z, sqri(gel(x,i)));
     677     6115247 :   return gerepileuptoint(av,z);
     678             : }
     679             : 
     680             : GEN
     681    11433219 : ZV_dotproduct(GEN x,GEN y)
     682             : {
     683             :   long lx;
     684    11433219 :   if (x == y) return ZV_dotsquare(x);
     685     7705555 :   lx = lg(x);
     686     7705555 :   if (lx == 1) return gen_0;
     687     7705555 :   return ZV_dotproduct_i(x, y, lx);
     688             : }
     689             : 
     690             : static GEN
     691         273 : _ZM_mul(void *data /*ignored*/, GEN x, GEN y)
     692         273 : { (void)data; return ZM_mul(x,y); }
     693             : static GEN
     694        2891 : _ZM_sqr(void *data /*ignored*/, GEN x)
     695        2891 : { (void)data; return ZM_sqr(x); }
     696             : GEN
     697           0 : ZM_pow(GEN x, GEN n)
     698             : {
     699           0 :   pari_sp av = avma;
     700           0 :   if (!signe(n)) return matid(lg(x)-1);
     701           0 :   return gerepileupto(av, gen_pow(x, n, NULL, &_ZM_sqr, &_ZM_mul));
     702             : }
     703             : GEN
     704        2450 : ZM_powu(GEN x, ulong n)
     705             : {
     706        2450 :   pari_sp av = avma;
     707        2450 :   if (!n) return matid(lg(x)-1);
     708        2450 :   return gerepileupto(av, gen_powu(x, n, NULL, &_ZM_sqr, &_ZM_mul));
     709             : }
     710             : /********************************************************************/
     711             : /**                                                                **/
     712             : /**                           ADD, SUB                             **/
     713             : /**                                                                **/
     714             : /********************************************************************/
     715             : static GEN
     716    16364831 : ZC_add_i(GEN x, GEN y, long lx)
     717             : {
     718    16364831 :   GEN A = cgetg(lx, t_COL);
     719             :   long i;
     720    16365194 :   for (i=1; i<lx; i++) gel(A,i) = addii(gel(x,i), gel(y,i));
     721    16364694 :   return A;
     722             : }
     723             : GEN
     724     7469129 : ZC_add(GEN x, GEN y) { return ZC_add_i(x, y, lg(x)); }
     725             : GEN
     726      129778 : ZC_Z_add(GEN x, GEN y)
     727             : {
     728      129778 :   long k, lx = lg(x);
     729      129778 :   GEN z = cgetg(lx, t_COL);
     730      129778 :   if (lx == 1) pari_err_TYPE2("+",x,y);
     731      129778 :   gel(z,1) = addii(y,gel(x,1));
     732      129778 :   for (k = 2; k < lx; k++) gel(z,k) = icopy(gel(x,k));
     733      129778 :   return z;
     734             : }
     735             : 
     736             : static GEN
     737     2094015 : ZC_sub_i(GEN x, GEN y, long lx)
     738             : {
     739             :   long i;
     740     2094015 :   GEN A = cgetg(lx, t_COL);
     741     2094015 :   for (i=1; i<lx; i++) gel(A,i) = subii(gel(x,i), gel(y,i));
     742     2094015 :   return A;
     743             : }
     744             : GEN
     745     1329317 : ZC_sub(GEN x, GEN y) { return ZC_sub_i(x, y, lg(x)); }
     746             : GEN
     747           0 : ZC_Z_sub(GEN x, GEN y)
     748             : {
     749           0 :   long k, lx = lg(x);
     750           0 :   GEN z = cgetg(lx, t_COL);
     751           0 :   if (lx == 1) pari_err_TYPE2("+",x,y);
     752           0 :   gel(z,1) = subii(gel(x,1), y);
     753           0 :   for (k = 2; k < lx; k++) gel(z,k) = icopy(gel(x,k));
     754           0 :   return z;
     755             : }
     756             : GEN
     757      107806 : Z_ZC_sub(GEN a, GEN x)
     758             : {
     759      107806 :   long k, lx = lg(x);
     760      107806 :   GEN z = cgetg(lx, t_COL);
     761      107821 :   if (lx == 1) pari_err_TYPE2("-",a,x);
     762      107821 :   gel(z,1) = subii(a, gel(x,1));
     763      107771 :   for (k = 2; k < lx; k++) gel(z,k) = negi(gel(x,k));
     764      107820 :   return z;
     765             : }
     766             : 
     767             : GEN
     768     1277172 : ZM_add(GEN x, GEN y)
     769             : {
     770     1277172 :   long lx = lg(x), l, j;
     771             :   GEN z;
     772     1277172 :   if (lx == 1) return cgetg(1, t_MAT);
     773     1269402 :   z = cgetg(lx, t_MAT); l = lgcols(x);
     774     1269402 :   for (j = 1; j < lx; j++) gel(z,j) = ZC_add_i(gel(x,j), gel(y,j), l);
     775     1269402 :   return z;
     776             : }
     777             : GEN
     778      692439 : ZM_sub(GEN x, GEN y)
     779             : {
     780      692439 :   long lx = lg(x), l, j;
     781             :   GEN z;
     782      692439 :   if (lx == 1) return cgetg(1, t_MAT);
     783      692439 :   z = cgetg(lx, t_MAT); l = lgcols(x);
     784      692439 :   for (j = 1; j < lx; j++) gel(z,j) = ZC_sub_i(gel(x,j), gel(y,j), l);
     785      692439 :   return z;
     786             : }
     787             : /********************************************************************/
     788             : /**                                                                **/
     789             : /**                         LINEAR COMBINATION                     **/
     790             : /**                                                                **/
     791             : /********************************************************************/
     792             : /* return X/c assuming division is exact */
     793             : GEN
     794     3169232 : ZC_Z_divexact(GEN x, GEN c)
     795     3169232 : { pari_APPLY_type(t_COL, diviiexact(gel(x,i), c)) }
     796             : 
     797             : GEN
     798      868656 : ZM_Z_divexact(GEN x, GEN c)
     799      868656 : { pari_APPLY_same(ZC_Z_divexact(gel(x,i), c)) }
     800             : 
     801             : GEN
     802    14127388 : ZC_Z_mul(GEN x, GEN c)
     803             : {
     804    14127388 :   if (!signe(c)) return zerocol(lg(x)-1);
     805    13706435 :   if (is_pm1(c)) return (signe(c) > 0)? ZC_copy(x): ZC_neg(x);
     806    11322402 :   pari_APPLY_type(t_COL, mulii(gel(x,i), c))
     807             : }
     808             : 
     809             : GEN
     810       15658 : ZC_z_mul(GEN x, long c)
     811             : {
     812       15658 :   if (!c) return zerocol(lg(x)-1);
     813        8511 :   if (c == 1) return ZC_copy(x);
     814        4031 :   if (c ==-1) return ZC_neg(x);
     815        2344 :   pari_APPLY_type(t_COL, mulsi(c, gel(x,i)))
     816             : }
     817             : 
     818             : GEN
     819       16603 : zv_z_mul(GEN x, long n)
     820       16603 : { pari_APPLY_long(x[i]*n) }
     821             : 
     822             : /* return a ZM */
     823             : GEN
     824      671839 : nm_Z_mul(GEN X, GEN c)
     825             : {
     826      671839 :   long i, j, h, l = lg(X), s = signe(c);
     827             :   GEN A;
     828      671839 :   if (l == 1) return cgetg(1, t_MAT);
     829      671839 :   h = lgcols(X);
     830      671839 :   if (!s) return zeromat(h-1, l-1);
     831      671839 :   if (is_pm1(c)) {
     832           0 :     if (s > 0) return Flm_to_ZM(X);
     833           0 :     X = Flm_to_ZM(X); ZM_togglesign(X); return X;
     834             :   }
     835      671839 :   A = cgetg(l, t_MAT);
     836     1360750 :   for (j = 1; j < l; j++)
     837             :   {
     838      688911 :     GEN a = cgetg(h, t_COL), x = gel(X, j);
     839      688911 :     for (i = 1; i < h; i++) gel(a,i) = muliu(c, x[i]);
     840      688911 :     gel(A,j) = a;
     841             :   }
     842      671839 :   return A;
     843             : }
     844             : GEN
     845     1054608 : ZM_Z_mul(GEN X, GEN c)
     846             : {
     847     1054608 :   long i, j, h, l = lg(X);
     848             :   GEN A;
     849     1054608 :   if (l == 1) return cgetg(1, t_MAT);
     850     1039600 :   h = lgcols(X);
     851     1039600 :   if (!signe(c)) return zeromat(h-1, l-1);
     852     1038130 :   if (is_pm1(c)) return (signe(c) > 0)? ZM_copy(X): ZM_neg(X);
     853      801582 :   A = cgetg(l, t_MAT);
     854     9403256 :   for (j = 1; j < l; j++)
     855             :   {
     856     8601674 :     GEN a = cgetg(h, t_COL), x = gel(X, j);
     857     8601674 :     for (i = 1; i < h; i++) gel(a,i) = mulii(c, gel(x,i));
     858     8601674 :     gel(A,j) = a;
     859             :   }
     860      801582 :   return A;
     861             : }
     862             : void
     863    47273197 : ZC_lincomb1_inplace_i(GEN X, GEN Y, GEN v, long n)
     864             : {
     865             :   long i;
     866    47273197 :   for (i = n; i; i--) gel(X,i) = addmulii_inplace(gel(X,i), gel(Y,i), v);
     867    47266482 : }
     868             : /* X <- X + v Y (elementary col operation) */
     869             : void
     870    41640512 : ZC_lincomb1_inplace(GEN X, GEN Y, GEN v)
     871             : {
     872    41640512 :   if (lgefint(v) != 2) return ZC_lincomb1_inplace_i(X, Y, v, lg(X)-1);
     873             : }
     874             : void
     875     8903725 : Flc_lincomb1_inplace(GEN X, GEN Y, ulong v, ulong q)
     876             : {
     877             :   long i;
     878     8903725 :   if (!v) return; /* v = 0 */
     879     8903725 :   for (i = lg(X)-1; i; i--) X[i] = Fl_add(X[i], Fl_mul(Y[i], v, q), q);
     880             : }
     881             : 
     882             : /* X + v Y, wasteful if (v = 0) */
     883             : static GEN
     884     1934810 : ZC_lincomb1(GEN v, GEN x, GEN y)
     885     1934810 : { pari_APPLY_type(t_COL, addmulii(gel(x,i), gel(y,i), v)) }
     886             : 
     887             : /* -X + vY */
     888             : static GEN
     889      250275 : ZC_lincomb_1(GEN v, GEN x, GEN y)
     890      250275 : { pari_APPLY_type(t_COL, mulsubii(gel(y,i), v, gel(x,i))) }
     891             : 
     892             : /* X,Y compatible ZV; u,v in Z. Returns A = u*X + v*Y */
     893             : GEN
     894     9260809 : ZC_lincomb(GEN u, GEN v, GEN X, GEN Y)
     895             : {
     896             :   long su, sv;
     897             :   GEN A;
     898             : 
     899     9260809 :   su = signe(u); if (!su) return ZC_Z_mul(Y, v);
     900     9260515 :   sv = signe(v); if (!sv) return ZC_Z_mul(X, u);
     901     9260193 :   if (is_pm1(v))
     902             :   {
     903     1598707 :     if (is_pm1(u))
     904             :     {
     905     1108061 :       if (su != sv) A = ZC_sub(X, Y);
     906      341200 :       else          A = ZC_add(X, Y);
     907     1108061 :       if (su < 0) ZV_togglesign(A); /* in place but was created above */
     908             :     }
     909             :     else
     910             :     {
     911      490646 :       if (sv > 0) A = ZC_lincomb1 (u, Y, X);
     912      202447 :       else        A = ZC_lincomb_1(u, Y, X);
     913             :     }
     914             :   }
     915     7661533 :   else if (is_pm1(u))
     916             :   {
     917     1694439 :     if (su > 0) A = ZC_lincomb1 (v, X, Y);
     918       47828 :     else        A = ZC_lincomb_1(v, X, Y);
     919             :   }
     920             :   else
     921             :   { /* not cgetg_copy: x may be a t_VEC */
     922     5967153 :     long i, lx = lg(X);
     923     5967153 :     A = cgetg(lx,t_COL);
     924     5967233 :     for (i=1; i<lx; i++) gel(A,i) = lincombii(u,v,gel(X,i),gel(Y,i));
     925             :   }
     926     9260058 :   return A;
     927             : }
     928             : 
     929             : /********************************************************************/
     930             : /**                                                                **/
     931             : /**                           CONVERSIONS                          **/
     932             : /**                                                                **/
     933             : /********************************************************************/
     934             : GEN
     935       75803 : ZV_to_nv(GEN x)
     936       75803 : { pari_APPLY_ulong(itou(gel(x,i))) }
     937             : 
     938             : GEN
     939       29854 : zm_to_ZM(GEN x)
     940       29854 : { pari_APPLY_type(t_MAT, zc_to_ZC(gel(x,i))) }
     941             : 
     942             : GEN
     943         112 : zmV_to_ZMV(GEN x)
     944         112 : { pari_APPLY_type(t_VEC, zm_to_ZM(gel(x,i))) }
     945             : 
     946             : /* same as Flm_to_ZM but do not assume positivity */
     947             : GEN
     948        1120 : ZM_to_zm(GEN x)
     949        1120 : { pari_APPLY_same(ZV_to_zv(gel(x,i))) }
     950             : 
     951             : GEN
     952      284480 : zv_to_Flv(GEN x, ulong p)
     953      284480 : { pari_APPLY_ulong(umodsu(x[i], p)) }
     954             : 
     955             : GEN
     956       19894 : zm_to_Flm(GEN x, ulong p)
     957       19894 : { pari_APPLY_same(zv_to_Flv(gel(x,i),p)) }
     958             : 
     959             : GEN
     960          35 : ZMV_to_zmV(GEN x)
     961          35 : { pari_APPLY_type(t_VEC, ZM_to_zm(gel(x,i))) }
     962             : 
     963             : /********************************************************************/
     964             : /**                                                                **/
     965             : /**                         COPY, NEGATION                         **/
     966             : /**                                                                **/
     967             : /********************************************************************/
     968             : GEN
     969     4809568 : ZC_copy(GEN x)
     970             : {
     971     4809568 :   long i, lx = lg(x);
     972     4809568 :   GEN y = cgetg(lx, t_COL);
     973    43917201 :   for (i=1; i<lx; i++)
     974             :   {
     975    39107527 :     GEN c = gel(x,i);
     976    39107527 :     gel(y,i) = lgefint(c) == 2? gen_0: icopy(c);
     977             :   }
     978     4809674 :   return y;
     979             : }
     980             : 
     981             : GEN
     982      281104 : ZM_copy(GEN x)
     983      281104 : { pari_APPLY_same(ZC_copy(gel(x,i))) }
     984             : 
     985             : void
     986       86547 : ZV_neg_inplace(GEN M)
     987             : {
     988       86547 :   long l = lg(M);
     989       86547 :   while (--l > 0) gel(M,l) = negi(gel(M,l));
     990       86547 : }
     991             : GEN
     992     1837671 : ZC_neg(GEN x)
     993     1837671 : { pari_APPLY_type(t_COL, negi(gel(x,i))) }
     994             : GEN
     995       31561 : zv_neg(GEN x)
     996       31561 : { pari_APPLY_long(-x[i]) }
     997             : GEN
     998         227 : zv_neg_inplace(GEN M)
     999             : {
    1000         227 :   long l = lg(M);
    1001         227 :   while (--l > 0) M[l] = -M[l];
    1002         227 :   return M;
    1003             : }
    1004             : GEN
    1005      278075 : ZM_neg(GEN x)
    1006      278075 : { pari_APPLY_same(ZC_neg(gel(x,i))) }
    1007             : 
    1008             : void
    1009     1424670 : ZV_togglesign(GEN M)
    1010             : {
    1011     1424670 :   long l = lg(M);
    1012     1424670 :   while (--l > 0) togglesign_safe(&gel(M,l));
    1013     1424670 : }
    1014             : void
    1015           0 : ZM_togglesign(GEN M)
    1016             : {
    1017           0 :   long l = lg(M);
    1018           0 :   while (--l > 0) ZV_togglesign(gel(M,l));
    1019           0 : }
    1020             : 
    1021             : /********************************************************************/
    1022             : /**                                                                **/
    1023             : /**                        "DIVISION" mod HNF                      **/
    1024             : /**                                                                **/
    1025             : /********************************************************************/
    1026             : /* Reduce ZC x modulo ZM y in HNF, may return x itself (not a copy) */
    1027             : GEN
    1028      920571 : ZC_hnfremdiv(GEN x, GEN y, GEN *Q)
    1029             : {
    1030      920571 :   long i, l = lg(x);
    1031             :   GEN q;
    1032             : 
    1033      920571 :   if (Q) *Q = cgetg(l,t_COL);
    1034      920571 :   if (l == 1) return cgetg(1,t_COL);
    1035     5299977 :   for (i = l-1; i>0; i--)
    1036             :   {
    1037     4379406 :     q = diviiround(gel(x,i), gcoeff(y,i,i));
    1038     4379406 :     if (signe(q)) {
    1039     2359855 :       togglesign(q);
    1040     2359855 :       x = ZC_lincomb(gen_1, q, x, gel(y,i));
    1041             :     }
    1042     4379406 :     if (Q) gel(*Q, i) = q;
    1043             :   }
    1044      920571 :   return x;
    1045             : }
    1046             : 
    1047             : /* x = y Q + R, may return some columns of x (not copies) */
    1048             : GEN
    1049       71050 : ZM_hnfdivrem(GEN x, GEN y, GEN *Q)
    1050             : {
    1051       71050 :   long lx = lg(x), i;
    1052       71050 :   GEN R = cgetg(lx, t_MAT);
    1053       71050 :   if (Q)
    1054             :   {
    1055       22184 :     GEN q = cgetg(lx, t_MAT); *Q = q;
    1056       22184 :     for (i=1; i<lx; i++) gel(R,i) = ZC_hnfremdiv(gel(x,i),y,(GEN*)(q+i));
    1057             :   }
    1058             :   else
    1059      150665 :     for (i=1; i<lx; i++)
    1060             :     {
    1061      101799 :       pari_sp av = avma;
    1062      101799 :       GEN z = ZC_hnfrem(gel(x,i),y);
    1063      101799 :       gel(R,i) = (avma == av)? ZC_copy(z): gerepileupto(av, z);
    1064             :     }
    1065       71050 :   return R;
    1066             : }
    1067             : 
    1068             : 
    1069             : /********************************************************************/
    1070             : /**                                                                **/
    1071             : /**                               TESTS                            **/
    1072             : /**                                                                **/
    1073             : /********************************************************************/
    1074             : int
    1075    21181440 : zv_equal0(GEN V)
    1076             : {
    1077    21181440 :   long l = lg(V);
    1078    55541739 :   while (--l > 0)
    1079    28035402 :     if (V[l]) return 0;
    1080     6324897 :   return 1;
    1081             : }
    1082             : 
    1083             : int
    1084      956044 : ZV_equal0(GEN V)
    1085             : {
    1086      956044 :   long l = lg(V);
    1087     3819166 :   while (--l > 0)
    1088     2630890 :     if (signe(gel(V,l))) return 0;
    1089      232232 :   return 1;
    1090             : }
    1091             : 
    1092             : static int
    1093    15956004 : ZV_equal_lg(GEN V, GEN W, long l)
    1094             : {
    1095    40688350 :   while (--l > 0)
    1096    22663292 :     if (!equalii(gel(V,l), gel(W,l))) return 0;
    1097     2069054 :   return 1;
    1098             : }
    1099             : int
    1100    14612760 : ZV_equal(GEN V, GEN W)
    1101             : {
    1102    14612760 :   long l = lg(V);
    1103    14612760 :   if (lg(W) != l) return 0;
    1104    14612760 :   return ZV_equal_lg(V, W, l);
    1105             : }
    1106             : int
    1107      497530 : ZM_equal(GEN A, GEN B)
    1108             : {
    1109      497530 :   long i, m, l = lg(A);
    1110      497530 :   if (lg(B) != l) return 0;
    1111      496648 :   if (l == 1) return 1;
    1112      496648 :   m = lgcols(A);
    1113      496648 :   if (lgcols(B) != m) return 0;
    1114     1732147 :   for (i = 1; i < l; i++)
    1115     1343222 :     if (!ZV_equal_lg(gel(A,i), gel(B,i), m)) return 0;
    1116      388925 :   return 1;
    1117             : }
    1118             : int
    1119          77 : ZM_equal0(GEN A)
    1120             : {
    1121          77 :   long i, j, m, l = lg(A);
    1122          77 :   if (l == 1) return 1;
    1123          77 :   m = lgcols(A);
    1124          84 :   for (j = 1; j < l; j++)
    1125         308 :     for (i = 1; i < m; i++)
    1126         301 :       if (signe(gcoeff(A,i,j))) return 0;
    1127           0 :   return 1;
    1128             : }
    1129             : int
    1130    21644542 : zv_equal(GEN V, GEN W)
    1131             : {
    1132    21644542 :   long l = lg(V);
    1133    21644542 :   if (lg(W) != l) return 0;
    1134   203149942 :   while (--l > 0)
    1135   160634517 :     if (V[l] != W[l]) return 0;
    1136    21235528 :   return 1;
    1137             : }
    1138             : 
    1139             : int
    1140        1337 : zvV_equal(GEN V, GEN W)
    1141             : {
    1142        1337 :   long l = lg(V);
    1143        1337 :   if (lg(W) != l) return 0;
    1144       79037 :   while (--l > 0)
    1145       77301 :     if (!zv_equal(gel(V,l),gel(W,l))) return 0;
    1146         406 :   return 1;
    1147             : }
    1148             : 
    1149             : int
    1150      161359 : ZM_ishnf(GEN x)
    1151             : {
    1152      161359 :   long i,j, lx = lg(x);
    1153      482755 :   for (i=1; i<lx; i++)
    1154             :   {
    1155      367365 :     GEN xii = gcoeff(x,i,i);
    1156      367365 :     if (signe(xii) <= 0) return 0;
    1157      800451 :     for (j=1; j<i; j++)
    1158      459483 :       if (signe(gcoeff(x,i,j))) return 0;
    1159      855709 :     for (j=i+1; j<lx; j++)
    1160             :     {
    1161      534313 :       GEN xij = gcoeff(x,i,j);
    1162      534313 :       if (signe(xij)<0 || cmpii(xij,xii)>=0) return 0;
    1163             :     }
    1164             :   }
    1165      115390 :   return 1;
    1166             : }
    1167             : int
    1168      121254 : ZM_isidentity(GEN x)
    1169             : {
    1170      121254 :   long i,j, lx = lg(x);
    1171             : 
    1172      121254 :   if (lx == 1) return 1;
    1173      121254 :   if (lx != lgcols(x)) return 0;
    1174      874148 :   for (j=1; j<lx; j++)
    1175             :   {
    1176      752915 :     GEN c = gel(x,j);
    1177     4234431 :     for (i=1; i<j; )
    1178     2728601 :       if (signe(gel(c,i++))) return 0;
    1179             :     /* i = j */
    1180      752915 :       if (!equali1(gel(c,i++))) return 0;
    1181     4234403 :     for (   ; i<lx; )
    1182     2728608 :       if (signe(gel(c,i++))) return 0;
    1183             :   }
    1184      121233 :   return 1;
    1185             : }
    1186             : int
    1187       81765 : ZM_isdiagonal(GEN x)
    1188             : {
    1189       81765 :   long i,j, lx = lg(x);
    1190       81765 :   if (lx == 1) return 1;
    1191       81765 :   if (lx != lgcols(x)) return 0;
    1192             : 
    1193      239957 :   for (j=1; j<lx; j++)
    1194             :   {
    1195      207058 :     GEN c = gel(x,j);
    1196      414853 :     for (i=1; i<j; i++)
    1197      256654 :       if (signe(gel(c,i))) return 0;
    1198      586258 :     for (i++; i<lx; i++)
    1199      428066 :       if (signe(gel(c,i))) return 0;
    1200             :   }
    1201       32899 :   return 1;
    1202             : }
    1203             : int
    1204      105623 : ZM_isscalar(GEN x, GEN s)
    1205             : {
    1206      105623 :   long i, j, lx = lg(x);
    1207             : 
    1208      105623 :   if (lx == 1) return 1;
    1209      105623 :   if (!s) s = gcoeff(x,1,1);
    1210      105623 :   if (equali1(s)) return ZM_isidentity(x);
    1211       21070 :   if (lx != lgcols(x)) return 0;
    1212      195166 :   for (j=1; j<lx; j++)
    1213             :   {
    1214      175824 :     GEN c = gel(x,j);
    1215     1765984 :     for (i=1; i<j; )
    1216     1415244 :       if (signe(gel(c,i++))) return 0;
    1217             :     /* i = j */
    1218      174916 :       if (!equalii(gel(c,i++), s)) return 0;
    1219     1772368 :     for (   ; i<lx; )
    1220     1424121 :       if (signe(gel(c,i++))) return 0;
    1221             :   }
    1222       19342 :   return 1;
    1223             : }
    1224             : 
    1225             : long
    1226      105694 : ZC_is_ei(GEN x)
    1227             : {
    1228      105694 :   long i, j = 0, l = lg(x);
    1229     1001705 :   for (i = 1; i < l; i++)
    1230             :   {
    1231      896012 :     GEN c = gel(x,i);
    1232      896012 :     long s = signe(c);
    1233      896012 :     if (!s) continue;
    1234      105681 :     if (s < 0 || !is_pm1(c) || j) return 0;
    1235      105680 :     j = i;
    1236             :   }
    1237      105693 :   return j;
    1238             : }
    1239             : 
    1240             : /********************************************************************/
    1241             : /**                                                                **/
    1242             : /**                       MISCELLANEOUS                            **/
    1243             : /**                                                                **/
    1244             : /********************************************************************/
    1245             : /* assume lg(x) = lg(y), x,y in Z^n */
    1246             : int
    1247      730289 : ZV_cmp(GEN x, GEN y)
    1248             : {
    1249      730289 :   long fl,i, lx = lg(x);
    1250     1444746 :   for (i=1; i<lx; i++)
    1251     1162269 :     if (( fl = cmpii(gel(x,i), gel(y,i)) )) return fl;
    1252      282477 :   return 0;
    1253             : }
    1254             : /* assume lg(x) = lg(y), x,y in Z^n */
    1255             : int
    1256        3459 : ZV_abscmp(GEN x, GEN y)
    1257             : {
    1258        3459 :   long fl,i, lx = lg(x);
    1259       18761 :   for (i=1; i<lx; i++)
    1260       18681 :     if (( fl = abscmpii(gel(x,i), gel(y,i)) )) return fl;
    1261          80 :   return 0;
    1262             : }
    1263             : 
    1264             : long
    1265     3581538 : zv_content(GEN x)
    1266             : {
    1267     3581538 :   long i, s, l = lg(x);
    1268     3581538 :   if (l == 1) return 0;
    1269     3581531 :   s = labs(x[1]);
    1270     3581531 :   for (i = 2; i < l && s != 1; i++) s = ugcd(s, labs(x[i]));
    1271     3581531 :   return s;
    1272             : }
    1273             : GEN
    1274      194173 : ZV_content(GEN x)
    1275             : {
    1276      194173 :   long i, l = lg(x);
    1277      194173 :   pari_sp av = avma;
    1278             :   GEN c;
    1279      194173 :   if (l == 1) return gen_0;
    1280      194173 :   if (l == 2) return absi(gel(x,1));
    1281      129150 :   c = gel(x,1);
    1282      310940 :   for (i = 2; i < l; i++)
    1283             :   {
    1284      216601 :     c = gcdii(c, gel(x,i));
    1285      216601 :     if (is_pm1(c)) { set_avma(av); return gen_1; }
    1286             :   }
    1287       94339 :   return gerepileuptoint(av, c);
    1288             : }
    1289             : 
    1290             : GEN
    1291     1682020 : ZM_det_triangular(GEN mat)
    1292             : {
    1293             :   pari_sp av;
    1294     1682020 :   long i,l = lg(mat);
    1295             :   GEN s;
    1296             : 
    1297     1682020 :   if (l<3) return l<2? gen_1: icopy(gcoeff(mat,1,1));
    1298     1569000 :   av = avma; s = gcoeff(mat,1,1);
    1299     1569000 :   for (i=2; i<l; i++) s = mulii(s,gcoeff(mat,i,i));
    1300     1569000 :   return gerepileuptoint(av,s);
    1301             : }
    1302             : 
    1303             : /* assumes no overflow */
    1304             : long
    1305      380136 : zv_prod(GEN v)
    1306             : {
    1307      380136 :   long n, i, l = lg(v);
    1308      380136 :   if (l == 1) return 1;
    1309      286938 :   n = v[1]; for (i = 2; i < l; i++) n *= v[i];
    1310      286938 :   return n;
    1311             : }
    1312             : 
    1313             : static GEN
    1314   278298503 : _mulii(void *E, GEN a, GEN b)
    1315   278298503 : { (void) E; return mulii(a, b); }
    1316             : 
    1317             : /* product of ulongs */
    1318             : GEN
    1319      365084 : zv_prod_Z(GEN v)
    1320             : {
    1321      365084 :   pari_sp av = avma;
    1322      365084 :   long k, m, n = lg(v)-1;
    1323             :   GEN V;
    1324      365084 :   switch(n) {
    1325       19502 :     case 0: return gen_1;
    1326       56623 :     case 1: return utoi(v[1]);
    1327      125468 :     case 2: return muluu(v[1], v[2]);
    1328             :   }
    1329      163491 :   m = n >> 1;
    1330      163491 :   V = cgetg(m + (odd(n)? 2: 1), t_VEC);
    1331      163491 :   for (k = 1; k <= m; k++) gel(V,k) = muluu(v[k<<1], v[(k<<1)-1]);
    1332      163491 :   if (odd(n)) gel(V,k) = utoipos(v[n]);
    1333      163491 :   return gerepileuptoint(av, gen_product(V, NULL, &_mulii));
    1334             : }
    1335             : GEN
    1336    14694393 : vecsmall_prod(GEN v)
    1337             : {
    1338    14694393 :   pari_sp av = avma;
    1339    14694393 :   long k, m, n = lg(v)-1;
    1340             :   GEN V;
    1341    14694393 :   switch (n) {
    1342           0 :     case 0: return gen_1;
    1343           0 :     case 1: return stoi(v[1]);
    1344          21 :     case 2: return mulss(v[1], v[2]);
    1345             :   }
    1346    14694372 :   m = n >> 1;
    1347    14694372 :   V = cgetg(m + (odd(n)? 2: 1), t_VEC);
    1348    14694372 :   for (k = 1; k <= m; k++) gel(V,k) = mulss(v[k<<1], v[(k<<1)-1]);
    1349    14694372 :   if (odd(n)) gel(V,k) = stoi(v[n]);
    1350    14694372 :   return gerepileuptoint(av, gen_product(V, NULL, &_mulii));
    1351             : }
    1352             : 
    1353             : GEN
    1354     7530724 : ZV_prod(GEN v)
    1355             : {
    1356     7530724 :   pari_sp av = avma;
    1357     7530724 :   long i, l = lg(v);
    1358             :   GEN n;
    1359     7530724 :   if (l == 1) return gen_1;
    1360     7515282 :   if (l > 7) return gerepileuptoint(av, gen_product(v, NULL, _mulii));
    1361      174269 :   n = gel(v,1);
    1362      174269 :   if (l == 2) return icopy(n);
    1363       99971 :   for (i = 2; i < l; i++) n = mulii(n, gel(v,i));
    1364       99971 :   return gerepileuptoint(av, n);
    1365             : }
    1366             : /* assumes no overflow */
    1367             : long
    1368        1596 : zv_sum(GEN v)
    1369             : {
    1370        1596 :   long n, i, l = lg(v);
    1371        1596 :   if (l == 1) return 0;
    1372        1575 :   n = v[1]; for (i = 2; i < l; i++) n += v[i];
    1373        1575 :   return n;
    1374             : }
    1375             : GEN
    1376          56 : ZV_sum(GEN v)
    1377             : {
    1378          56 :   pari_sp av = avma;
    1379          56 :   long i, l = lg(v);
    1380             :   GEN n;
    1381          56 :   if (l == 1) return gen_0;
    1382          56 :   n = gel(v,1);
    1383          56 :   if (l == 2) return icopy(n);
    1384          56 :   for (i = 2; i < l; i++) n = addii(n, gel(v,i));
    1385          56 :   return gerepileuptoint(av, n);
    1386             : }
    1387             : 
    1388             : /********************************************************************/
    1389             : /**                                                                **/
    1390             : /**         GRAM SCHMIDT REDUCTION (integer matrices)              **/
    1391             : /**                                                                **/
    1392             : /********************************************************************/
    1393             : 
    1394             : /* L[k,] += q * L[l,], l < k. Inefficient if q = 0 */
    1395             : static void
    1396       75810 : Zupdate_row(long k, long l, GEN q, GEN L, GEN B)
    1397             : {
    1398       75810 :   long i, qq = itos_or_0(q);
    1399       75810 :   if (!qq)
    1400             :   {
    1401        2799 :     for(i=1;i<l;i++)  gcoeff(L,k,i) = addii(gcoeff(L,k,i),mulii(q,gcoeff(L,l,i)));
    1402        2799 :     gcoeff(L,k,l) = addii(gcoeff(L,k,l), mulii(q,B));
    1403        2799 :     return;
    1404             :   }
    1405       73011 :   if (qq == 1) {
    1406       19885 :     for (i=1;i<l; i++) gcoeff(L,k,i) = addii(gcoeff(L,k,i),gcoeff(L,l,i));
    1407       19885 :     gcoeff(L,k,l) = addii(gcoeff(L,k,l), B);
    1408       53126 :   } else if (qq == -1) {
    1409       18257 :     for (i=1;i<l; i++) gcoeff(L,k,i) = subii(gcoeff(L,k,i),gcoeff(L,l,i));
    1410       18257 :     gcoeff(L,k,l) = addii(gcoeff(L,k,l), negi(B));
    1411             :   } else {
    1412       34869 :     for(i=1;i<l;i++) gcoeff(L,k,i) = addii(gcoeff(L,k,i),mulsi(qq,gcoeff(L,l,i)));
    1413       34869 :     gcoeff(L,k,l) = addii(gcoeff(L,k,l), mulsi(qq,B));
    1414             :   }
    1415             : }
    1416             : 
    1417             : /* update L[k,] */
    1418             : static void
    1419      173611 : ZRED(long k, long l, GEN x, GEN L, GEN B)
    1420             : {
    1421      173611 :   GEN q = truedivii(addii(B,shifti(gcoeff(L,k,l),1)), shifti(B,1));
    1422      173611 :   if (!signe(q)) return;
    1423       75810 :   q = negi(q);
    1424       75810 :   Zupdate_row(k,l,q,L,B);
    1425       75810 :   gel(x,k) = ZC_lincomb(gen_1, q, gel(x,k), gel(x,l));
    1426             : }
    1427             : 
    1428             : /* Gram-Schmidt reduction, x a ZM */
    1429             : static void
    1430      202705 : ZincrementalGS(GEN x, GEN L, GEN B, long k)
    1431             : {
    1432             :   long i, j;
    1433      666642 :   for (j=1; j<=k; j++)
    1434             :   {
    1435      463937 :     pari_sp av = avma;
    1436      463937 :     GEN u = ZV_dotproduct(gel(x,k), gel(x,j));
    1437      974791 :     for (i=1; i<j; i++)
    1438             :     {
    1439      510854 :       u = subii(mulii(gel(B,i+1), u), mulii(gcoeff(L,k,i), gcoeff(L,j,i)));
    1440      510854 :       u = diviiexact(u, gel(B,i));
    1441             :     }
    1442      463937 :     gcoeff(L,k,j) = gerepileuptoint(av, u);
    1443             :   }
    1444      202705 :   gel(B,k+1) = gcoeff(L,k,k); gcoeff(L,k,k) = gen_1;
    1445      202705 : }
    1446             : 
    1447             : /* Variant reducemodinvertible(ZC v, ZM y), when y singular.
    1448             :  * Very inefficient if y is not LLL-reduced of maximal rank */
    1449             : static GEN
    1450       22493 : ZC_reducemodmatrix_i(GEN v, GEN y)
    1451             : {
    1452       22493 :   GEN B, L, x = shallowconcat(y, v);
    1453       22493 :   long k, lx = lg(x), nx = lx-1;
    1454             : 
    1455       22493 :   B = scalarcol_shallow(gen_1, lx);
    1456       22493 :   L = zeromatcopy(nx, nx);
    1457       22493 :   for (k=1; k <= nx; k++) ZincrementalGS(x, L, B, k);
    1458       22493 :   for (k = nx-1; k >= 1; k--) ZRED(nx,k, x,L,gel(B,k+1));
    1459       22493 :   return gel(x,nx);
    1460             : }
    1461             : GEN
    1462       22493 : ZC_reducemodmatrix(GEN v, GEN y) {
    1463       22493 :   pari_sp av = avma;
    1464       22493 :   return gerepilecopy(av, ZC_reducemodmatrix_i(v,y));
    1465             : }
    1466             : 
    1467             : /* Variant reducemodinvertible(ZM v, ZM y), when y singular.
    1468             :  * Very inefficient if y is not LLL-reduced of maximal rank */
    1469             : static GEN
    1470       36477 : ZM_reducemodmatrix_i(GEN v, GEN y)
    1471             : {
    1472             :   GEN B, L, V;
    1473       36477 :   long j, k, lv = lg(v), nx = lg(y), lx = nx+1;
    1474             : 
    1475       36477 :   V = cgetg(lv, t_MAT);
    1476       36477 :   B = scalarcol_shallow(gen_1, lx);
    1477       36477 :   L = zeromatcopy(nx, nx);
    1478       36477 :   for (k=1; k < nx; k++) ZincrementalGS(y, L, B, k);
    1479      105175 :   for (j = 1; j < lg(v); j++)
    1480             :   {
    1481       68698 :     GEN x = shallowconcat(y, gel(v,j));
    1482       68698 :     ZincrementalGS(x, L, B, nx); /* overwrite last */
    1483       68698 :     for (k = nx-1; k >= 1; k--) ZRED(nx,k, x,L,gel(B,k+1));
    1484       68698 :     gel(V,j) = gel(x,nx);
    1485             :   }
    1486       36477 :   return V;
    1487             : }
    1488             : GEN
    1489       36477 : ZM_reducemodmatrix(GEN v, GEN y) {
    1490       36477 :   pari_sp av = avma;
    1491       36477 :   return gerepilecopy(av, ZM_reducemodmatrix_i(v,y));
    1492             : }
    1493             : 
    1494             : GEN
    1495       12427 : ZC_reducemodlll(GEN x,GEN y)
    1496             : {
    1497       12427 :   pari_sp av = avma;
    1498       12427 :   GEN z = ZC_reducemodmatrix(x, ZM_lll(y, 0.75, LLL_INPLACE));
    1499       12427 :   return gerepilecopy(av, z);
    1500             : }
    1501             : GEN
    1502           0 : ZM_reducemodlll(GEN x,GEN y)
    1503             : {
    1504           0 :   pari_sp av = avma;
    1505           0 :   GEN z = ZM_reducemodmatrix(x, ZM_lll(y, 0.75, LLL_INPLACE));
    1506           0 :   return gerepilecopy(av, z);
    1507             : }

Generated by: LCOV version 1.13