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 - alglin2.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.13.0 lcov report (development 25825-edc109b529) Lines: 1018 1128 90.2 %
Date: 2020-09-21 06:08:33 Functions: 79 85 92.9 %
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             : /********************************************************************/
      15             : /**                                                                **/
      16             : /**                         LINEAR ALGEBRA                         **/
      17             : /**                         (second part)                          **/
      18             : /**                                                                **/
      19             : /********************************************************************/
      20             : #include "pari.h"
      21             : #include "paripriv.h"
      22             : /*******************************************************************/
      23             : /*                                                                 */
      24             : /*                   CHARACTERISTIC POLYNOMIAL                     */
      25             : /*                                                                 */
      26             : /*******************************************************************/
      27             : 
      28             : static GEN
      29       22908 : Flm_charpoly_i(GEN x, ulong p)
      30             : {
      31       22908 :   long lx = lg(x), r, i;
      32       22908 :   GEN H, y = cgetg(lx+1, t_VEC);
      33       22908 :   gel(y,1) = pol1_Flx(0); H = Flm_hess(x, p);
      34      148585 :   for (r = 1; r < lx; r++)
      35             :   {
      36      125678 :     pari_sp av2 = avma;
      37      125678 :     ulong a = 1;
      38      125678 :     GEN z, b = zero_Flx(0);
      39      375183 :     for (i = r-1; i; i--)
      40             :     {
      41      300128 :       a = Fl_mul(a, ucoeff(H,i+1,i), p);
      42      300099 :       if (!a) break;
      43      249518 :       b = Flx_add(b, Flx_Fl_mul(gel(y,i), Fl_mul(a,ucoeff(H,i,r),p), p), p);
      44             :     }
      45      125676 :     z = Flx_sub(Flx_shift(gel(y,r), 1),
      46      125636 :                 Flx_Fl_mul(gel(y,r), ucoeff(H,r,r), p), p);
      47             :     /* (X - H[r,r])y[r] - b */
      48      125678 :     gel(y,r+1) = gerepileuptoleaf(av2, Flx_sub(z, b, p));
      49             :   }
      50       22907 :   return gel(y,lx);
      51             : }
      52             : 
      53             : GEN
      54        1661 : Flm_charpoly(GEN x, ulong p)
      55             : {
      56        1661 :   pari_sp av = avma;
      57        1661 :   return gerepileuptoleaf(av, Flm_charpoly_i(x,p));
      58             : }
      59             : 
      60             : GEN
      61       19031 : FpM_charpoly(GEN x, GEN p)
      62             : {
      63       19031 :   pari_sp av = avma;
      64             :   long lx, r, i;
      65             :   GEN y, H;
      66             : 
      67       19031 :   if (lgefint(p) == 3)
      68             :   {
      69       18296 :     ulong pp = p[2];
      70       18296 :     y = Flx_to_ZX(Flm_charpoly_i(ZM_to_Flm(x,pp), pp));
      71       18296 :     return gerepileupto(av, y);
      72             :   }
      73         735 :   lx = lg(x); y = cgetg(lx+1, t_VEC);
      74         735 :   gel(y,1) = pol_1(0); H = FpM_hess(x, p);
      75        4046 :   for (r = 1; r < lx; r++)
      76             :   {
      77        4046 :     pari_sp av2 = avma;
      78        4046 :     GEN z, a = gen_1, b = pol_0(0);
      79        9226 :     for (i = r-1; i; i--)
      80             :     {
      81        7007 :       a = Fp_mul(a, gcoeff(H,i+1,i), p);
      82        7007 :       if (!signe(a)) break;
      83        5180 :       b = ZX_add(b, ZX_Z_mul(gel(y,i), Fp_mul(a,gcoeff(H,i,r),p)));
      84             :     }
      85        4046 :     b = FpX_red(b, p);
      86        4046 :     z = FpX_sub(RgX_shift_shallow(gel(y,r), 1),
      87        4046 :                 FpX_Fp_mul(gel(y,r), gcoeff(H,r,r), p), p);
      88        4046 :     z = FpX_sub(z,b,p);
      89        4046 :     if (r+1 == lx) { gel(y,lx) = z; break; }
      90        3311 :     gel(y,r+1) = gerepileupto(av2, z); /* (X - H[r,r])y[r] - b */
      91             :   }
      92         735 :   return gerepileupto(av, gel(y,lx));
      93             : }
      94             : 
      95             : GEN
      96         392 : charpoly0(GEN x, long v, long flag)
      97             : {
      98         392 :   if (v<0) v = 0;
      99         392 :   switch(flag)
     100             :   {
     101          14 :     case 0: return caradj(x,v,NULL);
     102          14 :     case 1: return caract(x,v);
     103          14 :     case 2: return carhess(x,v);
     104          14 :     case 3: return carberkowitz(x,v);
     105           7 :     case 4:
     106           7 :       if (typ(x) != t_MAT) pari_err_TYPE("charpoly",x);
     107           7 :       RgM_check_ZM(x, "charpoly");
     108           7 :       x = ZM_charpoly(x); setvarn(x, v); return x;
     109         329 :     case 5:
     110         329 :       return charpoly(x, v);
     111             :   }
     112           0 :   pari_err_FLAG("charpoly");
     113             :   return NULL; /* LCOV_EXCL_LINE */
     114             : }
     115             : 
     116             : /* characteristic pol. Easy cases. Return NULL in case it's not so easy. */
     117             : static GEN
     118        7690 : easychar(GEN x, long v)
     119             : {
     120             :   pari_sp av;
     121             :   long lx;
     122             :   GEN p1;
     123             : 
     124        7690 :   switch(typ(x))
     125             :   {
     126          35 :     case t_INT: case t_REAL: case t_INTMOD:
     127             :     case t_FRAC: case t_PADIC:
     128          35 :       p1=cgetg(4,t_POL);
     129          35 :       p1[1]=evalsigne(1) | evalvarn(v);
     130          35 :       gel(p1,2) = gneg(x); gel(p1,3) = gen_1;
     131          35 :       return p1;
     132             : 
     133          14 :     case t_COMPLEX: case t_QUAD:
     134          14 :       p1 = cgetg(5,t_POL);
     135          14 :       p1[1] = evalsigne(1) | evalvarn(v);
     136          14 :       gel(p1,2) = gnorm(x); av = avma;
     137          14 :       gel(p1,3) = gerepileupto(av, gneg(gtrace(x)));
     138          14 :       gel(p1,4) = gen_1; return p1;
     139             : 
     140          28 :     case t_FFELT: {
     141          28 :       pari_sp ltop=avma;
     142          28 :       p1 = FpX_to_mod(FF_charpoly(x), FF_p_i(x));
     143          28 :       setvarn(p1,v); return gerepileupto(ltop,p1);
     144             :     }
     145             : 
     146         840 :     case t_POLMOD:
     147             :     {
     148         840 :       GEN A = gel(x,2), T = gel(x,1);
     149         840 :       if (typ(A)==t_POL && RgX_is_QX(A) && RgX_is_ZX(T))
     150         812 :         return QXQ_charpoly(A, T, v);
     151             :       else
     152          28 :         return RgXQ_charpoly(A, T, v);
     153             :     }
     154        6773 :     case t_MAT:
     155        6773 :       lx=lg(x);
     156        6773 :       if (lx==1) return pol_1(v);
     157        6717 :       if (lgcols(x) != lx) break;
     158        6710 :       return NULL;
     159             :   }
     160           7 :   pari_err_TYPE("easychar",x);
     161             :   return NULL; /* LCOV_EXCL_LINE */
     162             : }
     163             : /* compute charpoly by mapping to Fp first, return lift to Z */
     164             : static GEN
     165          35 : RgM_Fp_charpoly(GEN x, GEN p, long v)
     166             : {
     167             :   GEN T;
     168          35 :   if (lgefint(p) == 3)
     169             :   {
     170          21 :     ulong pp = itou(p);
     171          21 :     T = Flm_charpoly_i(RgM_to_Flm(x, pp), pp);
     172          21 :     T = Flx_to_ZX(T);
     173             :   }
     174             :   else
     175          14 :     T = FpM_charpoly(RgM_to_FpM(x, p), p);
     176          35 :   setvarn(T, v); return T;
     177             : }
     178             : GEN
     179        2911 : charpoly(GEN x, long v)
     180             : {
     181        2911 :   GEN T, p = NULL;
     182             :   long prec;
     183        2911 :   if ((T = easychar(x,v))) return T;
     184        2771 :   switch(RgM_type(x, &p,&T,&prec))
     185             :   {
     186        1910 :     case t_INT:
     187        1910 :       T = ZM_charpoly(x); setvarn(T, v); break;
     188          35 :     case t_INTMOD:
     189          35 :       if (!BPSW_psp(p)) T = carberkowitz(x, v);
     190             :       else
     191             :       {
     192          35 :         pari_sp av = avma;
     193          35 :         T = RgM_Fp_charpoly(x,p,v);
     194          35 :         T = gerepileupto(av, FpX_to_mod(T,p));
     195             :       }
     196          35 :       break;
     197          49 :     case t_REAL:
     198             :     case t_COMPLEX:
     199             :     case t_PADIC:
     200          49 :       T = carhess(x, v);
     201          49 :       break;
     202         777 :     default:
     203         777 :       T = carberkowitz(x, v);
     204             :   }
     205        2764 :   return T;
     206             : }
     207             : 
     208             : /* We possibly worked with an "invalid" polynomial p, satisfying
     209             :  * varn(p) > gvar2(p). Fix this. */
     210             : static GEN
     211        3876 : fix_pol(pari_sp av, GEN p)
     212             : {
     213        3876 :   long w = gvar2(p), v = varn(p);
     214        3876 :   if (w == v) pari_err_PRIORITY("charpoly", p, "=", w);
     215        3869 :   if (varncmp(w,v) < 0) p = gerepileupto(av, poleval(p, pol_x(v)));
     216        3869 :   return p;
     217             : }
     218             : GEN
     219          14 : caract(GEN x, long v)
     220             : {
     221          14 :   pari_sp av = avma;
     222             :   GEN  T, C, x_k, Q;
     223             :   long k, n;
     224             : 
     225          14 :   if ((T = easychar(x,v))) return T;
     226             : 
     227          14 :   n = lg(x)-1;
     228          14 :   if (n == 1) return fix_pol(av, deg1pol(gen_1, gneg(gcoeff(x,1,1)), v));
     229             : 
     230          14 :   x_k = pol_x(v); /* to be modified in place */
     231          14 :   T = scalarpol(det(x), v); C = utoineg(n); Q = pol_x(v);
     232          28 :   for (k=1; k<=n; k++)
     233             :   {
     234          28 :     GEN mk = utoineg(k), d;
     235          28 :     gel(x_k,2) = mk;
     236          28 :     d = det(RgM_Rg_add_shallow(x, mk));
     237          28 :     T = RgX_add(RgX_mul(T, x_k), RgX_Rg_mul(Q, gmul(C, d)));
     238          28 :     if (k == n) break;
     239             : 
     240          14 :     Q = RgX_mul(Q, x_k);
     241          14 :     C = diviuexact(mulsi(k-n,C), k+1); /* (-1)^k binomial(n,k) */
     242             :   }
     243          14 :   return fix_pol(av, RgX_Rg_div(T, mpfact(n)));
     244             : }
     245             : 
     246             : /* C = charpoly(x, v) */
     247             : static GEN
     248          21 : RgM_adj_from_char(GEN x, long v, GEN C)
     249             : {
     250          21 :   if (varn(C) != v) /* problem with variable priorities */
     251             :   {
     252           7 :     C = gdiv(gsub(C, gsubst(C, v, gen_0)), pol_x(v));
     253           7 :     if (odd(lg(x))) C = RgX_neg(C); /* even dimension */
     254           7 :     return gsubst(C, v, x);
     255             :   }
     256             :   else
     257             :   {
     258          14 :     C = RgX_shift_shallow(C, -1);
     259          14 :     if (odd(lg(x))) C = RgX_neg(C); /* even dimension */
     260          14 :     return RgX_RgM_eval(C, x);
     261             :   }
     262             : }
     263             : /* assume x square matrice */
     264             : static GEN
     265        6359 : mattrace(GEN x)
     266             : {
     267        6359 :   long i, lx = lg(x);
     268             :   GEN t;
     269        6359 :   if (lx < 3) return lx == 1? gen_0: gcopy(gcoeff(x,1,1));
     270        6289 :   t = gcoeff(x,1,1);
     271       17747 :   for (i = 2; i < lx; i++) t = gadd(t, gcoeff(x,i,i));
     272        6289 :   return t;
     273             : }
     274             : static int
     275        2322 : bad_char(GEN q, long n)
     276             : {
     277             :   forprime_t S;
     278             :   ulong p;
     279        2322 :   if (!signe(q)) return 0;
     280          42 :   (void)u_forprime_init(&S, 2, n);
     281          98 :   while ((p = u_forprime_next(&S)))
     282          70 :     if (!umodiu(q, p)) return 1;
     283          28 :   return 0;
     284             : }
     285             : /* Using traces: return the characteristic polynomial of x (in variable v).
     286             :  * If py != NULL, the adjoint matrix is put there. */
     287             : GEN
     288        2385 : caradj(GEN x, long v, GEN *py)
     289             : {
     290             :   pari_sp av, av0;
     291             :   long i, k, n;
     292             :   GEN T, y, t;
     293             : 
     294        2385 :   if ((T = easychar(x, v)))
     295             :   {
     296          42 :     if (py)
     297             :     {
     298          42 :       if (typ(x) != t_MAT) pari_err_TYPE("matadjoint",x);
     299          42 :       *py = cgetg(1,t_MAT);
     300             :     }
     301          42 :     return T;
     302             :   }
     303             : 
     304        2343 :   n = lg(x)-1; av0 = avma;
     305        2343 :   T = cgetg(n+3,t_POL); T[1] = evalsigne(1) | evalvarn(v);
     306        2343 :   gel(T,n+2) = gen_1;
     307        2343 :   if (!n) { if (py) *py = cgetg(1,t_MAT); return T; }
     308        2343 :   av = avma; t = gerepileupto(av, gneg(mattrace(x)));
     309        2343 :   gel(T,n+1) = t;
     310        2343 :   if (n == 1) {
     311           7 :     T = fix_pol(av0, T);
     312           7 :     if (py) *py = matid(1); return T;
     313             :   }
     314        2336 :   if (n == 2) {
     315          14 :     GEN a = gcoeff(x,1,1), b = gcoeff(x,1,2);
     316          14 :     GEN c = gcoeff(x,2,1), d = gcoeff(x,2,2);
     317          14 :     av = avma;
     318          14 :     gel(T,2) = gerepileupto(av, gsub(gmul(a,d), gmul(b,c)));
     319          14 :     T = fix_pol(av0, T);
     320          14 :     if (py) {
     321           7 :       y = cgetg(3, t_MAT);
     322           7 :       gel(y,1) = mkcol2(gcopy(d), gneg(c));
     323           7 :       gel(y,2) = mkcol2(gneg(b), gcopy(a));
     324           7 :       *py = y;
     325             :     }
     326          14 :     return T;
     327             :   }
     328             :   /* l > 3 */
     329        2322 :   if (bad_char(residual_characteristic(x), n))
     330             :   { /* n! not invertible in base ring */
     331          14 :     T = charpoly(x, v);
     332          14 :     if (!py) return gerepileupto(av, T);
     333          14 :     *py = RgM_adj_from_char(x, v, T);
     334          14 :     gerepileall(av, 2, &T,py);
     335          14 :     return T;
     336             :   }
     337        2308 :   av = avma; y = RgM_shallowcopy(x);
     338        9239 :   for (i = 1; i <= n; i++) gcoeff(y,i,i) = gadd(gcoeff(y,i,i), t);
     339        4623 :   for (k = 2; k < n; k++)
     340             :   {
     341        2315 :     GEN y0 = y;
     342        2315 :     y = RgM_mul(y, x);
     343        2315 :     t = gdivgs(mattrace(y), -k);
     344        9274 :     for (i = 1; i <= n; i++) gcoeff(y,i,i) = gadd(gcoeff(y,i,i), t);
     345        2315 :     y = gclone(y);
     346        2315 :     gel(T,n-k+2) = gerepilecopy(av, t); av = avma;
     347        2315 :     if (k > 2) gunclone(y0);
     348             :   }
     349        2308 :   t = gmul(gcoeff(x,1,1),gcoeff(y,1,1));
     350        6931 :   for (i=2; i<=n; i++) t = gadd(t, gmul(gcoeff(x,1,i),gcoeff(y,i,1)));
     351        2308 :   gel(T,2) = gerepileupto(av, gneg(t));
     352        2308 :   T = fix_pol(av0, T);
     353        2308 :   if (py) *py = odd(n)? gcopy(y): RgM_neg(y);
     354        2308 :   gunclone(y); return T;
     355             : }
     356             : 
     357             : GEN
     358        2371 : adj(GEN x)
     359             : {
     360             :   GEN y;
     361        2371 :   (void)caradj(x, fetch_var(), &y);
     362        2371 :   (void)delete_var(); return y;
     363             : }
     364             : 
     365             : GEN
     366           7 : adjsafe(GEN x)
     367             : {
     368           7 :   const long v = fetch_var();
     369           7 :   pari_sp av = avma;
     370             :   GEN C, A;
     371           7 :   if (typ(x) != t_MAT) pari_err_TYPE("matadjoint",x);
     372           7 :   if (lg(x) < 3) return gcopy(x);
     373           7 :   C = charpoly(x,v);
     374           7 :   A = RgM_adj_from_char(x, v, C);
     375           7 :   (void)delete_var(); return gerepileupto(av, A);
     376             : }
     377             : 
     378             : GEN
     379         112 : matadjoint0(GEN x, long flag)
     380             : {
     381         112 :   switch(flag)
     382             :   {
     383         105 :     case 0: return adj(x);
     384           7 :     case 1: return adjsafe(x);
     385             :   }
     386           0 :   pari_err_FLAG("matadjoint");
     387             :   return NULL; /* LCOV_EXCL_LINE */
     388             : }
     389             : 
     390             : /*******************************************************************/
     391             : /*                                                                 */
     392             : /*                       Frobenius form                            */
     393             : /*                                                                 */
     394             : /*******************************************************************/
     395             : 
     396             : /* The following section implement a mix of Ozello and Storjohann algorithms
     397             : 
     398             : P. Ozello, doctoral thesis (in French):
     399             : Calcul exact des formes de Jordan et de Frobenius d'une matrice, Chapitre 2
     400             : http://tel.archives-ouvertes.fr/tel-00323705
     401             : 
     402             : A. Storjohann,  Diss. ETH No. 13922
     403             : Algorithms for Matrix Canonical Forms, Chapter 9
     404             : https://cs.uwaterloo.ca/~astorjoh/diss2up.pdf
     405             : 
     406             : We use Storjohann Lemma 9.14 (step1, step2, step3) Ozello theorem 4,
     407             : and Storjohann Lemma 9.18
     408             : */
     409             : 
     410             : /* Elementary transforms */
     411             : 
     412             : /* M <- U^(-1) M U, U = E_{i,j}(k) => U^(-1) = E{i,j}(-k)
     413             :  * P = U * P */
     414             : static void
     415       14182 : transL(GEN M, GEN P, GEN k, long i, long j)
     416             : {
     417       14182 :   long l, n = lg(M)-1;
     418      173628 :   for(l=1; l<=n; l++) /* M[,j]-=k*M[,i] */
     419      159446 :     gcoeff(M,l,j) = gsub(gcoeff(M,l,j), gmul(gcoeff(M,l,i), k));
     420      173628 :   for(l=1; l<=n; l++) /* M[i,]+=k*M[j,] */
     421      159446 :     gcoeff(M,i,l) = gadd(gcoeff(M,i,l), gmul(gcoeff(M,j,l), k));
     422       14182 :   if (P)
     423      163653 :     for(l=1; l<=n; l++)
     424      150570 :       gcoeff(P,i,l) = gadd(gcoeff(P,i,l), gmul(gcoeff(P,j,l), k));
     425       14182 : }
     426             : 
     427             : /* j = a or b */
     428             : static void
     429        2261 : transD(GEN M, GEN P, long a, long b, long j)
     430             : {
     431             :   long l, n;
     432        2261 :   GEN k = gcoeff(M,a,b), ki;
     433             : 
     434        2261 :   if (gequal1(k)) return;
     435        1127 :   ki = ginv(k); n = lg(M)-1;
     436       11802 :   for(l=1; l<=n; l++)
     437       10675 :     if (l!=j)
     438             :     {
     439        9548 :       gcoeff(M,l,j) = gmul(gcoeff(M,l,j), k);
     440        9548 :       gcoeff(M,j,l) = (j==a && l==b)? gen_1: gmul(gcoeff(M,j,l), ki);
     441             :     }
     442        1127 :   if (P)
     443        9457 :     for(l=1; l<=n; l++)
     444        8603 :       gcoeff(P,j,l) = gmul(gcoeff(P,j,l), ki);
     445             : }
     446             : 
     447             : static void
     448         378 : transS(GEN M, GEN P, long i, long j)
     449             : {
     450         378 :   long l, n = lg(M)-1;
     451         378 :   swap(gel(M,i), gel(M,j));
     452        5495 :   for (l=1; l<=n; l++)
     453        5117 :     swap(gcoeff(M,i,l), gcoeff(M,j,l));
     454         378 :   if (P)
     455        4060 :     for (l=1; l<=n; l++)
     456        3843 :       swap(gcoeff(P,i,l), gcoeff(P,j,l));
     457         378 : }
     458             : 
     459             : /* Convert companion matrix to polynomial*/
     460             : static GEN
     461         280 : minpoly_polslice(GEN M, long i, long j, long v)
     462             : {
     463         280 :   long k, d = j+1-i;
     464         280 :   GEN P = cgetg(d+3,t_POL);
     465         280 :   P[1] = evalsigne(1)|evalvarn(v);
     466        1379 :   for (k=0; k<d; k++)
     467        1099 :     gel(P,k+2) = gneg(gcoeff(M,i+k, j));
     468         280 :   gel(P,d+2) = gen_1;
     469         280 :   return P;
     470             : }
     471             : 
     472             : static GEN
     473          49 : minpoly_listpolslice(GEN M, GEN V, long v)
     474             : {
     475          49 :   long i, n = lg(M)-1, nb = lg(V)-1;
     476          49 :   GEN W = cgetg(nb+1, t_VEC);
     477         147 :   for (i=1; i<=nb; i++)
     478          98 :     gel(W,i) = minpoly_polslice(M, V[i], i < nb? V[i+1]-1: n, v);
     479          49 :   return W;
     480             : }
     481             : 
     482             : static int
     483          91 : minpoly_dvdslice(GEN M, long i, long j, long k)
     484             : {
     485          91 :   pari_sp av = avma;
     486          91 :   long r = signe(RgX_rem(minpoly_polslice(M, i, j-1, 0),
     487             :                         minpoly_polslice(M, j, k, 0)));
     488          91 :   return gc_bool(av, r == 0);
     489             : }
     490             : 
     491             : static void
     492           0 : RgM_replace(GEN M, GEN M2)
     493             : {
     494           0 :   long n = lg(M)-1, m = nbrows(M), i, j;
     495           0 :   for(i=1; i<=n; i++)
     496           0 :     for(j=1; j<=m; j++)
     497           0 :       gcoeff(M, i, j) = gcoeff(M2, i, j);
     498           0 : }
     499             : 
     500             : static void
     501           0 : gerepilemat2_inplace(pari_sp av, GEN M, GEN P)
     502             : {
     503           0 :   GEN M2 = M, P2 = P;
     504           0 :   gerepileall(av, P ? 2: 1, &M2, &P2);
     505           0 :   RgM_replace(M, M2);
     506           0 :   if (P) RgM_replace(P, P2);
     507           0 : }
     508             : 
     509             : /* Lemma 9.14 */
     510             : static long
     511         581 : weakfrobenius_step1(GEN M, GEN P, long j0)
     512             : {
     513         581 :   pari_sp av = avma;
     514         581 :   long n = lg(M)-1, k, j;
     515        2821 :   for (j = j0; j < n; ++j)
     516             :   {
     517        2401 :     if (gequal0(gcoeff(M, j+1, j)))
     518             :     {
     519        1456 :       for (k = j+2; k <= n; ++k)
     520        1295 :         if (!gequal0(gcoeff(M,k,j))) break;
     521         518 :       if (k > n) return j;
     522         357 :       transS(M, P, k, j+1);
     523             :     }
     524        2240 :     transD(M, P, j+1, j, j+1);
     525             :     /* Now M[j+1,j] = 1 */
     526       25452 :     for (k = 1; k <= n; ++k)
     527       23212 :       if (k != j+1 && !gequal0(gcoeff(M,k,j))) /* zero M[k,j] */
     528             :       {
     529       13559 :         transL(M, P, gneg(gcoeff(M,k,j)), k, j+1);
     530       13559 :         gcoeff(M,k,j) = gen_0; /* avoid approximate 0 */
     531             :       }
     532        2240 :     if (gc_needed(av,1))
     533             :     {
     534           0 :       if (DEBUGMEM > 1)
     535           0 :         pari_warn(warnmem,"RgM_minpoly stage 1: j0=%ld, j=%ld", j0, j);
     536           0 :       gerepilemat2_inplace(av, M, P);
     537             :     }
     538             :   }
     539         420 :   return n;
     540             : }
     541             : 
     542             : static void
     543         581 : weakfrobenius_step2(GEN M, GEN P, long j)
     544             : {
     545         581 :   pari_sp av = avma;
     546         581 :   long i, k, n = lg(M)-1;
     547        3437 :   for(i=j; i>=2; i--)
     548             :   {
     549        6237 :     for(k=j+1; k<=n; k++)
     550        3381 :       if (!gequal0(gcoeff(M,i,k)))
     551         623 :         transL(M, P, gcoeff(M,i,k), i-1, k);
     552        2856 :     if (gc_needed(av,1))
     553             :     {
     554           0 :       if (DEBUGMEM > 1)
     555           0 :         pari_warn(warnmem,"RgM_minpoly stage 2: j=%ld, i=%ld", j, i);
     556           0 :       gerepilemat2_inplace(av, M, P);
     557             :     }
     558             :   }
     559         581 : }
     560             : 
     561             : static long
     562         581 : weakfrobenius_step3(GEN M, GEN P, long j0, long j)
     563             : {
     564         581 :   long i, k, n = lg(M)-1;
     565         581 :   if (j == n) return 0;
     566         161 :   if (gequal0(gcoeff(M, j0, j+1)))
     567             :   {
     568         588 :     for (k=j+2; k<=n; k++)
     569         448 :       if (!gequal0(gcoeff(M, j0, k))) break;
     570         140 :     if (k > n) return 0;
     571           0 :     transS(M, P, k, j+1);
     572             :   }
     573          21 :   transD(M, P, j0, j+1, j+1);
     574          21 :   for (i=j+2; i<=n; i++)
     575           0 :     if (!gequal0(gcoeff(M, j0, i)))
     576           0 :       transL(M, P, gcoeff(M, j0, i),j+1, i);
     577          21 :   return 1;
     578             : }
     579             : 
     580             : /* flag: 0 -> full Frobenius from , 1 -> weak Frobenius form */
     581             : static GEN
     582         420 : RgM_Frobenius(GEN M, long flag, GEN *pt_P, GEN *pt_v)
     583             : {
     584         420 :   pari_sp av = avma, av2, ltop;
     585         420 :   long n = lg(M)-1, eps, j0 = 1, nb = 0;
     586             :   GEN v, P;
     587         420 :   v = cgetg(n+1, t_VECSMALL);
     588         420 :   ltop = avma;
     589         420 :   P = pt_P ? matid(n): NULL;
     590         420 :   M = RgM_shallowcopy(M);
     591         420 :   av2 = avma;
     592        1001 :   while (j0 <= n)
     593             :   {
     594         581 :     long j = weakfrobenius_step1(M, P, j0);
     595         581 :     weakfrobenius_step2(M, P, j);
     596         581 :     eps = weakfrobenius_step3(M, P, j0, j);
     597         581 :     if (eps == 0)
     598             :     {
     599         560 :       v[++nb] = j0;
     600         560 :       if (flag == 0 && nb > 1 && !minpoly_dvdslice(M, v[nb-1], j0, j))
     601             :       {
     602           0 :         j = j0; j0 = v[nb-1]; nb -= 2;
     603           0 :         transL(M, P, gen_1, j, j0); /*lemma 9.18*/
     604             :       } else
     605         560 :         j0 = j+1;
     606             :     }
     607             :     else
     608          21 :       transS(M, P, j0, j+1); /*theorem 4*/
     609         581 :     if (gc_needed(av,1))
     610             :     {
     611           0 :       if (DEBUGMEM > 1)
     612           0 :         pari_warn(warnmem,"weakfrobenius j0=%ld",j0);
     613           0 :       gerepilemat2_inplace(av2, M, P);
     614             :     }
     615             :   }
     616         420 :   fixlg(v, nb+1);
     617         420 :   if (pt_v) *pt_v = v;
     618         420 :   gerepileall(pt_v ? ltop: av, P? 2: 1, &M, &P);
     619         420 :   if (pt_P) *pt_P = P;
     620         420 :   return M;
     621             : }
     622             : 
     623             : static GEN
     624          49 : RgM_minpoly(GEN M, long v)
     625             : {
     626          49 :   pari_sp av = avma;
     627             :   GEN V, W;
     628          49 :   if (lg(M) == 1) return pol_1(v);
     629          49 :   M = RgM_Frobenius(M, 1, NULL, &V);
     630          49 :   W = minpoly_listpolslice(M, V, v);
     631          49 :   if (varncmp(v,gvar2(W)) >= 0)
     632           0 :     pari_err_PRIORITY("matfrobenius", M, "<=", v);
     633          49 :   return gerepileupto(av, RgX_normalize(glcm0(W, NULL)));
     634             : }
     635             : 
     636             : GEN
     637           0 : Frobeniusform(GEN V, long n)
     638             : {
     639             :   long i, j, k;
     640           0 :   GEN M = zeromatcopy(n,n);
     641           0 :   for (k=1,i=1;i<lg(V);i++,k++)
     642             :   {
     643           0 :     GEN  P = gel(V,i);
     644           0 :     long d = degpol(P);
     645           0 :     if (k+d-1 > n) pari_err_PREC("matfrobenius");
     646           0 :     for (j=0; j<d-1; j++, k++) gcoeff(M,k+1,k) = gen_1;
     647           0 :     for (j=0; j<d; j++) gcoeff(M,k-j,k) = gneg(gel(P, 1+d-j));
     648             :   }
     649           0 :   return M;
     650             : }
     651             : 
     652             : GEN
     653         371 : matfrobenius(GEN M, long flag, long v)
     654             : {
     655             :   long n;
     656         371 :   if (typ(M)!=t_MAT) pari_err_TYPE("matfrobenius",M);
     657         371 :   if (v < 0) v = 0;
     658         371 :   n = lg(M)-1;
     659         371 :   if (n && lgcols(M)!=n+1) pari_err_DIM("matfrobenius");
     660         371 :   if (flag > 2) pari_err_FLAG("matfrobenius");
     661         371 :   switch (flag)
     662             :   {
     663           7 :   case 0:
     664           7 :     return RgM_Frobenius(M, 0, NULL, NULL);
     665           0 :   case 1:
     666             :     {
     667           0 :       pari_sp av = avma;
     668             :       GEN V, W, F;
     669           0 :       F = RgM_Frobenius(M, 0, NULL, &V);
     670           0 :       W = minpoly_listpolslice(F, V, v);
     671           0 :       if (varncmp(v, gvar2(W)) >= 0)
     672           0 :         pari_err_PRIORITY("matfrobenius", M, "<=", v);
     673           0 :       return gerepileupto(av, W);
     674             :     }
     675         364 :   case 2:
     676             :     {
     677         364 :       GEN P, F, R = cgetg(3, t_VEC);
     678         364 :       F = RgM_Frobenius(M, 0, &P, NULL);
     679         364 :       gel(R,1) = F; gel(R,2) = P;
     680         364 :       return R;
     681             :     }
     682           0 :   default:
     683           0 :     pari_err_FLAG("matfrobenius");
     684             :   }
     685             :   return NULL; /*LCOV_EXCL_LINE*/
     686             : }
     687             : 
     688             : /*******************************************************************/
     689             : /*                                                                 */
     690             : /*                       MINIMAL POLYNOMIAL                        */
     691             : /*                                                                 */
     692             : /*******************************************************************/
     693             : 
     694             : static GEN
     695          63 : RgXQ_minpoly_naive(GEN y, GEN P)
     696             : {
     697          63 :   long n = lgpol(P);
     698          63 :   GEN M = ker(RgXQ_matrix_pow(y,n,n,P));
     699          63 :   return content(RgM_to_RgXV(M,varn(P)));
     700             : }
     701             : 
     702             : static GEN
     703           0 : RgXQ_minpoly_FpXQ(GEN x, GEN y, GEN p)
     704             : {
     705           0 :   pari_sp av = avma;
     706             :   GEN r;
     707           0 :   if (lgefint(p) == 3)
     708             :   {
     709           0 :     ulong pp = uel(p, 2);
     710           0 :     r = Flx_to_ZX_inplace(Flxq_minpoly(RgX_to_Flx(x, pp),
     711             :                                    RgX_to_Flx(y, pp), pp));
     712             :   }
     713             :   else
     714           0 :     r = FpXQ_minpoly(RgX_to_FpX(x, p), RgX_to_FpX(y, p), p);
     715           0 :   return gerepileupto(av, FpX_to_mod(r, p));
     716             : }
     717             : 
     718             : static GEN
     719          21 : RgXQ_minpoly_FpXQXQ(GEN x, GEN S, GEN pol, GEN p)
     720             : {
     721          21 :   pari_sp av = avma;
     722             :   GEN r;
     723          21 :   GEN T = RgX_to_FpX(pol, p);
     724          21 :   if (signe(T)==0) pari_err_OP("minpoly",x,S);
     725          21 :   if (lgefint(p) == 3)
     726             :   {
     727          13 :     ulong pp = uel(p, 2);
     728          13 :     GEN Tp = ZX_to_Flx(T, pp);
     729          13 :     r = FlxX_to_ZXX(FlxqXQ_minpoly(RgX_to_FlxqX(x, Tp, pp),
     730             :                                    RgX_to_FlxqX(S, Tp, pp), Tp, pp));
     731             :   }
     732             :   else
     733           8 :     r = FpXQXQ_minpoly(RgX_to_FpXQX(x, T, p), RgX_to_FpXQX(S, T, p), T, p);
     734          21 :   return gerepileupto(av, FpXQX_to_mod(r, T, p));
     735             : }
     736             : 
     737             : #define code(t1,t2) ((t1 << 6) | t2)
     738             : 
     739             : static GEN
     740         959 : RgXQ_minpoly_fast(GEN x, GEN y)
     741             : {
     742             :   GEN p, pol;
     743             :   long pa;
     744         959 :   long t = RgX_type2(x,y, &p,&pol,&pa);
     745         959 :   switch(t)
     746             :   {
     747           0 :     case t_INTMOD: return RgXQ_minpoly_FpXQ(x, y, p);
     748          77 :     case t_FFELT:  return FFXQ_minpoly(x, y, pol);
     749          21 :     case code(t_POLMOD, t_INTMOD):
     750          21 :                    return RgXQ_minpoly_FpXQXQ(x, y, pol, p);
     751         861 :     default:       return NULL;
     752             :   }
     753             : }
     754             : 
     755             : #undef code
     756             : 
     757             : static GEN
     758        1239 : easymin(GEN x, long v)
     759             : {
     760             :   GEN G, R, dR;
     761        1239 :   long tx = typ(x);
     762        1239 :   if (tx == t_FFELT)
     763             :   {
     764         154 :     R = FpX_to_mod(FF_minpoly(x), FF_p_i(x));
     765         154 :     setvarn(R,v); return R;
     766             :   }
     767        1085 :   if (tx == t_POLMOD)
     768             :   {
     769        1036 :     GEN a = gel(x,2), b = gel(x,1);
     770        1036 :     if (typ(a) != t_POL || varn(a) != varn(b))
     771             :     {
     772          77 :       if (varncmp(gvar(a), v) <= 0) pari_err_PRIORITY("minpoly", x, "<", v);
     773          70 :       return deg1pol(gen_1, gneg_i(a), v);
     774             :     }
     775         959 :     R = RgXQ_minpoly_fast(a,b);
     776         959 :     if (R) { setvarn(R, v); return R; }
     777         861 :     if (!issquarefree(b))
     778             :     {
     779          63 :       R = RgXQ_minpoly_naive(a, b);
     780          63 :       setvarn(R,v); return R;
     781             :     }
     782             :   }
     783         847 :   R = easychar(x, v); if (!R) return NULL;
     784         798 :   dR = RgX_deriv(R);  if (!lgpol(dR)) return NULL;
     785         798 :   G = RgX_normalize(RgX_gcd(R,dR));
     786         798 :   return RgX_div(R,G);
     787             : }
     788             : GEN
     789        1239 : minpoly(GEN x, long v)
     790             : {
     791        1239 :   pari_sp av = avma;
     792             :   GEN P;
     793        1239 :   if (v < 0) v = 0;
     794        1239 :   P = easymin(x,v);
     795        1232 :   if (P) return gerepileupto(av,P);
     796             :   /* typ(x) = t_MAT */
     797          49 :   set_avma(av); return RgM_minpoly(x,v);
     798             : }
     799             : 
     800             : /*******************************************************************/
     801             : /*                                                                 */
     802             : /*                       HESSENBERG FORM                           */
     803             : /*                                                                 */
     804             : /*******************************************************************/
     805             : static int
     806         462 : relative0(GEN a, GEN a0, long bit)
     807         462 : { return gequal0(a) || (bit && gexpo(a)-gexpo(a0) < bit); }
     808             : /* assume x a non-empty square t_MAT */
     809             : static GEN
     810          70 : RgM_hess(GEN x0, long prec)
     811             : {
     812             :   pari_sp av;
     813          70 :   long lx = lg(x0), bit = prec? 8 - bit_accuracy(prec): 0, m, i, j;
     814             :   GEN x;
     815             : 
     816          70 :   if (bit) x0 = RgM_shallowcopy(x0);
     817          70 :   av = avma; x = RgM_shallowcopy(x0);
     818         287 :   for (m=2; m<lx-1; m++)
     819             :   {
     820         217 :     GEN t = NULL;
     821         560 :     for (i=m+1; i<lx; i++)
     822             :     {
     823         462 :       t = gcoeff(x,i,m-1);
     824         462 :       if (!relative0(t, gcoeff(x0,i,m-1), bit)) break;
     825             :     }
     826         217 :     if (i == lx) continue;
     827         665 :     for (j=m-1; j<lx; j++) swap(gcoeff(x,i,j), gcoeff(x,m,j));
     828         119 :     swap(gel(x,i), gel(x,m));
     829         119 :     if (bit)
     830             :     {
     831         329 :       for (j=m-1; j<lx; j++) swap(gcoeff(x0,i,j), gcoeff(x0,m,j));
     832          63 :       swap(gel(x0,i), gel(x0,m));
     833             :     }
     834         119 :     t = ginv(t);
     835             : 
     836         427 :     for (i=m+1; i<lx; i++)
     837             :     {
     838         308 :       GEN c = gcoeff(x,i,m-1);
     839         308 :       if (gequal0(c)) continue;
     840             : 
     841         266 :       c = gmul(c,t); gcoeff(x,i,m-1) = gen_0;
     842        1330 :       for (j=m; j<lx; j++)
     843        1064 :         gcoeff(x,i,j) = gsub(gcoeff(x,i,j), gmul(c,gcoeff(x,m,j)));
     844        1981 :       for (j=1; j<lx; j++)
     845        1715 :         gcoeff(x,j,m) = gadd(gcoeff(x,j,m), gmul(c,gcoeff(x,j,i)));
     846         266 :       if (gc_needed(av,2))
     847             :       {
     848           0 :         if (DEBUGMEM>1) pari_warn(warnmem,"hess, m = %ld", m);
     849           0 :         gerepileall(av,2, &x, &t);
     850             :       }
     851             :     }
     852             :   }
     853          70 :   return x;
     854             : }
     855             : 
     856             : GEN
     857          70 : hess(GEN x)
     858             : {
     859          70 :   pari_sp av = avma;
     860          70 :   GEN p = NULL, T = NULL;
     861          70 :   long prec, lx = lg(x);
     862          70 :   if (typ(x) != t_MAT) pari_err_TYPE("hess",x);
     863          70 :   if (lx == 1) return cgetg(1,t_MAT);
     864          70 :   if (lgcols(x) != lx) pari_err_DIM("hess");
     865          70 :   switch(RgM_type(x, &p, &T, &prec))
     866             :   {
     867          42 :     case t_REAL:
     868          42 :     case t_COMPLEX: break;
     869          28 :     default: prec = 0;
     870             :   }
     871          70 :   return gerepilecopy(av, RgM_hess(x,prec));
     872             : }
     873             : 
     874             : GEN
     875       22908 : Flm_hess(GEN x, ulong p)
     876             : {
     877       22908 :   long lx = lg(x), m, i, j;
     878       22908 :   if (lx == 1) return cgetg(1,t_MAT);
     879       22908 :   if (lgcols(x) != lx) pari_err_DIM("hess");
     880             : 
     881       22908 :   x = Flm_copy(x);
     882      103024 :   for (m=2; m<lx-1; m++)
     883             :   {
     884       80489 :     ulong t = 0;
     885      283977 :     for (i=m+1; i<lx; i++) { t = ucoeff(x,i,m-1); if (t) break; }
     886       80489 :     if (i == lx) continue;
     887      452553 :     for (j=m-1; j<lx; j++) lswap(ucoeff(x,i,j), ucoeff(x,m,j));
     888       50881 :     swap(gel(x,i), gel(x,m)); t = Fl_inv(t, p);
     889             : 
     890      349790 :     for (i=m+1; i<lx; i++)
     891             :     {
     892      299282 :       ulong c = ucoeff(x,i,m-1);
     893      299282 :       if (!c) continue;
     894             : 
     895      174384 :       c = Fl_mul(c,t,p); ucoeff(x,i,m-1) = 0;
     896     2012573 :       for (j=m; j<lx; j++)
     897     1838246 :         ucoeff(x,i,j) = Fl_sub(ucoeff(x,i,j), Fl_mul(c,ucoeff(x,m,j), p), p);
     898     3154663 :       for (j=1; j<lx; j++)
     899     2980666 :         ucoeff(x,j,m) = Fl_add(ucoeff(x,j,m), Fl_mul(c,ucoeff(x,j,i), p), p);
     900             :     }
     901             :   }
     902       22535 :   return x;
     903             : }
     904             : GEN
     905         735 : FpM_hess(GEN x, GEN p)
     906             : {
     907         735 :   pari_sp av = avma;
     908         735 :   long lx = lg(x), m, i, j;
     909         735 :   if (lx == 1) return cgetg(1,t_MAT);
     910         735 :   if (lgcols(x) != lx) pari_err_DIM("hess");
     911         735 :   if (lgefint(p) == 3)
     912             :   {
     913           0 :     ulong pp = p[2];
     914           0 :     x = Flm_hess(ZM_to_Flm(x, pp), pp);
     915           0 :     return gerepileupto(av, Flm_to_ZM(x));
     916             :   }
     917         735 :   x = RgM_shallowcopy(x);
     918        3311 :   for (m=2; m<lx-1; m++)
     919             :   {
     920        2576 :     GEN t = NULL;
     921        5481 :     for (i=m+1; i<lx; i++) { t = gcoeff(x,i,m-1); if (signe(t)) break; }
     922        2576 :     if (i == lx) continue;
     923       12019 :     for (j=m-1; j<lx; j++) swap(gcoeff(x,i,j), gcoeff(x,m,j));
     924        1890 :     swap(gel(x,i), gel(x,m)); t = Fp_inv(t, p);
     925             : 
     926        8239 :     for (i=m+1; i<lx; i++)
     927             :     {
     928        6349 :       GEN c = gcoeff(x,i,m-1);
     929        6349 :       if (!signe(c)) continue;
     930             : 
     931        5061 :       c = Fp_mul(c,t, p); gcoeff(x,i,m-1) = gen_0;
     932       29785 :       for (j=m; j<lx; j++)
     933       24724 :         gcoeff(x,i,j) = Fp_sub(gcoeff(x,i,j), Fp_mul(c,gcoeff(x,m,j),p), p);
     934       45465 :       for (j=1; j<lx; j++)
     935       40404 :         gcoeff(x,j,m) = Fp_add(gcoeff(x,j,m), Fp_mul(c,gcoeff(x,j,i),p), p);
     936        5061 :       if (gc_needed(av,2))
     937             :       {
     938           0 :         if (DEBUGMEM>1) pari_warn(warnmem,"hess, m = %ld", m);
     939           0 :         gerepileall(av,2, &x, &t);
     940             :       }
     941             :     }
     942             :   }
     943         735 :   return gerepilecopy(av,x);
     944             : }
     945             : /* H in Hessenberg form */
     946             : static GEN
     947          63 : RgM_hess_charpoly(GEN H, long v)
     948             : {
     949          63 :   long n = lg(H), r, i;
     950          63 :   GEN y = cgetg(n+1, t_VEC);
     951          63 :   gel(y,1) = pol_1(v);
     952         371 :   for (r = 1; r < n; r++)
     953             :   {
     954         308 :     pari_sp av2 = avma;
     955         308 :     GEN z, a = gen_1, b = pol_0(v);
     956         812 :     for (i = r-1; i; i--)
     957             :     {
     958         581 :       a = gmul(a, gcoeff(H,i+1,i));
     959         581 :       if (gequal0(a)) break;
     960         504 :       b = RgX_add(b, RgX_Rg_mul(gel(y,i), gmul(a,gcoeff(H,i,r))));
     961             :     }
     962         308 :     z = RgX_sub(RgX_shift_shallow(gel(y,r), 1),
     963         308 :                 RgX_Rg_mul(gel(y,r), gcoeff(H,r,r)));
     964         308 :     gel(y,r+1) = gerepileupto(av2, RgX_sub(z, b)); /* (X - H[r,r])y[r] - b */
     965             :   }
     966          63 :   return gel(y,n);
     967             : }
     968             : GEN
     969          63 : carhess(GEN x, long v)
     970             : {
     971             :   pari_sp av;
     972             :   GEN y;
     973          63 :   if ((y = easychar(x,v))) return y;
     974          63 :   av = avma; y = RgM_hess_charpoly(hess(x), v);
     975          63 :   return fix_pol(av, y);
     976             : }
     977             : 
     978             : /* Bound for sup norm of charpoly(M/dM), M integral: let B = |M|oo / |dM|,
     979             :  *   s = max_k binomial(n,k) (kB^2)^(k/2),
     980             :  * return ceil(log2(s)) */
     981             : static long
     982        3156 : charpoly_bound(GEN M, GEN dM, GEN N)
     983             : {
     984        3156 :   pari_sp av = avma;
     985        3156 :   GEN B = itor(N, LOWDEFAULTPREC);
     986        3156 :   GEN s = real_0(LOWDEFAULTPREC), bin, B2;
     987        3156 :   long n = lg(M)-1, k;
     988        3156 :   bin = gen_1;
     989        3156 :   if (dM) B = divri(B, dM);
     990        3156 :   B2 = sqrr(B);
     991       14017 :   for (k = n; k >= (n+1)>>1; k--)
     992             :   {
     993       10861 :     GEN t = mulri(powruhalf(mulur(k, B2), k), bin);
     994       10861 :     if (abscmprr(t, s) > 0) s = t;
     995       10861 :     bin = diviuexact(muliu(bin, k), n-k+1);
     996             :   }
     997        3156 :   return gc_long(av, ceil(dbllog2(s)));
     998             : }
     999             : 
    1000             : static GEN
    1001        3522 : QM_charpoly_ZX_slice(GEN A, GEN dM, GEN P, GEN *mod)
    1002             : {
    1003        3522 :   pari_sp av = avma;
    1004        3522 :   long i, n = lg(P)-1;
    1005             :   GEN H, T;
    1006        3522 :   if (n == 1)
    1007             :   {
    1008        2930 :     ulong p = uel(P,1), dp = dM ? umodiu(dM, p): 1;
    1009        2930 :     GEN Hp, a = ZM_to_Flm(A, p);
    1010        2930 :     Hp = Flm_charpoly_i(a, p);
    1011        2930 :     if (dp != 1) Hp = Flx_rescale(Hp, Fl_inv(dp, p), p);
    1012        2930 :     Hp = gerepileupto(av, Flx_to_ZX(Hp));
    1013        2930 :     *mod = utoi(p); return Hp;
    1014             :   }
    1015         592 :   T = ZV_producttree(P);
    1016         592 :   A = ZM_nv_mod_tree(A, P, T);
    1017         592 :   H = cgetg(n+1, t_VEC);
    1018        2253 :   for(i=1; i <= n; i++)
    1019             :   {
    1020        1661 :     ulong p = uel(P,i), dp = dM ? umodiu(dM, p): 1;
    1021        1661 :     gel(H,i) = Flm_charpoly(gel(A, i), p);
    1022        1661 :     if (dp != 1) gel(H,i) = Flx_rescale(gel(H,i), Fl_inv(dp, p), p);
    1023             :   }
    1024         592 :   H = nxV_chinese_center_tree(H, P, T, ZV_chinesetree(P,T));
    1025         592 :   *mod = gmael(T, lg(T)-1, 1);
    1026         592 :   gerepileall(av, 2, &H, mod);
    1027         592 :   return H;
    1028             : }
    1029             : 
    1030             : GEN
    1031        3522 : QM_charpoly_ZX_worker(GEN P, GEN M, GEN dM)
    1032             : {
    1033        3522 :   GEN V = cgetg(3, t_VEC);
    1034        3522 :   gel(V,1) = QM_charpoly_ZX_slice(M, equali1(dM) ? NULL:dM, P, &gel(V,2));
    1035        3522 :   return V;
    1036             : }
    1037             : 
    1038             : /* Assume M a square ZM, dM integer. Return charpoly(M / dM) in Z[X] */
    1039             : static GEN
    1040        3660 : QM_charpoly_ZX_i(GEN M, GEN dM, long bound)
    1041             : {
    1042        3660 :   long n = lg(M)-1;
    1043             :   forprime_t S;
    1044        3660 :   GEN worker = snm_closure(is_entry("_QM_charpoly_ZX_worker"),
    1045             :                            mkvec2(M, dM? dM: gen_1));
    1046        3660 :   if (!n) return pol_1(0);
    1047        3660 :   if (bound < 0)
    1048             :   {
    1049        3394 :     GEN N = ZM_supnorm(M);
    1050        3394 :     if (signe(N) == 0) return monomial(gen_1, n, 0);
    1051        3156 :     bound = charpoly_bound(M, dM, N) + 1;
    1052             :   }
    1053        3422 :   if (DEBUGLEVEL>5) err_printf("ZM_charpoly: bound 2^%ld\n", bound);
    1054        3422 :   init_modular_big(&S);
    1055        3422 :   return gen_crt("QM_charpoly_ZX", worker, &S, dM, bound, 0, NULL,
    1056             :                  nxV_chinese_center, FpX_center);
    1057             : }
    1058             : 
    1059             : GEN
    1060         266 : QM_charpoly_ZX_bound(GEN M, long bit)
    1061             : {
    1062         266 :   pari_sp av = avma;
    1063         266 :   GEN dM; M = Q_remove_denom(M, &dM);
    1064         266 :   return gerepilecopy(av, QM_charpoly_ZX_i(M, dM, bit));
    1065             : }
    1066             : GEN
    1067        1477 : QM_charpoly_ZX(GEN M)
    1068             : {
    1069        1477 :   pari_sp av = avma;
    1070        1477 :   GEN dM; M = Q_remove_denom(M, &dM);
    1071        1477 :   return gerepilecopy(av, QM_charpoly_ZX_i(M, dM, -1));
    1072             : }
    1073             : GEN
    1074        1917 : ZM_charpoly(GEN M)
    1075             : {
    1076        1917 :   pari_sp av = avma;
    1077        1917 :   return gerepilecopy(av, QM_charpoly_ZX_i(M, NULL, -1));
    1078             : }
    1079             : 
    1080             : /*******************************************************************/
    1081             : /*                                                                 */
    1082             : /*        CHARACTERISTIC POLYNOMIAL (BERKOWITZ'S ALGORITHM)        */
    1083             : /*                                                                 */
    1084             : /*******************************************************************/
    1085             : GEN
    1086        1470 : carberkowitz(GEN x, long v)
    1087             : {
    1088             :   long lx, i, j, k, r;
    1089             :   GEN V, S, C, Q;
    1090             :   pari_sp av0, av;
    1091        1470 :   if ((V = easychar(x,v))) return V;
    1092        1470 :   lx = lg(x); av0 = avma;
    1093        1470 :   V = cgetg(lx+1, t_VEC);
    1094        1470 :   S = cgetg(lx+1, t_VEC);
    1095        1470 :   C = cgetg(lx+1, t_VEC);
    1096        1470 :   Q = cgetg(lx+1, t_VEC);
    1097        1470 :   av = avma;
    1098        1470 :   gel(C,1) = gen_m1;
    1099        1470 :   gel(V,1) = gen_m1;
    1100       11431 :   for (i=2;i<=lx; i++) gel(C,i) = gel(Q,i) = gel(S,i) = gel(V,i) = gen_0;
    1101        1470 :   gel(V,2) = gcoeff(x,1,1);
    1102        9961 :   for (r = 2; r < lx; r++)
    1103             :   {
    1104             :     pari_sp av2;
    1105             :     GEN t;
    1106             : 
    1107       68628 :     for (i = 1; i < r; i++) gel(S,i) = gcoeff(x,i,r);
    1108        8491 :     gel(C,2) = gcoeff(x,r,r);
    1109       60137 :     for (i = 1; i < r-1; i++)
    1110             :     {
    1111       51646 :       av2 = avma; t = gmul(gcoeff(x,r,1), gel(S,1));
    1112      759304 :       for (j = 2; j < r; j++) t = gadd(t, gmul(gcoeff(x,r,j), gel(S,j)));
    1113       51646 :       gel(C,i+2) = gerepileupto(av2, t);
    1114      810950 :       for (j = 1; j < r; j++)
    1115             :       {
    1116      759304 :         av2 = avma; t = gmul(gcoeff(x,j,1), gel(S,1));
    1117    15028454 :         for (k = 2; k < r; k++) t = gadd(t, gmul(gcoeff(x,j,k), gel(S,k)));
    1118      759304 :         gel(Q,j) = gerepileupto(av2, t);
    1119             :       }
    1120      810950 :       for (j = 1; j < r; j++) gel(S,j) = gel(Q,j);
    1121             :     }
    1122        8491 :     av2 = avma; t = gmul(gcoeff(x,r,1), gel(S,1));
    1123       60137 :     for (j = 2; j < r; j++) t = gadd(t, gmul(gcoeff(x,r,j), gel(S,j)));
    1124        8491 :     gel(C,r+1) = gerepileupto(av2, t);
    1125        8491 :     if (gc_needed(av0,1))
    1126             :     {
    1127           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"carberkowitz");
    1128           0 :       gerepileall(av, 2, &C, &V);
    1129             :     }
    1130       85610 :     for (i = 1; i <= r+1; i++)
    1131             :     {
    1132       77119 :       av2 = avma; t = gmul(gel(C,i), gel(V,1));
    1133      577045 :       for (j = 2; j <= minss(r,i); j++)
    1134      499926 :         t = gadd(t, gmul(gel(C,i+1-j), gel(V,j)));
    1135       77119 :       gel(Q,i) = gerepileupto(av2, t);
    1136             :     }
    1137       85610 :     for (i = 1; i <= r+1; i++) gel(V,i) = gel(Q,i);
    1138             :   }
    1139        1470 :   V = RgV_to_RgX(vecreverse(V), v); /* not gtopoly: fail if v > gvar(V) */
    1140        1470 :   V = odd(lx)? gcopy(V): RgX_neg(V);
    1141        1470 :   return fix_pol(av0, V);
    1142             : }
    1143             : 
    1144             : /*******************************************************************/
    1145             : /*                                                                 */
    1146             : /*                            NORMS                                */
    1147             : /*                                                                 */
    1148             : /*******************************************************************/
    1149             : GEN
    1150      805355 : gnorm(GEN x)
    1151             : {
    1152             :   pari_sp av;
    1153             :   long lx, i;
    1154             :   GEN y;
    1155             : 
    1156      805355 :   switch(typ(x))
    1157             :   {
    1158       17840 :     case t_INT:  return sqri(x);
    1159      180853 :     case t_REAL: return sqrr(x);
    1160        1078 :     case t_FRAC: return sqrfrac(x);
    1161      601125 :     case t_COMPLEX: av = avma; return gerepileupto(av, cxnorm(x));
    1162         315 :     case t_QUAD:    av = avma; return gerepileupto(av, quadnorm(x));
    1163             : 
    1164          14 :     case t_POL: case t_SER: case t_RFRAC: av = avma;
    1165          14 :       return gerepileupto(av, greal(gmul(conj_i(x),x)));
    1166             : 
    1167          28 :     case t_FFELT:
    1168          28 :       y = cgetg(3, t_INTMOD);
    1169          28 :       gel(y,1) = FF_p(x);
    1170          28 :       gel(y,2) = FF_norm(x); return y;
    1171             : 
    1172        4102 :     case t_POLMOD:
    1173             :     {
    1174        4102 :       GEN T = gel(x,1), a = gel(x,2);
    1175        4102 :       if (typ(a) != t_POL || varn(a) != varn(T)) return gpowgs(a, degpol(T));
    1176        3934 :       return RgXQ_norm(a, T);
    1177             :     }
    1178           0 :     case t_VEC: case t_COL: case t_MAT:
    1179           0 :       y = cgetg_copy(x, &lx);
    1180           0 :       for (i=1; i<lx; i++) gel(y,i) = gnorm(gel(x,i));
    1181           0 :       return y;
    1182             :   }
    1183           0 :   pari_err_TYPE("gnorm",x);
    1184             :   return NULL; /* LCOV_EXCL_LINE */
    1185             : }
    1186             : 
    1187             : /* return |q|^2, complex modulus */
    1188             : static GEN
    1189          28 : cxquadnorm(GEN q, long prec)
    1190             : {
    1191          28 :   GEN X = gel(q,1), c = gel(X,2); /* (1-D)/4, -D/4 */
    1192          28 :   if (signe(c) > 0) return quadnorm(q); /* imaginary */
    1193          21 :   if (!prec) pari_err_TYPE("gnorml2", q);
    1194           7 :   return sqrr(quadtofp(q, prec));
    1195             : }
    1196             : 
    1197             : static GEN
    1198    26350695 : gnorml2_i(GEN x, long prec)
    1199             : {
    1200             :   pari_sp av;
    1201             :   long i, lx;
    1202             :   GEN s;
    1203             : 
    1204    26350695 :   switch(typ(x))
    1205             :   {
    1206     7562041 :     case t_INT:  return sqri(x);
    1207    14550162 :     case t_REAL: return sqrr(x);
    1208           7 :     case t_FRAC: return sqrfrac(x);
    1209     1039446 :     case t_COMPLEX: av = avma; return gerepileupto(av, cxnorm(x));
    1210          21 :     case t_QUAD:    av = avma; return gerepileupto(av, cxquadnorm(x,prec));
    1211             : 
    1212       11739 :     case t_POL: lx = lg(x)-1; x++; break;
    1213             : 
    1214     3187285 :     case t_VEC:
    1215             :     case t_COL:
    1216     3187285 :     case t_MAT: lx = lg(x); break;
    1217             : 
    1218           0 :     default: pari_err_TYPE("gnorml2",x);
    1219             :       return NULL; /* LCOV_EXCL_LINE */
    1220             :   }
    1221     3199024 :   if (lx == 1) return gen_0;
    1222     3199024 :   av = avma;
    1223     3199024 :   s = gnorml2(gel(x,1));
    1224    23123981 :   for (i=2; i<lx; i++)
    1225             :   {
    1226    19924964 :     s = gadd(s, gnorml2(gel(x,i)));
    1227    19924960 :     if (gc_needed(av,1))
    1228             :     {
    1229           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"gnorml2");
    1230           0 :       s = gerepileupto(av, s);
    1231             :     }
    1232             :   }
    1233     3199017 :   return gerepileupto(av,s);
    1234             : }
    1235             : GEN
    1236    26350653 : gnorml2(GEN x) { return gnorml2_i(x, 0); }
    1237             : 
    1238             : static GEN pnormlp(GEN,GEN,long);
    1239             : static GEN
    1240          63 : pnormlpvec(long i0, GEN x, GEN p, long prec)
    1241             : {
    1242          63 :   pari_sp av = avma;
    1243          63 :   long i, lx = lg(x);
    1244          63 :   GEN s = gen_0;
    1245         224 :   for (i=i0; i<lx; i++)
    1246             :   {
    1247         161 :     s = gadd(s, pnormlp(gel(x,i),p,prec));
    1248         161 :     if (gc_needed(av,1))
    1249             :     {
    1250           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"gnormlp, i = %ld", i);
    1251           0 :       s = gerepileupto(av, s);
    1252             :     }
    1253             :   }
    1254          63 :   return s;
    1255             : }
    1256             : /* (||x||_p)^p */
    1257             : static GEN
    1258         196 : pnormlp(GEN x, GEN p, long prec)
    1259             : {
    1260         196 :   switch(typ(x))
    1261             :   {
    1262         119 :     case t_INT: case t_REAL: x = mpabs(x); break;
    1263           0 :     case t_FRAC: x = absfrac(x); break;
    1264          14 :     case t_COMPLEX: case t_QUAD: x = gabs(x,prec); break;
    1265           7 :     case t_POL: return pnormlpvec(2, x, p, prec);
    1266          56 :     case t_VEC: case t_COL: case t_MAT: return pnormlpvec(1, x, p, prec);
    1267           0 :     default: pari_err_TYPE("gnormlp",x);
    1268             :   }
    1269         133 :   return gpow(x, p, prec);
    1270             : }
    1271             : 
    1272             : GEN
    1273         847 : gnormlp(GEN x, GEN p, long prec)
    1274             : {
    1275         847 :   pari_sp av = avma;
    1276         847 :   if (!p || (typ(p) == t_INFINITY && inf_get_sign(p) > 0))
    1277         686 :     return gsupnorm(x, prec);
    1278         161 :   if (gsigne(p) <= 0) pari_err_DOMAIN("normlp", "p", "<=", gen_0, p);
    1279         154 :   if (is_scalar_t(typ(x))) return gabs(x, prec);
    1280          91 :   if (typ(p) == t_INT)
    1281             :   {
    1282          63 :     ulong pp = itou_or_0(p);
    1283          63 :     switch(pp)
    1284             :     {
    1285          28 :       case 1: return gnorml1(x, prec);
    1286          28 :       case 2: x = gnorml2_i(x, prec); break;
    1287           7 :       default: x = pnormlp(x, p, prec); break;
    1288             :     }
    1289          35 :     if (pp && typ(x) == t_INT && Z_ispowerall(x, pp, &x))
    1290           7 :       return gerepileuptoleaf(av, x);
    1291          28 :     if (pp == 2) return gerepileupto(av, gsqrt(x, prec));
    1292             :   }
    1293             :   else
    1294          28 :     x = pnormlp(x, p, prec);
    1295          28 :   x = gpow(x, ginv(p), prec);
    1296          28 :   return gerepileupto(av, x);
    1297             : }
    1298             : 
    1299             : GEN
    1300        7077 : gnorml1(GEN x,long prec)
    1301             : {
    1302        7077 :   pari_sp av = avma;
    1303             :   long lx,i;
    1304             :   GEN s;
    1305        7077 :   switch(typ(x))
    1306             :   {
    1307        5005 :     case t_INT: case t_REAL: return mpabs(x);
    1308           0 :     case t_FRAC: return absfrac(x);
    1309             : 
    1310         462 :     case t_COMPLEX: case t_QUAD:
    1311         462 :       return gabs(x,prec);
    1312             : 
    1313        1561 :     case t_POL:
    1314        1561 :       lx = lg(x); s = gen_0;
    1315        6937 :       for (i=2; i<lx; i++) s = gadd(s, gnorml1(gel(x,i),prec));
    1316        1561 :       break;
    1317             : 
    1318          49 :     case t_VEC: case t_COL: case t_MAT:
    1319          49 :       lx = lg(x); s = gen_0;
    1320         168 :       for (i=1; i<lx; i++) s = gadd(s, gnorml1(gel(x,i),prec));
    1321          49 :       break;
    1322             : 
    1323           0 :     default: pari_err_TYPE("gnorml1",x);
    1324             :       return NULL; /* LCOV_EXCL_LINE */
    1325             :   }
    1326        1610 :   return gerepileupto(av, s);
    1327             : }
    1328             : /* As gnorml1, except for t_QUAD and t_COMPLEX: |x + wy| := |x| + |y|
    1329             :  * Still a norm of R-vector spaces, and can be cheaply computed without
    1330             :  * square roots */
    1331             : GEN
    1332      345065 : gnorml1_fake(GEN x)
    1333             : {
    1334      345065 :   pari_sp av = avma;
    1335             :   long lx, i;
    1336             :   GEN s;
    1337      345065 :   switch(typ(x))
    1338             :   {
    1339      256158 :     case t_INT: case t_REAL: return mpabs(x);
    1340           0 :     case t_FRAC: return absfrac(x);
    1341             : 
    1342           0 :     case t_COMPLEX:
    1343           0 :       s = gadd(gnorml1_fake(gel(x,1)), gnorml1_fake(gel(x,2)));
    1344           0 :       break;
    1345           0 :     case t_QUAD:
    1346           0 :       s = gadd(gnorml1_fake(gel(x,2)), gnorml1_fake(gel(x,3)));
    1347           0 :       break;
    1348             : 
    1349       88907 :     case t_POL:
    1350       88907 :       lx = lg(x); s = gen_0;
    1351      305823 :       for (i=2; i<lx; i++) s = gadd(s, gnorml1_fake(gel(x,i)));
    1352       88907 :       break;
    1353             : 
    1354           0 :     case t_VEC: case t_COL: case t_MAT:
    1355           0 :       lx = lg(x); s = gen_0;
    1356           0 :       for (i=1; i<lx; i++) s = gadd(s, gnorml1_fake(gel(x,i)));
    1357           0 :       break;
    1358             : 
    1359           0 :     default: pari_err_TYPE("gnorml1_fake",x);
    1360             :       return NULL; /* LCOV_EXCL_LINE */
    1361             :   }
    1362       88907 :   return gerepileupto(av, s);
    1363             : }
    1364             : 
    1365             : static void
    1366       73430 : store(GEN z, GEN *m) { if (!*m || gcmp(z, *m) > 0) *m = z; }
    1367             : /* compare |x| to *m or |x|^2 to *msq, whichever is easiest, and update
    1368             :  * the pointed value if x is larger */
    1369             : void
    1370       83223 : gsupnorm_aux(GEN x, GEN *m, GEN *msq, long prec)
    1371             : {
    1372             :   long i, lx;
    1373             :   GEN z;
    1374       83223 :   switch(typ(x))
    1375             :   {
    1376        9079 :     case t_COMPLEX: z = cxnorm(x); store(z, msq); return;
    1377           7 :     case t_QUAD:  z = cxquadnorm(x,prec); store(z, msq); return;
    1378       64344 :     case t_INT: case t_REAL: z = mpabs(x); store(z,m); return;
    1379           0 :     case t_FRAC: z = absfrac(x); store(z,m); return;
    1380             : 
    1381        5677 :     case t_POL: lx = lg(x)-1; x++; break;
    1382             : 
    1383        4116 :     case t_VEC:
    1384             :     case t_COL:
    1385        4116 :     case t_MAT: lx = lg(x); break;
    1386             : 
    1387           0 :     default: pari_err_TYPE("gsupnorm",x);
    1388             :       return; /* LCOV_EXCL_LINE */
    1389             :   }
    1390       83230 :   for (i=1; i<lx; i++) gsupnorm_aux(gel(x,i), m, msq, prec);
    1391             : }
    1392             : GEN
    1393        9786 : gsupnorm(GEN x, long prec)
    1394             : {
    1395        9786 :   GEN m = NULL, msq = NULL;
    1396        9786 :   pari_sp av = avma;
    1397        9786 :   gsupnorm_aux(x, &m, &msq, prec);
    1398             :   /* now set m = max (m, sqrt(msq)) */
    1399        9786 :   if (msq) {
    1400        1239 :     msq = gsqrt(msq, prec);
    1401        1239 :     if (!m || gcmp(m, msq) < 0) m = msq;
    1402        8547 :   } else if (!m) m = gen_0;
    1403        9786 :   return gerepilecopy(av, m);
    1404             : }
    1405             : 
    1406             : /*******************************************************************/
    1407             : /*                                                                 */
    1408             : /*                            TRACES                               */
    1409             : /*                                                                 */
    1410             : /*******************************************************************/
    1411             : GEN
    1412          35 : matcompanion(GEN x)
    1413             : {
    1414          35 :   long n = degpol(x), j;
    1415             :   GEN y, c;
    1416             : 
    1417          35 :   if (typ(x)!=t_POL) pari_err_TYPE("matcompanion",x);
    1418          35 :   if (!signe(x)) pari_err_DOMAIN("matcompanion","polynomial","=",gen_0,x);
    1419          28 :   if (n == 0) return cgetg(1, t_MAT);
    1420             : 
    1421          28 :   y = cgetg(n+1,t_MAT);
    1422         105 :   for (j=1; j < n; j++) gel(y,j) = col_ei(n, j+1);
    1423          28 :   c = cgetg(n+1,t_COL); gel(y,n) = c;
    1424          28 :   if (gequal1(gel(x, n+2)))
    1425         112 :     for (j=1; j<=n; j++) gel(c,j) = gneg(gel(x,j+1));
    1426             :   else
    1427             :   { /* not monic. Hardly ever used */
    1428           7 :     pari_sp av = avma;
    1429           7 :     GEN d = gclone(gneg(gel(x,n+2)));
    1430           7 :     set_avma(av);
    1431          21 :     for (j=1; j<=n; j++) gel(c,j) = gdiv(gel(x,j+1), d);
    1432           7 :     gunclone(d);
    1433             :   }
    1434          28 :   return y;
    1435             : }
    1436             : 
    1437             : GEN
    1438      272551 : gtrace(GEN x)
    1439             : {
    1440             :   pari_sp av;
    1441      272551 :   long i, lx, tx = typ(x);
    1442             :   GEN y, z;
    1443             : 
    1444      272551 :   switch(tx)
    1445             :   {
    1446        1822 :     case t_INT: case t_REAL: case t_FRAC:
    1447        1822 :       return gmul2n(x,1);
    1448             : 
    1449      268307 :     case t_COMPLEX:
    1450      268307 :       return gmul2n(gel(x,1),1);
    1451             : 
    1452          28 :     case t_QUAD:
    1453          28 :       y = gel(x,1);
    1454          28 :       if (!gequal0(gel(y,3)))
    1455             :       { /* assume quad. polynomial is either x^2 + d or x^2 - x + d */
    1456          28 :         av = avma;
    1457          28 :         return gerepileupto(av, gadd(gel(x,3), gmul2n(gel(x,2),1)));
    1458             :       }
    1459           0 :       return gmul2n(gel(x,2),1);
    1460             : 
    1461           7 :     case t_POL:
    1462           7 :       y = cgetg_copy(x, &lx); y[1] = x[1];
    1463          21 :       for (i=2; i<lx; i++) gel(y,i) = gtrace(gel(x,i));
    1464           7 :       return normalizepol_lg(y, lx);
    1465             : 
    1466          14 :     case t_SER:
    1467          14 :       if (ser_isexactzero(x)) return gcopy(x);
    1468           7 :       y = cgetg_copy(x, &lx); y[1] = x[1];
    1469          21 :       for (i=2; i<lx; i++) gel(y,i) = gtrace(gel(x,i));
    1470           7 :       return normalize(y);
    1471             : 
    1472         623 :     case t_POLMOD:
    1473         623 :       y = gel(x,1); z = gel(x,2);
    1474         623 :       if (typ(z) != t_POL || varn(y) != varn(z)) return gmulsg(degpol(y), z);
    1475         406 :       av = avma;
    1476         406 :       return gerepileupto(av, quicktrace(z, polsym(y, degpol(y)-1)));
    1477             : 
    1478          28 :     case t_FFELT:
    1479          28 :       y=cgetg(3, t_INTMOD);
    1480          28 :       gel(y,1) = FF_p(x);
    1481          28 :       gel(y,2) = FF_trace(x);
    1482          28 :       return y;
    1483             : 
    1484           7 :     case t_RFRAC:
    1485           7 :       av = avma; return gerepileupto(av, gadd(x, conj_i(x)));
    1486             : 
    1487           0 :     case t_VEC: case t_COL:
    1488           0 :       y = cgetg_copy(x, &lx);
    1489           0 :       for (i=1; i<lx; i++) gel(y,i) = gtrace(gel(x,i));
    1490           0 :       return y;
    1491             : 
    1492        1715 :     case t_MAT:
    1493        1715 :       lx = lg(x); if (lx == 1) return gen_0;
    1494             :       /*now lx >= 2*/
    1495        1708 :       if (lx != lgcols(x)) pari_err_DIM("gtrace");
    1496        1701 :       av = avma; return gerepileupto(av, mattrace(x));
    1497             :   }
    1498           0 :   pari_err_TYPE("gtrace",x);
    1499             :   return NULL; /* LCOV_EXCL_LINE */
    1500             : }
    1501             : 
    1502             : /* Cholesky decomposition for positive definite matrix a
    1503             :  * [GTM138, Algo 2.7.6, matrix Q]
    1504             :  * If a is not positive definite return NULL. */
    1505             : GEN
    1506        1363 : qfgaussred_positive(GEN a)
    1507             : {
    1508        1363 :   pari_sp av = avma;
    1509             :   GEN b;
    1510        1363 :   long i,j,k, n = lg(a);
    1511             : 
    1512        1363 :   if (typ(a)!=t_MAT) pari_err_TYPE("qfgaussred_positive",a);
    1513        1363 :   if (n == 1) return cgetg(1, t_MAT);
    1514        1356 :   if (lgcols(a)!=n) pari_err_DIM("qfgaussred_positive");
    1515        1356 :   b = cgetg(n,t_MAT);
    1516        8288 :   for (j=1; j<n; j++)
    1517             :   {
    1518        6932 :     GEN p1=cgetg(n,t_COL), p2=gel(a,j);
    1519             : 
    1520        6932 :     gel(b,j) = p1;
    1521       41349 :     for (i=1; i<=j; i++) gel(p1,i) = gel(p2,i);
    1522       34417 :     for (   ; i<n ; i++) gel(p1,i) = gen_0;
    1523             :   }
    1524        8288 :   for (k=1; k<n; k++)
    1525             :   {
    1526        6932 :     GEN bk, p = gcoeff(b,k,k), invp;
    1527        6932 :     if (gsigne(p)<=0) return gc_NULL(av); /* not positive definite */
    1528        6932 :     invp = ginv(p);
    1529        6932 :     bk = row(b, k);
    1530       34417 :     for (i=k+1; i<n; i++) gcoeff(b,k,i) = gmul(gel(bk,i), invp);
    1531       34417 :     for (i=k+1; i<n; i++)
    1532             :     {
    1533       27485 :       GEN c = gel(bk, i);
    1534      155433 :       for (j=i; j<n; j++)
    1535      127948 :         gcoeff(b,i,j) = gsub(gcoeff(b,i,j), gmul(c,gcoeff(b,k,j)));
    1536             :     }
    1537        6932 :     if (gc_needed(av,1))
    1538             :     {
    1539           0 :       if (DEBUGMEM>1) pari_warn(warnmem,"qfgaussred_positive");
    1540           0 :       b=gerepilecopy(av,b);
    1541             :     }
    1542             :   }
    1543        1356 :   return gerepilecopy(av,b);
    1544             : }
    1545             : 
    1546             : /* Maximal pivot strategy: x is a suitable pivot if it is non zero and either
    1547             :  * - an exact type, or
    1548             :  * - it is maximal among remaining non-zero (t_REAL) pivots */
    1549             : static int
    1550       23534 : suitable(GEN x, long k, GEN *pp, long *pi)
    1551             : {
    1552       23534 :   long t = typ(x);
    1553       23534 :   switch(t)
    1554             :   {
    1555        4522 :     case t_INT: return signe(x) != 0;
    1556       18858 :     case t_FRAC: return 1;
    1557         154 :     case t_REAL: {
    1558         154 :       GEN p = *pp;
    1559         154 :       if (signe(x) && (!p || abscmprr(p, x) < 0)) { *pp = x; *pi = k; }
    1560         154 :       return 0;
    1561             :     }
    1562           0 :     default: return !gequal0(x);
    1563             :   }
    1564             : }
    1565             : 
    1566             : /* Gauss reduction (arbitrary symetric matrix, only the part above the
    1567             :  * diagonal is considered). If signature is non-zero, return only the
    1568             :  * signature, in which case gsigne() should be defined for elements of a. */
    1569             : static GEN
    1570        4263 : gaussred(GEN a, long signature)
    1571             : {
    1572             :   GEN r, ak, al;
    1573             :   pari_sp av, av1;
    1574        4263 :   long n = lg(a), i, j, k, l, sp, sn, t;
    1575             : 
    1576        4263 :   if (typ(a) != t_MAT) pari_err_TYPE("gaussred",a);
    1577        4263 :   if (n == 1) return signature? mkvec2(gen_0, gen_0): cgetg(1, t_MAT);
    1578        4263 :   if (lgcols(a) != n) pari_err_DIM("gaussred");
    1579        4263 :   n--;
    1580             : 
    1581        4263 :   av = avma;
    1582        4263 :   r = const_vecsmall(n, 1);
    1583        4263 :   av1= avma;
    1584        4263 :   a = RgM_shallowcopy(a);
    1585        4263 :   t = n; sp = sn = 0;
    1586       27657 :   while (t)
    1587             :   {
    1588       23394 :     long pind = 0;
    1589       23394 :     GEN invp, p = NULL;
    1590       82544 :     k=1; while (k<=n && (!r[k] || !suitable(gcoeff(a,k,k), k, &p, &pind))) k++;
    1591       23394 :     if (k > n && p) k = pind;
    1592       23394 :     if (k <= n)
    1593             :     {
    1594       23387 :       p = gcoeff(a,k,k); invp = ginv(p); /* != 0 */
    1595       23387 :       if (signature) { /* skip if (!signature): gsigne may fail ! */
    1596       23338 :         if (gsigne(p) > 0) sp++; else sn++;
    1597             :       }
    1598       23387 :       r[k] = 0; t--;
    1599       23387 :       ak = row(a, k);
    1600      164682 :       for (i=1; i<=n; i++)
    1601      141295 :         gcoeff(a,k,i) = r[i]? gmul(gcoeff(a,k,i), invp): gen_0;
    1602             : 
    1603      164682 :       for (i=1; i<=n; i++) if (r[i])
    1604             :       {
    1605       58940 :         GEN c = gel(ak,i); /* - p * a[k,i] */
    1606       58940 :         if (gequal0(c)) continue;
    1607      445130 :         for (j=1; j<=n; j++) if (r[j])
    1608      237965 :           gcoeff(a,i,j) = gsub(gcoeff(a,i,j), gmul(c,gcoeff(a,k,j)));
    1609             :       }
    1610       23387 :       gcoeff(a,k,k) = p;
    1611       23387 :       if (gc_needed(av1,1))
    1612             :       {
    1613           0 :         if(DEBUGMEM>1) pari_warn(warnmem,"gaussred (t = %ld)", t);
    1614           0 :         a = gerepilecopy(av1, a);
    1615             :       }
    1616             :     }
    1617             :     else
    1618             :     { /* all remaining diagonal coeffs are currently 0 */
    1619           7 :       for (k=1; k<=n; k++) if (r[k])
    1620             :       {
    1621           7 :         l=k+1; while (l<=n && (!r[l] || !suitable(gcoeff(a,k,l), l, &p, &pind))) l++;
    1622           7 :         if (l > n && p) l = pind;
    1623           7 :         if (l > n) continue;
    1624             : 
    1625           7 :         p = gcoeff(a,k,l); invp = ginv(p);
    1626           7 :         sp++; sn++;
    1627           7 :         r[k] = r[l] = 0; t -= 2;
    1628           7 :         ak = row(a, k);
    1629           7 :         al = row(a, l);
    1630          35 :         for (i=1; i<=n; i++) if (r[i])
    1631             :         {
    1632          14 :           gcoeff(a,k,i) = gmul(gcoeff(a,k,i), invp);
    1633          14 :           gcoeff(a,l,i) = gmul(gcoeff(a,l,i), invp);
    1634             :         } else {
    1635          14 :           gcoeff(a,k,i) = gen_0;
    1636          14 :           gcoeff(a,l,i) = gen_0;
    1637             :         }
    1638             : 
    1639          35 :         for (i=1; i<=n; i++) if (r[i])
    1640             :         { /* c = a[k,i] * p, d = a[l,i] * p; */
    1641          14 :           GEN c = gel(ak,i), d = gel(al,i);
    1642          70 :           for (j=1; j<=n; j++) if (r[j])
    1643          28 :             gcoeff(a,i,j) = gsub(gcoeff(a,i,j),
    1644          28 :                                  gadd(gmul(gcoeff(a,l,j), c),
    1645          28 :                                       gmul(gcoeff(a,k,j), d)));
    1646             :         }
    1647          35 :         for (i=1; i<=n; i++) if (r[i])
    1648             :         {
    1649          14 :           GEN c = gcoeff(a,k,i), d = gcoeff(a,l,i);
    1650          14 :           gcoeff(a,k,i) = gadd(c, d);
    1651          14 :           gcoeff(a,l,i) = gsub(c, d);
    1652             :         }
    1653           7 :         gcoeff(a,k,l) = gen_1;
    1654           7 :         gcoeff(a,l,k) = gen_m1;
    1655           7 :         gcoeff(a,k,k) = gmul2n(p,-1);
    1656           7 :         gcoeff(a,l,l) = gneg(gcoeff(a,k,k));
    1657           7 :         if (gc_needed(av1,1))
    1658             :         {
    1659           0 :           if(DEBUGMEM>1) pari_warn(warnmem,"gaussred");
    1660           0 :           a = gerepilecopy(av1, a);
    1661             :         }
    1662           7 :         break;
    1663             :       }
    1664           7 :       if (k > n) break;
    1665             :     }
    1666             :   }
    1667        4263 :   if (!signature) return gerepilecopy(av, a);
    1668        4249 :   set_avma(av); return mkvec2s(sp, sn);
    1669             : }
    1670             : 
    1671             : GEN
    1672          14 : qfgaussred(GEN a) { return gaussred(a,0); }
    1673             : 
    1674             : GEN
    1675        4249 : qfsign(GEN a) { return gaussred(a,1); }
    1676             : 
    1677             : /* x -= s(y+u*x) */
    1678             : /* y += s(x-u*y), simultaneously */
    1679             : static void
    1680       19180 : rot(GEN x, GEN y, GEN s, GEN u) {
    1681       19180 :   GEN x1 = subrr(x, mulrr(s,addrr(y,mulrr(u,x))));
    1682       19180 :   GEN y1 = addrr(y, mulrr(s,subrr(x,mulrr(u,y))));
    1683       19180 :   affrr(x1,x);
    1684       19180 :   affrr(y1,y);
    1685       19180 : }
    1686             : 
    1687             : /* Diagonalization of a REAL symetric matrix. Return a vector [L, r]:
    1688             :  * L = vector of eigenvalues
    1689             :  * r = matrix of eigenvectors */
    1690             : GEN
    1691          28 : jacobi(GEN a, long prec)
    1692             : {
    1693             :   pari_sp av1;
    1694          28 :   long de, e, e1, e2, i, j, p, q, l = lg(a);
    1695             :   GEN c, ja, L, r, L2, r2, unr;
    1696             : 
    1697          28 :   if (typ(a) != t_MAT) pari_err_TYPE("jacobi",a);
    1698          28 :   ja = cgetg(3,t_VEC);
    1699          28 :   L = cgetg(l,t_COL); gel(ja,1) = L;
    1700          28 :   r = cgetg(l,t_MAT); gel(ja,2) = r;
    1701          28 :   if (l == 1) return ja;
    1702          28 :   if (lgcols(a) != l) pari_err_DIM("jacobi");
    1703             : 
    1704          28 :   e1 = HIGHEXPOBIT-1;
    1705         224 :   for (j=1; j<l; j++)
    1706             :   {
    1707         196 :     GEN z = gtofp(gcoeff(a,j,j), prec);
    1708         196 :     gel(L,j) = z;
    1709         196 :     e = expo(z); if (e < e1) e1 = e;
    1710             :   }
    1711         224 :   for (j=1; j<l; j++)
    1712             :   {
    1713         196 :     gel(r,j) = cgetg(l,t_COL);
    1714        1582 :     for (i=1; i<l; i++) gcoeff(r,i,j) = utor(i==j? 1: 0, prec);
    1715             :   }
    1716          28 :   av1 = avma;
    1717             : 
    1718          28 :   e2 = -(long)HIGHEXPOBIT; p = q = 1;
    1719          28 :   c = cgetg(l,t_MAT);
    1720         224 :   for (j=1; j<l; j++)
    1721             :   {
    1722         196 :     gel(c,j) = cgetg(j,t_COL);
    1723         791 :     for (i=1; i<j; i++)
    1724             :     {
    1725         595 :       GEN z = gtofp(gcoeff(a,i,j), prec);
    1726         595 :       gcoeff(c,i,j) = z;
    1727         595 :       if (!signe(z)) continue;
    1728         308 :       e = expo(z); if (e > e2) { e2 = e; p = i; q = j; }
    1729             :     }
    1730             :   }
    1731          28 :   a = c; unr = real_1(prec);
    1732          28 :   de = prec2nbits(prec);
    1733             : 
    1734             :  /* e1 = min expo(a[i,i])
    1735             :   * e2 = max expo(a[i,j]), i != j */
    1736        1568 :   while (e1-e2 < de)
    1737             :   {
    1738        1540 :     pari_sp av2 = avma;
    1739             :     GEN x, y, t, c, s, u;
    1740             :     /* compute attached rotation in the plane formed by basis vectors number
    1741             :      * p and q */
    1742        1540 :     x = subrr(gel(L,q),gel(L,p));
    1743        1540 :     if (signe(x))
    1744             :     {
    1745        1512 :       x = divrr(x, shiftr(gcoeff(a,p,q),1));
    1746        1512 :       y = sqrtr(addrr(unr, sqrr(x)));
    1747        1512 :       t = invr((signe(x)>0)? addrr(x,y): subrr(x,y));
    1748             :     }
    1749             :     else
    1750          28 :       y = t = unr;
    1751        1540 :     c = sqrtr(addrr(unr,sqrr(t)));
    1752        1540 :     s = divrr(t,c);
    1753        1540 :     u = divrr(t,addrr(unr,c));
    1754             : 
    1755             :     /* compute successive transforms of a and the matrix of accumulated
    1756             :      * rotations (r) */
    1757        4144 :     for (i=1;   i<p; i++) rot(gcoeff(a,i,p), gcoeff(a,i,q), s,u);
    1758        4039 :     for (i=p+1; i<q; i++) rot(gcoeff(a,p,i), gcoeff(a,i,q), s,u);
    1759        4487 :     for (i=q+1; i<l; i++) rot(gcoeff(a,p,i), gcoeff(a,q,i), s,u);
    1760        1540 :     y = gcoeff(a,p,q);
    1761        1540 :     t = mulrr(t, y); shiftr_inplace(y, -de - 1);
    1762        1540 :     x = gel(L,p); subrrz(x,t, x);
    1763        1540 :     y = gel(L,q); addrrz(y,t, y);
    1764       12670 :     for (i=1; i<l; i++) rot(gcoeff(r,i,p), gcoeff(r,i,q), s,u);
    1765             : 
    1766        1540 :     e2 = -(long)HIGHEXPOBIT; p = q = 1;
    1767       12670 :     for (j=1; j<l; j++)
    1768             :     {
    1769       46396 :       for (i=1; i<j; i++)
    1770             :       {
    1771       35266 :         GEN z = gcoeff(a,i,j);
    1772       35266 :         if (!signe(z)) continue;
    1773       31080 :         e = expo(z); if (e > e2) { e2=e; p=i; q=j; }
    1774             :       }
    1775       46396 :       for (i=j+1; i<l; i++)
    1776             :       {
    1777       35266 :         GEN z = gcoeff(a,j,i);
    1778       35266 :         if (!signe(z)) continue;
    1779       31080 :         e = expo(z); if (e > e2) { e2=e; p=j; q=i; }
    1780             :       }
    1781             :     }
    1782        1540 :     set_avma(av2);
    1783             :   }
    1784             :   /* sort eigenvalues from smallest to largest */
    1785          28 :   c = indexsort(L);
    1786         224 :   r2 = vecpermute(r, c); for (i=1; i<l; i++) gel(r,i) = gel(r2,i);
    1787         224 :   L2 = vecpermute(L, c); for (i=1; i<l; i++) gel(L,i) = gel(L2,i);
    1788          28 :   set_avma(av1); return ja;
    1789             : }
    1790             : 
    1791             : /*************************************************************************/
    1792             : /**                                                                     **/
    1793             : /**                   Q-vector space -> Z-modules                       **/
    1794             : /**                                                                     **/
    1795             : /*************************************************************************/
    1796             : 
    1797             : GEN
    1798         126 : matrixqz0(GEN x,GEN p)
    1799             : {
    1800         126 :   if (typ(x) != t_MAT) pari_err_TYPE("matrixqz",x);
    1801         126 :   if (!p) return QM_minors_coprime(x,NULL);
    1802          98 :   if (typ(p) != t_INT) pari_err_TYPE("matrixqz",p);
    1803          98 :   if (signe(p)>=0) return QM_minors_coprime(x,p);
    1804          91 :   if (!RgM_is_QM(x)) pari_err_TYPE("matrixqz", x);
    1805          91 :   if (absequaliu(p,1)) return QM_ImZ(x); /* p = -1 */
    1806          63 :   if (absequaliu(p,2)) return QM_ImQ(x); /* p = -2 */
    1807           7 :   pari_err_FLAG("QM_minors_coprime");
    1808             :   return NULL; /* LCOV_EXCL_LINE */
    1809             : }
    1810             : 
    1811             : GEN
    1812          35 : QM_minors_coprime(GEN x, GEN D)
    1813             : {
    1814          35 :   pari_sp av = avma, av1;
    1815             :   long i, j, m, n, lP;
    1816             :   GEN P, y;
    1817             : 
    1818          35 :   n = lg(x)-1; if (!n) return gcopy(x);
    1819          35 :   m = nbrows(x);
    1820          35 :   if (n > m) pari_err_DOMAIN("QM_minors_coprime","n",">",strtoGENstr("m"),x);
    1821          28 :   y = x; x = cgetg(n+1,t_MAT);
    1822          77 :   for (j=1; j<=n; j++)
    1823             :   {
    1824          49 :     gel(x,j) = Q_primpart(gel(y,j));
    1825          49 :     RgV_check_ZV(gel(x,j), "QM_minors_coprime");
    1826             :   }
    1827             :   /* x now a ZM */
    1828          28 :   if (n==m)
    1829             :   {
    1830          21 :     if (gequal0(ZM_det(x)))
    1831          14 :       pari_err_DOMAIN("QM_minors_coprime", "rank(A)", "<",stoi(n),x);
    1832           7 :     set_avma(av); return matid(n);
    1833             :   }
    1834             :   /* m > n */
    1835           7 :   if (!D || gequal0(D))
    1836             :   {
    1837           7 :     pari_sp av2 = avma;
    1838           7 :     D = ZM_detmult(shallowtrans(x));
    1839           7 :     if (is_pm1(D)) { set_avma(av2); return ZM_copy(x); }
    1840             :   }
    1841           7 :   P = gel(Z_factor(D), 1); lP = lg(P);
    1842           7 :   av1 = avma;
    1843          14 :   for (i=1; i < lP; i++)
    1844             :   {
    1845           7 :     GEN p = gel(P,i), pov2 = shifti(p, -1);
    1846             :     for(;;)
    1847          14 :     {
    1848          21 :       GEN N, M = FpM_ker(x, p);
    1849          21 :       long lM = lg(M);
    1850          21 :       if (lM==1) break;
    1851             : 
    1852          14 :       FpM_center_inplace(M, p, pov2);
    1853          14 :       N = ZM_Z_divexact(ZM_mul(x,M), p);
    1854          28 :       for (j=1; j<lM; j++)
    1855             :       {
    1856          14 :         long k = n; while (!signe(gcoeff(M,k,j))) k--;
    1857          14 :         gel(x,k) = gel(N,j);
    1858             :       }
    1859          14 :       if (gc_needed(av1,1))
    1860             :       {
    1861           0 :         if (DEBUGMEM>1) pari_warn(warnmem,"QM_minors_coprime, p = %Ps", p);
    1862           0 :         x = gerepilecopy(av1, x); pov2 = shifti(p, -1);
    1863             :       }
    1864             :     }
    1865             :   }
    1866           7 :   return gerepilecopy(av, x);
    1867             : }
    1868             : 
    1869             : static GEN
    1870        1162 : QM_ImZ_all_i(GEN A, GEN *U, long remove, long hnf, long linindep)
    1871             : {
    1872        1162 :   GEN V = NULL, D;
    1873        1162 :   A = Q_remove_denom(A,&D);
    1874        1162 :   if (D)
    1875             :   {
    1876             :     long l, lA;
    1877         126 :     V = matkermod(A,D,NULL);
    1878         126 :     l = lg(V); lA = lg(A);
    1879         126 :     if (l == 1) V = scalarmat_shallow(D, lA-1);
    1880             :     else
    1881             :     {
    1882          84 :       if (l < lA) V = hnfmodid(V,D);
    1883          84 :       A = ZM_Z_divexact(ZM_mul(A, V), D);
    1884             :     }
    1885             :   }
    1886        1162 :   if (!linindep && ZM_rank(A)==lg(A)-1) linindep = 1;
    1887        1162 :   if (hnf || !linindep) A = ZM_hnflll(A, U, remove);
    1888        1162 :   if (U && V)
    1889             :   {
    1890          35 :     if (hnf)    *U = ZM_mul(V,*U);
    1891           0 :     else        *U = V;
    1892             :   }
    1893        1162 :   return A;
    1894             : }
    1895             : GEN
    1896          28 : QM_ImZ_all(GEN x, GEN *U, long remove, long hnf)
    1897             : {
    1898          28 :   pari_sp av = avma;
    1899          28 :   x = QM_ImZ_all_i(x, U, remove, hnf, 0);
    1900          28 :   if (U) gerepileall(av, 2, &x, &U); else x = gerepilecopy(av, x);
    1901          28 :   return x;
    1902             : }
    1903             : GEN
    1904           0 : QM_ImZ_hnfall(GEN x, GEN *U, long remove) { return QM_ImZ_all(x, U, remove, 1); }
    1905             : GEN
    1906           0 : QM_ImZ_hnf(GEN x) { return QM_ImZ_hnfall(x, NULL, 1); }
    1907             : GEN
    1908          28 : QM_ImZ(GEN x) { return QM_ImZ_all(x, NULL, 1, 0); }
    1909             : 
    1910             : GEN
    1911        1141 : QM_ImQ_all(GEN x, GEN *U, long remove, long hnf)
    1912             : {
    1913        1141 :   pari_sp av = avma;
    1914        1141 :   long i, n = lg(x), m;
    1915        1141 :   GEN ir, V, D, c, K = NULL;
    1916             : 
    1917        1141 :   if (U) *U = matid(n-1);
    1918        1141 :   if (n==1) return gcopy(x);
    1919        1134 :   m = lg(gel(x,1));
    1920             : 
    1921        1134 :   x = RgM_shallowcopy(x);
    1922        7287 :   for (i=1; i<n; i++)
    1923             :   {
    1924        6153 :     gel(x,i) = Q_primitive_part(gel(x,i), &c);
    1925        6153 :     if (U && c && signe(c)) gcoeff(*U,i,i) = ginv(c);
    1926             :   }
    1927             : 
    1928        1134 :   ir = ZM_indexrank(x);
    1929        1134 :   if (U)
    1930             :   {
    1931          77 :     *U = vecpermute(*U, gel(ir,2));
    1932          77 :     if (remove < 2) K = ZM_ker(x);
    1933             :   }
    1934        1134 :   x = vecpermute(x, gel(ir,2));
    1935             : 
    1936        1134 :   D = absi(ZM_det(rowpermute(x,gel(ir,1))));
    1937        1134 :   x = RgM_Rg_div(x, D);
    1938        1134 :   x = QM_ImZ_all_i(x, U? &V: NULL, remove, hnf, 1);
    1939             : 
    1940        1134 :   if (U)
    1941             :   {
    1942          77 :     *U = RgM_Rg_div(RgM_mul(*U,V),D);
    1943          77 :     if (remove < 2) *U = shallowconcat(K,*U);
    1944          77 :     if (!remove) x = shallowconcat(zeromatcopy(m-1,lg(K)-1), x);
    1945          77 :     gerepileall(av, 2, &x, U);
    1946             :   }
    1947        1057 :   else x = gerepilecopy(av,x);
    1948        1134 :   return x;
    1949             : }
    1950             : GEN
    1951        1085 : QM_ImQ_hnfall(GEN x, GEN *U, long remove) { return QM_ImQ_all(x, U, remove, 1); }
    1952             : GEN
    1953        1008 : QM_ImQ_hnf(GEN x) { return QM_ImQ_hnfall(x, NULL, 1); }
    1954             : GEN
    1955          56 : QM_ImQ(GEN x) { return QM_ImQ_all(x, NULL, 1, 0); }
    1956             : 
    1957             : GEN
    1958        3612 : intersect(GEN x, GEN y)
    1959             : {
    1960        3612 :   long j, lx = lg(x);
    1961             :   pari_sp av;
    1962             :   GEN z;
    1963             : 
    1964        3612 :   if (typ(x)!=t_MAT) pari_err_TYPE("intersect",x);
    1965        3612 :   if (typ(y)!=t_MAT) pari_err_TYPE("intersect",y);
    1966        3612 :   if (lx==1 || lg(y)==1) return cgetg(1,t_MAT);
    1967             : 
    1968        3612 :   av = avma; z = ker(shallowconcat(x,y));
    1969       15946 :   for (j=lg(z)-1; j; j--) setlg(z[j], lx);
    1970        3612 :   return gerepileupto(av, image(RgM_mul(x,z)));
    1971             : }

Generated by: LCOV version 1.13