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-bordeaux1.fr machine (x86_64 architecture), and agregate them in the final report:

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

LCOV - code coverage report
Current view: top level - basemath - qfsolve.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.8.0 lcov report (development 16624-25b9976) Lines: 581 590 98.5 %
Date: 2014-06-24 Functions: 29 29 100.0 %
Legend: Lines: hit not hit | Branches: + taken - not taken # not executed Branches: 348 383 90.9 %

           Branch data     Line data    Source code
       1                 :            : /* Copyright (C) 2000-2004  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                 :            : /* Copyright (C) 2014 Denis Simon
      15                 :            :  * Adapted from qfsolve.gp v. 09/01/2014
      16                 :            :  *   http://www.math.unicaen.fr/~simon/qfsolve.gp
      17                 :            :  *
      18                 :            :  * Author: Denis SIMON <simon@math.unicaen.fr> */
      19                 :            : 
      20                 :            : #include "pari.h"
      21                 :            : #include "paripriv.h"
      22                 :            : 
      23                 :            : /* LINEAR ALGEBRA */
      24                 :            : /* complete by 0s, assume l-1 <= n */
      25                 :            : static GEN
      26                 :      24030 : vecextend(GEN v, long n)
      27                 :            : {
      28                 :      24030 :   long i, l = lg(v);
      29                 :      24030 :   GEN w = cgetg(n+1, t_COL);
      30         [ +  + ]:      79228 :   for (i = 1; i < l; i++) gel(w,i) = gel(v,i);
      31         [ +  + ]:     157135 :   for (     ; i <=n; i++) gel(w,i) = gen_0;
      32                 :      24030 :   return w;
      33                 :            : }
      34                 :            : 
      35                 :            : /* Gives a unimodular matrix with the last column(s) equal to Mv.
      36                 :            :  * Mv can be a column vector or a rectangular matrix.
      37                 :            :  * redflag = 0 or 1. If redflag = 1, LLL-reduce the n-#v first columns. */
      38                 :            : static GEN
      39                 :      59541 : completebasis(GEN Mv, long redflag)
      40                 :            : {
      41                 :            :   GEN U;
      42                 :            :   long m, n;
      43                 :            : 
      44         [ +  + ]:      59541 :   if (typ(Mv) == t_COL) Mv = mkmat(Mv);
      45                 :      59541 :   n = lg(Mv)-1;
      46                 :      59541 :   m = nbrows(Mv); /* m x n */
      47         [ +  + ]:      59541 :   if (m == n) return Mv;
      48                 :      59491 :   (void)ZM_hnfall(shallowtrans(Mv), &U, 0);
      49                 :      59491 :   U = ZM_inv(shallowtrans(U), gen_1);
      50 [ +  - ][ +  + ]:      59491 :   if (m==1 || !redflag) return U;
      51                 :            :   /* LLL-reduce the m-n first columns */
      52                 :      59541 :   return shallowconcat(ZM_lll(vecslice(U,1,m-n), 0.99, LLL_INPLACE),
      53                 :      23872 :                        vecslice(U, m-n+1,m));
      54                 :            : }
      55                 :            : 
      56                 :            : /* Compute the kernel of M mod p.
      57                 :            :  * returns [d,U], where
      58                 :            :  * d = dim (ker M mod p)
      59                 :            :  * U in GLn(Z), and its first d columns span the kernel. */
      60                 :            : static GEN
      61                 :      35522 : kermodp(GEN M, GEN p, long *d)
      62                 :            : {
      63                 :            :   long j, l;
      64                 :            :   GEN K, B, U;
      65                 :            : 
      66                 :      35522 :   K = FpM_center(FpM_ker(M, p), p, shifti(p,-1));
      67                 :      35522 :   B = completebasis(K,0);
      68                 :      35522 :   l = lg(M); U = cgetg(l, t_MAT);
      69         [ +  + ]:     299960 :   for (j =  1; j < l; j++) gel(U,j) = gel(B,l-j);
      70                 :      35522 :   *d = lg(K)-1; return U;
      71                 :            : }
      72                 :            : 
      73                 :            : /* INVARIANTS COMPUTATIONS */
      74                 :            : 
      75                 :            : static GEN
      76                 :      20597 : principal_minor(GEN G, long  i)
      77                 :      20597 : { return rowslice(vecslice(G,1,i), 1,i); }
      78                 :            : static GEN
      79                 :        560 : det_minors(GEN G)
      80                 :            : {
      81                 :        560 :   long i, l = lg(G);
      82                 :        560 :   GEN v = cgetg(l+1, t_VEC);
      83                 :        560 :   gel(v,1) = gen_1;
      84         [ +  + ]:       2800 :   for (i = 2; i <= l; i++) gel(v,i) = ZM_det(principal_minor(G,i-1));
      85                 :        560 :   return v;
      86                 :            : }
      87                 :            : 
      88                 :            : /* Given a symmetric matrix G over Z, compute the Witt invariant
      89                 :            :  *  of G at the prime p (at real place if p = NULL)
      90                 :            :  * Assume that none of the determinant G[1..i,1..i] is 0. */
      91                 :            : static long
      92                 :         85 : qflocalinvariant(GEN G, GEN p)
      93                 :            : {
      94                 :         85 :   long i, j, c, l = lg(G);
      95                 :         85 :   GEN diag, v = det_minors(G);
      96                 :            :   /* Diagonalize G first. */
      97                 :         85 :   diag = cgetg(l, t_VEC);
      98         [ +  + ]:        425 :   for (i = 1; i < l; i++) gel(diag,i) = mulii(gel(v,i+1), gel(v,i));
      99                 :            : 
     100                 :            :   /* Then compute the product of the Hilbert symbols */
     101                 :            :   /* (diag[i],diag[j])_p for i < j */
     102                 :         85 :   c = 1;
     103         [ +  + ]:        340 :   for (i = 1; i < l-1; i++)
     104         [ +  + ]:        765 :     for (j = i+1; j < l; j++)
     105         [ +  + ]:        510 :       if (hilbertii(gel(diag,i), gel(diag,j), p) < 0) c = -c;
     106                 :         85 :   return c;
     107                 :            : }
     108                 :            : 
     109                 :            : static GEN
     110                 :       5179 : hilberts(GEN a, GEN b, GEN P, long lP)
     111                 :            : {
     112                 :       5179 :   GEN v = cgetg(lP, t_VECSMALL);
     113                 :            :   long i;
     114         [ +  + ]:      30055 :   for (i = 1; i < lP; i++) v[i] = hilbertii(a, b, gel(P,i)) < 0;
     115                 :       5179 :   return v;
     116                 :            : }
     117                 :            : 
     118                 :            : /* G symmetrix matrix or qfb or list of quadratic forms with same discriminant.
     119                 :            :  * P must be equal to factor(-abs(2*matdet(G)))[,1]. */
     120                 :            : static GEN
     121                 :       3286 : qflocalinvariants(GEN G, GEN P)
     122                 :            : {
     123                 :            :   GEN sol;
     124                 :       3286 :   long i, j, l, lP = lg(P);
     125                 :            : 
     126                 :            :   /* convert G into a vector of symmetric matrices */
     127         [ +  + ]:       3286 :   G = (typ(G) == t_VEC)? shallowcopy(G): mkvec(G);
     128                 :       3286 :   l = lg(G);
     129         [ +  + ]:       8555 :   for (j = 1; j < l; j++)
     130                 :            :   {
     131                 :       5269 :     GEN g = gel(G,j);
     132 [ +  + ][ +  + ]:       5269 :     if (typ(g) == t_QFI || typ(g) == t_QFR) gel(G,j) = gtomat(g);
     133                 :            :   }
     134                 :       3286 :   sol = cgetg(l, t_MAT);
     135         [ +  + ]:       3286 :   if (lg(gel(G,1)) == 3)
     136                 :            :   { /* in dimension 2, each invariant is a single Hilbert symbol. */
     137                 :       2811 :     GEN d = negi(ZM_det(gel(G,1)));
     138         [ +  + ]:       7605 :     for (j = 1; j < l; j++)
     139                 :            :     {
     140                 :       4794 :       GEN a = gcoeff(gel(G,j),1,1);
     141                 :       4794 :       gel(sol,j) = hilberts(a, d, P, lP);
     142                 :            :     }
     143                 :            :   }
     144                 :            :   else /* in dimension n > 2, we compute a product of n Hilbert symbols. */
     145         [ +  + ]:        950 :     for (j = 1; j <l; j++)
     146                 :            :     {
     147                 :        475 :       GEN g = gel(G,j), v = det_minors(g), w = cgetg(lP, t_VECSMALL);
     148                 :        475 :       long n = lg(v);
     149                 :        475 :       gel(sol,j) = w;
     150         [ +  + ]:       2855 :       for (i = 1; i < lP; i++)
     151                 :            :       {
     152                 :       2380 :         GEN p = gel(P,i);
     153                 :       2380 :         long k = n-2, h = hilbertii(gel(v,k), gel(v,k+1),p);
     154         [ +  + ]:       9520 :         for (k--; k >= 1; k--)
     155         [ +  + ]:       7140 :           if (hilbertii(negi(gel(v,k)), gel(v,k+1),p) < 0) h = -h;
     156                 :       2380 :         w[i] = h < 0;
     157                 :            :       }
     158                 :            :     }
     159                 :       3286 :   return sol;
     160                 :            : }
     161                 :            : 
     162                 :            : /* QUADRATIC FORM REDUCTION */
     163                 :            : static GEN
     164                 :       6031 : qfb(GEN D, GEN a, GEN b, GEN c)
     165                 :            : {
     166         [ +  + ]:       6031 :   if (signe(D) < 0) return mkqfi(a,b,c);
     167                 :       6031 :   retmkqfr(a,b,c,real_0(DEFAULTPREC));
     168                 :            : }
     169                 :            : 
     170                 :            : /* Gauss reduction of the binary quadratic form
     171                 :            :  * Q = a*X^2+2*b*X*Y+c*Y^2 of discriminant D (divisible by 4)
     172                 :            :  * returns the reduced form */
     173                 :            : static GEN
     174                 :       1237 : qfbreduce(GEN D, GEN Q)
     175                 :            : {
     176                 :       1237 :   GEN a = gel(Q,1), b = shifti(gel(Q,2),-1), c = gel(Q,3);
     177         [ +  - ]:       3789 :   while (signe(a))
     178                 :            :   {
     179                 :            :     GEN r, q, nexta, nextc;
     180                 :       3789 :     q = dvmdii(b,a, &r); /* FIXME: export as dvmdiiround ? */
     181 [ +  + ][ +  + ]:       3789 :     if (signe(r) > 0 && absi_cmp(shifti(r,1), a) > 0) {
     182                 :        777 :       r = subii(r, absi(a)); q = addis(q, signe(a));
     183                 :            :     }
     184                 :       3789 :     nextc = a; nexta = subii(c, mulii(q, addii(r,b)));
     185         [ +  + ]:       3789 :     if (absi_cmp(nexta, a) >= 0) break;
     186                 :       3789 :     c = nextc; b = negi(r); a = nexta;
     187                 :            :   }
     188                 :       1237 :   return qfb(D,a,shifti(b,1),c);
     189                 :            : }
     190                 :            : 
     191                 :            : /* private version of qfgaussred:
     192                 :            :  * - early abort if k-th principal minor is singular, return stoi(k)
     193                 :            :  * - else return a matrix whose upper triangular part is qfgaussred(a) */
     194                 :            : static GEN
     195                 :      15090 : partialgaussred(GEN a)
     196                 :            : {
     197                 :      15090 :   long n = lg(a)-1, k;
     198                 :      15090 :   a = RgM_shallowcopy(a);
     199         [ +  + ]:      64204 :   for(k = 1; k < n; k++)
     200                 :            :   {
     201                 :      53500 :     GEN ak, p = gcoeff(a,k,k);
     202                 :            :     long i, j;
     203         [ +  + ]:      53500 :     if (isintzero(p)) return stoi(k);
     204                 :      49114 :     ak = row(a, k);
     205         [ +  + ]:     212796 :     for (i=k+1; i<=n; i++) gcoeff(a,k,i) = gdiv(gcoeff(a,k,i), p);
     206         [ +  + ]:     212796 :     for (i=k+1; i<=n; i++)
     207                 :            :     {
     208                 :     163682 :       GEN c = gel(ak,i);
     209         [ +  + ]:     163682 :       if (gequal0(c)) continue;
     210         [ +  + ]:     552627 :       for (j=i; j<=n; j++)
     211                 :     410043 :         gcoeff(a,i,j) = gsub(gcoeff(a,i,j), gmul(c,gcoeff(a,k,j)));
     212                 :            :     }
     213                 :            :   }
     214                 :      15090 :   return a;
     215                 :            : }
     216                 :            : 
     217                 :            : /* LLL-reduce a positive definite qf QD bounding the indefinite G, dim G > 1.
     218                 :            :  * Then finishes by looking for trivial solution */
     219                 :            : static GEN qftriv(GEN G, GEN z, long base);
     220                 :            : static GEN
     221                 :      15090 : qflllgram_indef(GEN G, long base)
     222                 :            : {
     223                 :            :   GEN M, R, g, DM, S, dR;
     224                 :      15090 :   long i, j, n = lg(G)-1;
     225                 :            : 
     226                 :      15090 :   R = partialgaussred(G);
     227         [ +  + ]:      15090 :   if (typ(R) == t_INT) return qftriv(G, R, base);
     228                 :      10704 :   R = Q_remove_denom(R, &dR); /* avoid rational arithmetic */
     229                 :      10704 :   M = zeromatcopy(n,n);
     230                 :      10704 :   DM = zeromatcopy(n,n);
     231         [ +  + ]:      67359 :   for (i = 1; i <= n; i++)
     232                 :            :   {
     233                 :      56655 :     GEN d = absi_shallow(gcoeff(R,i,i));
     234         [ +  + ]:      56655 :     if (dR) {
     235                 :      50618 :       gcoeff(M,i,i) = dR;
     236                 :      50618 :       gcoeff(DM,i,i) = mulii(d,dR);
     237                 :            :     } else {
     238                 :       6037 :       gcoeff(M,i,i) = gen_1;
     239                 :       6037 :       gcoeff(DM,i,i) = d;
     240                 :            :     }
     241         [ +  + ]:     206094 :     for (j = i+1; j <= n; j++)
     242                 :            :     {
     243                 :     149439 :       gcoeff(M,i,j) = gcoeff(R,i,j);
     244                 :     149439 :       gcoeff(DM,i,j) = mulii(d, gcoeff(R,i,j));
     245                 :            :     }
     246                 :            :   }
     247                 :            :   /* G = M~*D*M, D diagonal, DM=|D|*M, g =  M~*|D|*M */
     248                 :      10704 :   g = ZM_transmultosym(M,DM);
     249                 :      10704 :   S = lllgramint(Q_primpart(g));
     250                 :      10704 :   R = qftriv(qf_apply_ZM(G,S), NULL, base);
     251      [ +  +  + ]:      10704 :   switch(typ(R))
     252                 :            :   {
     253                 :        154 :     case t_COL: return ZM_ZC_mul(S,R);
     254                 :       5130 :     case t_MAT: return mkvec2(R, S);
     255                 :            :     default:
     256                 :       5420 :       gel(R,2) = ZM_mul(S, gel(R,2));
     257                 :      15090 :       return R;
     258                 :            :   }
     259                 :            : }
     260                 :            : 
     261                 :            : /* G symetric, i < j, let E = E_{i,j}(a), G <- E~*G*E,  U <- U*E.
     262                 :            :  * Everybody integral */
     263                 :            : static void
     264                 :       2957 : qf_apply_transvect_Z(GEN G, GEN U, long i, long j, GEN a)
     265                 :            : {
     266                 :       2957 :   long k, n = lg(G)-1;
     267                 :       2957 :   gel(G, j) =  ZC_lincomb(gen_1, a, gel(G,j), gel(G,i));
     268         [ +  + ]:      16401 :   for (k = 1; k < n; k++) gcoeff(G, j, k) = gcoeff(G, k, j);
     269                 :       2957 :   gcoeff(G,j,j) = addmulii(gcoeff(G,j,j), a,
     270                 :       5914 :                            addmulii(gcoeff(G,i,j), a,gcoeff(G,i,i)));
     271                 :       2957 :   gel(U, j) =  ZC_lincomb(gen_1, a, gel(U,j), gel(U,i));
     272                 :       2957 : }
     273                 :            : 
     274                 :            : /* LLL reduction of the quadratic form G (Gram matrix)
     275                 :            :  * where we go on, even if an isotropic vector is found. */
     276                 :            : static GEN
     277                 :       7824 : qflllgram_indefgoon(GEN G)
     278                 :            : {
     279                 :            :   GEN red, U, A, U1,U2,U3,U5,U6, V, B, G2,G3,G4,G5, G6, a, g;
     280                 :       7824 :   long i, j, n = lg(G)-1;
     281                 :            : 
     282                 :       7824 :   red = qflllgram_indef(G,1);
     283         [ -  + ]:       7824 :   if (typ(red)==t_MAT) return red; /*no isotropic vector found: nothing to do*/
     284                 :            :   /* otherwise a solution is found: */
     285                 :       7824 :   U1 = gel(red,2);
     286                 :       7824 :   G2 = gel(red,1); /* G2[1,1] = 0 */
     287                 :       7824 :   U2 = gel(mathnf0(row(G2,1), 4), 2);
     288                 :       7824 :   G3 = qf_apply_ZM(G2,U2);
     289                 :       7824 :   U = ZM_mul(U1,U2); /* qf_apply(G,U) = G3 */
     290                 :            :   /* G3[1,] = [0,...,0,g], g^2 | det G */
     291                 :       7824 :   g = gcoeff(G3,1,n);
     292                 :       7824 :   a = diviiround(negi(gcoeff(G3,n,n)), shifti(g,1));
     293         [ +  + ]:       7824 :   if (signe(a)) qf_apply_transvect_Z(G3,U,1,n,a);
     294                 :            :   /* G3[n,n] reduced mod 2g */
     295         [ +  + ]:       7824 :   if (n == 2) return mkvec2(G3,U);
     296                 :       7349 :   V = rowpermute(vecslice(G3, 2,n-1), mkvecsmall2(1,n));
     297                 :       7349 :   A = mkmat2(mkcol2(gcoeff(G3,1,1),gcoeff(G3,1,n)),
     298                 :      14698 :              mkcol2(gcoeff(G3,1,n),gcoeff(G3,2,2)));
     299                 :       7349 :   B = ground(RgM_neg(RgM_mul(RgM_inv(A), V)));
     300                 :       7349 :   U3 = matid(n);
     301         [ +  + ]:      36562 :   for (j = 2; j < n; j++)
     302                 :            :   {
     303                 :      29213 :     gcoeff(U3,1,j) = gcoeff(B,1,j-1);
     304                 :      29213 :     gcoeff(U3,n,j) = gcoeff(B,2,j-1);
     305                 :            :   }
     306                 :       7349 :   G4 = qf_apply_ZM(G3,U3); /* the last column of G4 is reduced */
     307                 :       7349 :   U = ZM_mul(U,U3);
     308         [ +  + ]:       7349 :   if (n == 3) return mkvec2(G4,U);
     309                 :            : 
     310                 :       5744 :   red = qflllgram_indefgoon(rowslice(vecslice(G4,2,n-1),2,n-1));
     311         [ -  + ]:       5744 :   if (typ(red) == t_MAT) return mkvec2(G4,U);
     312                 :            :   /* Let U5:=matconcat(diagonal[1,red[2],1])
     313                 :            :    * return [qf_apply_ZM(G5, U5), U*U5] */
     314                 :       5744 :   G5 = gel(red,1);
     315                 :       5744 :   U5 = gel(red,2);
     316                 :       5744 :   G6 = cgetg(n+1,t_MAT);
     317                 :       5744 :   gel(G6,1) = gel(G4,1);
     318                 :       5744 :   gel(G6,n) = gel(G4,n);
     319         [ +  + ]:      33352 :   for (j=2; j<n; j++)
     320                 :            :   {
     321                 :      27608 :     gel(G6,j) = cgetg(n+1,t_COL);
     322                 :      27608 :     gcoeff(G6,1,j) = gcoeff(G4,j,1);
     323                 :      27608 :     gcoeff(G6,n,j) = gcoeff(G4,j,n);
     324         [ +  + ]:     182574 :     for (i=2; i<n; i++) gcoeff(G6,i,j) = gcoeff(G5,i-1,j-1);
     325                 :            :   }
     326                 :       5744 :   U6 = mkvec3(mkmat(gel(U,1)), ZM_mul(vecslice(U,2,n-1),U5), mkmat(gel(U,n)));
     327                 :       7824 :   return mkvec2(G6, shallowconcat1(U6));
     328                 :            : }
     329                 :            : 
     330                 :            : /* qf_apply_ZM(G,H),  where H = matrix of \tau_{i,j}, i != j */
     331                 :            : static GEN
     332                 :       8794 : qf_apply_tau(GEN G, long i, long j)
     333                 :            : {
     334                 :       8794 :   long l = lg(G), k;
     335                 :       8794 :   G = RgM_shallowcopy(G);
     336                 :       8794 :   swap(gel(G,i), gel(G,j));
     337         [ +  + ]:      55275 :   for (k = 1; k < l; k++) swap(gcoeff(G,i,k), gcoeff(G,j,k));
     338                 :       8794 :   return G;
     339                 :            : }
     340                 :            : 
     341                 :            : /* LLL reduction of the quadratic form G (Gram matrix)
     342                 :            :  * in dim 3 only, with detG = -1 and sign(G) = [2,1]; */
     343                 :            : static GEN
     344                 :       1611 : qflllgram_indefgoon2(GEN G)
     345                 :            : {
     346                 :            :   GEN red, G2, a, b, c, d, e, f, u, v, r, r3, U2, G3;
     347                 :            : 
     348                 :       1611 :   red = qflllgram_indef(G,1); /* always find an isotropic vector. */
     349                 :       1611 :   G2 = qf_apply_tau(gel(red,1),1,3); /* G2[3,3] = 0 */
     350                 :       1611 :   r = row(gel(red,2), 3);
     351                 :       1611 :   swap(gel(r,1), gel(r,3)); /* apply tau_{1,3} */
     352                 :       1611 :   a = gcoeff(G2,3,1);
     353                 :       1611 :   b = gcoeff(G2,3,2);
     354                 :       1611 :   d = bezout(a,b, &u,&v);
     355         [ -  + ]:       1611 :   if (!is_pm1(d))
     356                 :            :   {
     357                 :          0 :     a = diviiexact(a,d);
     358                 :          0 :     b = diviiexact(b,d);
     359                 :            :   }
     360                 :            :   /* for U2 = [-u,-b,0;-v,a,0;0,0,1]
     361                 :            :    * G3 = qf_apply_ZM(G2,U2) has known last row (-d, 0, 0),
     362                 :            :    * so apply to principal_minor(G3,2), instead */
     363                 :       1611 :   U2 = mkmat2(mkcol2(negi(u),negi(v)), mkcol2(negi(b),a));
     364                 :       1611 :   G3 = qf_apply_ZM(principal_minor(G2,2),U2);
     365                 :       1611 :   r3 = gel(r,3);
     366                 :       1611 :   r = ZV_ZM_mul(mkvec2(gel(r,1),gel(r,2)),U2);
     367                 :            : 
     368                 :       1611 :   a = gcoeff(G3,1,1);
     369                 :       1611 :   b = gcoeff(G3,1,2);
     370                 :       1611 :   c = negi(d); /* G3[1,3] */
     371                 :       1611 :   d = gcoeff(G3,2,2);
     372         [ +  + ]:       1611 :   if (mpodd(a))
     373                 :            :   {
     374                 :        897 :     e = addii(b,d);
     375                 :        897 :     a = addii(a, addii(b,e));
     376                 :        897 :     e = diviiround(negi(e),c);
     377                 :        897 :     f = diviiround(negi(a), shifti(c,1));
     378                 :        897 :     a = addmulii(addii(gel(r,1),gel(r,2)), f,r3);
     379                 :            :   }
     380                 :            :   else
     381                 :            :   {
     382                 :        714 :     e = diviiround(negi(b),c);
     383                 :        714 :     f = diviiround(negi(shifti(a,-1)), c);
     384                 :        714 :     a = addmulii(gel(r,1), f, r3);
     385                 :            :   }
     386                 :       1611 :   b = addmulii(gel(r,2), e, r3);
     387                 :       1611 :   return mkvec3(a,b, r3);
     388                 :            : }
     389                 :            : 
     390                 :            : /* QUADRATIC FORM MINIMIZATION */
     391                 :            : /* G symetric, return ZM_Z_divexact(G,d) */
     392                 :            : static GEN
     393                 :      43536 : ZsymM_Z_divexact(GEN G, GEN d)
     394                 :            : {
     395                 :      43536 :   long i,j,l = lg(G);
     396                 :      43536 :   GEN H = cgetg(l, t_MAT);
     397         [ +  + ]:     314666 :   for(j=1; j<l; j++)
     398                 :            :   {
     399                 :     271130 :     GEN c = cgetg(l, t_COL), b = gel(G,j);
     400         [ +  + ]:    1141454 :     for(i=1; i<j; i++) gcoeff(H,j,i) = gel(c,i) = diviiexact(gel(b,i),d);
     401                 :     271130 :     gel(c,j) = diviiexact(gel(b,j),d);
     402                 :     271130 :     gel(H,j) = c;
     403                 :            :   }
     404                 :      43536 :   return H;
     405                 :            : }
     406                 :            : 
     407                 :            : /* write symetric G as [A,B;B~,C], A dxd, C (n-d)x(n-d) */
     408                 :            : static void
     409                 :        355 : blocks4(GEN G, long d, long n, GEN *A, GEN *B, GEN *C)
     410                 :            : {
     411                 :        355 :   GEN G2 = vecslice(G,d+1,n);
     412                 :        355 :   *A = principal_minor(G, d);
     413                 :        355 :   *B = rowslice(G2, 1, d);
     414                 :        355 :   *C = rowslice(G2, d+1, n);
     415                 :        355 : }
     416                 :            : /* Minimization of the quadratic form G, deg G != 0, dim n >= 2
     417                 :            :  * G symmetric integral
     418                 :            :  * Returns [G',U,factd] with U in GLn(Q) such that G'=U~*G*U*constant
     419                 :            :  * is integral and has minimal determinant.
     420                 :            :  * In dimension 3 or 4, may return a prime p if the reduction at p is
     421                 :            :  * impossible because of local non-solvability.
     422                 :            :  * P,E = factor(+/- det(G)), "prime" -1 is ignored. Destroy E. */
     423                 :            : static GEN qfsolvemodp(GEN G, GEN p);
     424                 :            : static GEN
     425                 :       5850 : qfminimize(GEN G, GEN P, GEN E)
     426                 :            : {
     427                 :            :   GEN d, U, Ker, sol, aux, faE, faP;
     428                 :       5850 :   long n = lg(G)-1, lP = lg(P), i, dimKer, m;
     429                 :            : 
     430                 :       5850 :   faP = vectrunc_init(lP);
     431                 :       5850 :   faE = vecsmalltrunc_init(lP);
     432                 :       5850 :   U = NULL;
     433         [ +  + ]:      58200 :   for (i = 1; i < lP; i++)
     434                 :            :   {
     435                 :      52690 :     GEN p = gel(P,i);
     436                 :      52690 :     long vp = E[i];
     437 [ +  + ][ -  + ]:      52690 :     if (!vp || !p) continue;
     438                 :            : 
     439         [ -  + ]:      40567 :     if (DEBUGLEVEL >= 4) err_printf("    p^v = %Ps^%ld\n", p,vp);
     440                 :            :     /* The case vp = 1 can be minimized only if n is odd. */
     441 [ +  + ][ +  + ]:      40567 :     if (vp == 1 && n%2 == 0) {
     442                 :       5400 :       vectrunc_append(faP, p);
     443                 :       5400 :       vecsmalltrunc_append(faE, 1);
     444                 :       5400 :       continue;
     445                 :            :     }
     446                 :      35167 :     Ker = kermodp(G,p, &dimKer); /* dimKer <= vp */
     447         [ -  + ]:      35167 :     if (DEBUGLEVEL >= 4) err_printf("    dimKer = %ld\n",dimKer);
     448         [ -  + ]:      35167 :     if (dimKer == n)
     449                 :            :     { /* trivial case: dimKer = n */
     450         [ #  # ]:          0 :       if (DEBUGLEVEL >= 4) err_printf("     case 0: dimKer = n\n");
     451                 :          0 :       G = ZsymM_Z_divexact(G, p);
     452                 :          0 :       E[i] -= n;
     453                 :          0 :       i--; continue; /* same p */
     454                 :            :     }
     455                 :      35167 :     G = qf_apply_ZM(G, Ker);
     456         [ +  + ]:      35167 :     U = U? RgM_mul(U,Ker): Ker;
     457                 :            : 
     458                 :            :     /* 1st case: dimKer < vp */
     459                 :            :     /* then the kernel mod p contains a kernel mod p^2 */
     460         [ +  + ]:      35167 :     if (dimKer < vp)
     461                 :            :     {
     462         [ -  + ]:       1660 :       if (DEBUGLEVEL >= 4) err_printf("    case 1: dimker < vp\n");
     463         [ +  + ]:       1660 :       if (dimKer == 1)
     464                 :            :       {
     465                 :            :         long j;
     466                 :       1305 :         gel(G,1) = ZC_Z_divexact(gel(G,1), p);
     467         [ +  + ]:       8420 :         for (j = 1; j<=n; j++) gcoeff(G,1,j) = diviiexact(gcoeff(G,1,j), p);
     468                 :       1305 :         gel(U,1) = RgC_Rg_div(gel(U,1), p);
     469                 :       1305 :         E[i] -= 2;
     470                 :            :       }
     471                 :            :       else
     472                 :            :       {
     473                 :        355 :         GEN A,B,C, K2 = ZsymM_Z_divexact(principal_minor(G,dimKer),p);
     474                 :            :         long j, dimKer2;
     475                 :        355 :         K2 = kermodp(K2, p, &dimKer2);
     476         [ +  + ]:        730 :         for (j = dimKer2+1; j <= dimKer; j++) gel(K2,j) = ZC_Z_mul(gel(K2,j),p);
     477                 :            :         /* Write G = [A,B;B~,C] and apply [K2,0;0,p*Id]/p by blocks */
     478                 :        355 :         blocks4(G, dimKer,n, &A,&B,&C);
     479                 :        355 :         A = ZsymM_Z_divexact(qf_apply_ZM(A,K2), sqri(p));
     480                 :        355 :         B = ZM_Z_divexact(ZM_transmul(B,K2), p);
     481                 :        355 :         G = shallowmatconcat(mkmat2(mkcol2(A,B),
     482                 :            :                                     mkcol2(shallowtrans(B), C)));
     483                 :            :         /* U *= [K2,0;0,Id] */
     484                 :        355 :         U = shallowconcat(RgM_Rg_div(RgM_mul(vecslice(U,1,dimKer),K2), p),
     485                 :            :                           vecslice(U,dimKer+1,n));
     486                 :        355 :         E[i] -= 2*dimKer2;
     487                 :            :       }
     488                 :       1660 :       i--; continue; /* same p */
     489                 :            :     }
     490                 :            : 
     491                 :            :    /* vp = dimKer
     492                 :            :     * 2nd case: kernel has dim >= 2 and contains an element of norm 0 mod p^2
     493                 :            :     * search for an element of norm p^2... in the kernel */
     494                 :      33507 :     sol = NULL;
     495         [ +  + ]:      33507 :     if (dimKer > 2) {
     496         [ -  + ]:      12819 :       if (DEBUGLEVEL >= 4) err_printf("    case 2.1\n");
     497                 :      12819 :       dimKer = 3;
     498                 :      12819 :       sol = qfsolvemodp(ZsymM_Z_divexact(principal_minor(G,3),p),  p);
     499                 :      12819 :       sol = FpC_red(sol, p);
     500                 :            :     }
     501         [ +  + ]:      20688 :     else if (dimKer == 2)
     502                 :            :     {
     503                 :      11073 :       GEN a = modii(diviiexact(gcoeff(G,1,1),p), p);
     504                 :      11073 :       GEN b = modii(diviiexact(gcoeff(G,1,2),p), p);
     505                 :      11073 :       GEN c = diviiexact(gcoeff(G,2,2),p);
     506                 :      11073 :       GEN di= modii(subii(sqri(b), mulii(a,c)), p);
     507         [ +  + ]:      11073 :       if (kronecker(di,p) >= 0)
     508                 :            :       {
     509         [ -  + ]:      11038 :         if (DEBUGLEVEL >= 4) err_printf("    case 2.2\n");
     510         [ +  + ]:      11038 :         sol = signe(a)? mkcol2(Fp_sub(Fp_sqrt(di,p), b, p), a): vec_ei(2,1);
     511                 :            :       }
     512                 :            :     }
     513         [ +  + ]:      33507 :     if (sol)
     514                 :            :     {
     515                 :            :       long j;
     516                 :      23857 :       sol = FpC_center(sol, p, shifti(p,-1));
     517                 :      23857 :       sol = Q_primpart(sol);
     518         [ -  + ]:      23857 :       if (DEBUGLEVEL >= 4) err_printf("    sol = %Ps\n", sol);
     519                 :      23857 :       Ker = completebasis(vecextend(sol,n), 1);
     520         [ +  + ]:     187130 :       for(j=1; j<n; j++) gel(Ker,j) = ZC_Z_mul(gel(Ker,j), p);
     521                 :      23857 :       G = ZsymM_Z_divexact(qf_apply_ZM(G, Ker), sqri(p));
     522                 :      23857 :       U = RgM_Rg_div(RgM_mul(U,Ker), p);
     523                 :      23857 :       E[i] -= 2;
     524                 :      23857 :       i--; continue; /* same p */
     525                 :            :     }
     526                 :            :     /* Now 1 <= vp = dimKer <= 2 and kernel contains no vector with norm p^2 */
     527                 :            :     /* exchanging kernel and image makes minimization easier ? */
     528                 :       9650 :     m = (n-3)/2;
     529         [ +  + ]:       9650 :     d = ZM_det(G); if (odd(m)) d = negi(d);
     530 [ +  + ][ +  + ]:       9650 :     if ((vp==1 && kronecker(gmod(gdiv(negi(d), gcoeff(G,1,1)),p), p) >= 0)
     531 [ +  + ][ +  + ]:       3525 :      || (vp==2 && odd(n) && n >= 5)
                 [ -  + ]
     532 [ +  + ][ +  - ]:       3515 :      || (vp==2 && !odd(n) && kronecker(modii(diviiexact(d,sqri(p)), p),p) < 0))
                 [ +  + ]
     533                 :            :     {
     534                 :            :       long j;
     535         [ -  + ]:       6150 :       if (DEBUGLEVEL >= 4) err_printf("    case 3\n");
     536                 :       6150 :       Ker = matid(n);
     537         [ +  + ]:      43918 :       for (j = dimKer+1; j <= n; j++) gcoeff(Ker,j,j) = p;
     538                 :       6150 :       G = ZsymM_Z_divexact(qf_apply_ZM(G, Ker), p);
     539                 :       6150 :       U = RgM_mul(U,Ker);
     540                 :       6150 :       E[i] -= 2*dimKer-n;
     541                 :       6150 :       i--; continue; /* same p */
     542                 :            :     }
     543                 :            : 
     544                 :            :     /* Minimization was not possible so far. */
     545                 :            :     /* If n == 3 or 4, this proves the local non-solubility at p. */
     546 [ +  + ][ -  + ]:       3500 :     if (n == 3 || n == 4)
     547                 :            :     {
     548         [ -  + ]:        340 :       if (DEBUGLEVEL >= 1) err_printf(" no local solution at %Ps\n",p);
     549                 :        340 :       return(p);
     550                 :            :     }
     551                 :       3160 :     vectrunc_append(faP, p);
     552                 :       3160 :     vecsmalltrunc_append(faE, vp);
     553                 :            :   }
     554         [ +  + ]:       5510 :   if (!U) U = matid(n);
     555                 :            :   else
     556                 :            :   { /* apply LLL to avoid coefficient explosion */
     557                 :       4755 :     aux = lllint(Q_primpart(U));
     558                 :       4755 :     G = qf_apply_ZM(G,aux);
     559                 :       4755 :     U = RgM_mul(U,aux);
     560                 :            :   }
     561                 :       5850 :   return mkvec4(G, U, faP, faE);
     562                 :            : }
     563                 :            : 
     564                 :            : /* CLASS GROUP COMPUTATIONS */
     565                 :            : 
     566                 :            : /* Compute the square root of the quadratic form q of discriminant D. Not
     567                 :            :  * fully implemented; it only works for detqfb squarefree except at 2, where
     568                 :            :  * the valuation is 2 or 3.
     569                 :            :  * mkmat2(P,zv_to_ZV(E)) = factor(2*abs(det q)) */
     570                 :            : static GEN
     571                 :       1611 : qfbsqrt(GEN D, GEN q, GEN P, GEN E)
     572                 :            : {
     573                 :       1611 :   GEN a = gel(q,1), b = shifti(gel(q,2),-1), c = gel(q,3), mb = negi(b);
     574                 :            :   GEN m,n, aux,Q1,M, A,B,C;
     575                 :       1611 :   GEN d = subii(mulii(a,c), sqri(b));
     576                 :            :   long i;
     577                 :            : 
     578                 :            :   /* 1st step: solve m^2 = a (d), m*n = -b (d), n^2 = c (d) */
     579                 :       1611 :   m = n = mkintmod(gen_0,gen_1);
     580                 :       1611 :   E[1] -= 3;
     581         [ +  + ]:       6963 :   for (i = 1; i < lg(P); i++)
     582                 :            :   {
     583                 :       5352 :     GEN p = gel(P,i), N, M;
     584         [ +  + ]:       5352 :     if (!E[i]) continue;
     585         [ +  + ]:       5302 :     if (dvdii(a,p)) {
     586                 :       1279 :       aux = Fp_sqrt(c, p);
     587                 :       1279 :       N = aux;
     588                 :       1279 :       M = Fp_div(mb, aux, p);
     589                 :            :     } else {
     590                 :       4023 :       aux = Fp_sqrt(a, p);
     591                 :       4023 :       M = aux;
     592                 :       4023 :       N = Fp_div(mb, aux, p);
     593                 :            :     }
     594                 :       5302 :     n = chinese(n, mkintmod(N,p));
     595                 :       5302 :     m = chinese(m, mkintmod(M,p));
     596                 :            :   }
     597                 :       1611 :   m = centerlift(m);
     598                 :       1611 :   n = centerlift(n);
     599         [ -  + ]:       1611 :   if (DEBUGLEVEL >=4) err_printf("    [m,n] = [%Ps, %Ps]\n",m,n);
     600                 :            : 
     601                 :            :   /* 2nd step: build Q1, with det=-1 such that Q1(x,y,0) = G(x,y) */
     602                 :       1611 :   A = diviiexact(subii(sqri(n),c), d);
     603                 :       1611 :   B = diviiexact(addii(b, mulii(m,n)), d);
     604                 :       1611 :   C = diviiexact(subii(sqri(m), a), d);
     605                 :       1611 :   Q1 = mkmat3(mkcol3(A,B,n), mkcol3(B,C,m), mkcol3(n,m,d));
     606                 :       1611 :   Q1 = gneg(adj(Q1));
     607                 :            : 
     608                 :            :   /* 3rd step: reduce Q1 to [0,0,-1;0,1,0;-1,0,0] */
     609                 :       1611 :   M = qflllgram_indefgoon2(Q1);
     610         [ +  + ]:       1611 :   if (signe(gel(M,1)) < 0) M = ZC_neg(M);
     611                 :       1611 :   a = gel(M,1);
     612                 :       1611 :   b = gel(M,2);
     613                 :       1611 :   c = gel(M,3);
     614         [ +  + ]:       1611 :   if (mpodd(a))
     615                 :       1571 :     return qfb(D, a, shifti(b,1), shifti(c,1));
     616                 :            :   else
     617                 :       1611 :     return qfb(D, c, shifti(negi(b),1), shifti(a,1));
     618                 :            : }
     619                 :            : 
     620                 :            : /* \prod gen[i]^e[i] as a Qfb, e in {0,1}^n non-zero */
     621                 :            : static GEN
     622                 :       2811 : qfb_factorback(GEN D, GEN gen, GEN e)
     623                 :            : {
     624                 :       2811 :   GEN q = NULL;
     625                 :       2811 :   long j, l = lg(gen), n = 0;
     626         [ +  + ]:       9590 :   for (j = 1; j < l; j++)
     627 [ +  + ][ +  + ]:       6779 :     if (e[j]) { n++; q = q? qfbcompraw(q, gel(gen,j)): gel(gen,j); }
     628         [ +  + ]:       2811 :   return (n <= 1)? q: qfbreduce(D, q);
     629                 :            : }
     630                 :            : 
     631                 :            : /* unit form, assuming 4 | D */
     632                 :            : static GEN
     633                 :        695 : id(GEN D)
     634                 :        695 : { return mkmat2(mkcol2(gen_1,gen_0),mkcol2(gen_0,shifti(negi(D),-2))); }
     635                 :            : 
     636                 :            : /* Shanks/Bosma-Stevenhagen algorithm to compute the 2-Sylow of the class
     637                 :            :  * group of discriminant D. Only works for D = fundamental discriminant.
     638                 :            :  * When D = 1(4), work with 4D.
     639                 :            :  * P2D,E2D = factor(abs(2*D))
     640                 :            :  * Pm2D = factor(-abs(2*D))[,1].
     641                 :            :  * Return a form having Witt invariants W at Pm2D */
     642                 :            : static GEN
     643                 :       1895 : quadclass2(GEN D, GEN P2D, GEN E2D, GEN Pm2D, GEN W, int n_is_4)
     644                 :            : {
     645                 :            :   GEN gen, Wgen, U2;
     646                 :            :   long i, n, r, m, vD;
     647                 :            : 
     648         [ +  + ]:       1895 :   if (mpodd(D))
     649                 :            :   {
     650                 :        150 :     D = shifti(D,2);
     651                 :        150 :     E2D = shallowcopy(E2D);
     652                 :        150 :     E2D[1] = 3;
     653                 :            :   }
     654         [ +  + ]:       1895 :   if (zv_equal0(W)) return id(D);
     655                 :            : 
     656                 :       1285 :   n = lg(Pm2D)-1; /* >= 3 since W != 0 */
     657                 :       1285 :   r = n-3;
     658         [ +  + ]:       1285 :   m = (signe(D)>0)? r+1: r;
     659                 :            :   /* n=4: look among forms of type q or 2*q, since Q can be imprimitive */
     660         [ +  + ]:       1285 :   U2 = n_is_4? mkmat(hilberts(gen_2, D, Pm2D, lg(Pm2D))): NULL;
     661 [ +  + ][ +  + ]:       1285 :   if (U2 && zv_equal(gel(U2,1),W)) return gmul2n(id(D),1);
     662                 :            : 
     663                 :       1200 :   gen = cgetg(m+1, t_VEC);
     664         [ +  + ]:       3258 :   for (i = 1; i <= m; i++) /* no need to look at Pm2D[1]=-1, nor Pm2D[2]=2 */
     665                 :            :   {
     666                 :       2058 :     GEN p = gel(Pm2D,i+2), d;
     667                 :       2058 :     long vp = Z_pvalrem(D,p, &d);
     668                 :       2058 :     gel(gen,i) = qfb(D, powiu(p,vp), gen_0, negi(shifti(d,-2)));
     669                 :            :   }
     670                 :       1200 :   vD = Z_lval(D,2);  /* = 2 or 3 */
     671 [ +  + ][ +  + ]:       1200 :   if (vD == 2 && smodis(D,16) != 4)
     672                 :            :   {
     673                 :         85 :     GEN q2 = qfb(D, gen_2,gen_2, shifti(subsi(4,D),-3));
     674                 :         85 :     m++; r++; gen = shallowconcat(gen, mkvec(q2));
     675                 :            :   }
     676         [ +  + ]:       1200 :   if (vD == 3)
     677                 :            :   {
     678                 :       1040 :     GEN q2 = qfb(D, gen_2,gen_0, negi(shifti(D,-3)));
     679                 :       1040 :     m++; r++; gen = shallowconcat(gen, mkvec(q2));
     680                 :            :   }
     681         [ -  + ]:       1200 :   if (!r) return id(D);
     682                 :       1200 :   Wgen = qflocalinvariants(gen,Pm2D);
     683                 :            :   for(;;)
     684                 :            :   {
     685                 :       2661 :     GEN Wgen2, gen2, Ker, indexim = gel(Flm_indexrank(Wgen,2), 2);
     686                 :            :     long dKer;
     687         [ +  + ]:       2661 :     if (lg(indexim)-1 >= r)
     688                 :            :     {
     689                 :       1200 :       GEN W2 = Wgen, V;
     690         [ +  + ]:       1200 :       if (lg(indexim) < lg(Wgen)) W2 = vecpermute(Wgen,indexim);
     691         [ +  + ]:       1200 :       if (U2) W2 = shallowconcat(W2,U2);
     692                 :       1200 :       V = Flm_Flc_invimage(W2, W,2);
     693         [ +  - ]:       1200 :       if (V) {
     694                 :       1200 :         GEN Q = qfb_factorback(D, vecpermute(gen,indexim), V);
     695                 :       1200 :         Q = gtomat(Q);
     696 [ +  + ][ +  + ]:       1200 :         if (U2 && V[lg(V)-1]) Q = gmul2n(Q,1);
     697                 :       1200 :         return Q;
     698                 :            :       }
     699                 :            :     }
     700                 :       1461 :     Ker = Flm_ker(Wgen,2); dKer = lg(Ker)-1;
     701                 :       1461 :     gen2 = cgetg(m+1, t_VEC);
     702                 :       1461 :     Wgen2 = cgetg(m+1, t_MAT);
     703         [ +  + ]:       3072 :     for (i = 1; i <= dKer; i++)
     704                 :            :     {
     705                 :       1611 :       GEN q = qfb_factorback(D, gen, gel(Ker,i));
     706                 :       1611 :       q = qfbsqrt(D,q,P2D,E2D);
     707                 :       1611 :       gel(gen2,i) = q;
     708                 :       1611 :       gel(Wgen2,i) = gel(qflocalinvariants(q,Pm2D), 1);
     709                 :            :     }
     710         [ +  + ]:       3151 :     for (; i <=m; i++)
     711                 :            :     {
     712                 :       1690 :       long j = indexim[i-dKer];
     713                 :       1690 :       gel(gen2,i) = gel(gen,j);
     714                 :       1690 :       gel(Wgen2,i) = gel(Wgen,j);
     715                 :            :     }
     716                 :       1461 :     gen = gen2; Wgen = Wgen2;
     717                 :       3356 :   }
     718                 :            : }
     719                 :            : 
     720                 :            : /* QUADRATIC EQUATIONS */
     721                 :            : /* is x*y = -1 ? */
     722                 :            : static int
     723                 :       3230 : both_pm1(GEN x, GEN y)
     724 [ +  + ][ +  + ]:       3230 : { return is_pm1(x) && is_pm1(y) && signe(x) == -signe(y); }
                 [ +  + ]
     725                 :            : 
     726                 :            : /* Try to solve G = 0 with small coefficients. This is proved to work if
     727                 :            :  * -  det(G) = 1, dim <= 6 and G is LLL reduced
     728                 :            :  * Returns G if no solution is found.
     729                 :            :  * Exit with a norm 0 vector if one such is found.
     730                 :            :  * If base == 1 and norm 0 is obtained, returns [H~*G*H,H] where
     731                 :            :  * the 1st column of H is a norm 0 vector */
     732                 :            : static GEN
     733                 :      15090 : qftriv(GEN G, GEN R, long base)
     734                 :            : {
     735                 :      15090 :   long n = lg(G)-1, i;
     736                 :            :   GEN s, H;
     737                 :            : 
     738                 :            :   /* case 1: A basis vector is isotropic */
     739         [ +  + ]:      55655 :   for (i = 1; i <= n; i++)
     740         [ +  + ]:      48186 :     if (!signe(gcoeff(G,i,i)))
     741                 :            :     {
     742         [ +  + ]:       7621 :       if (!base) return col_ei(n,i);
     743                 :       7183 :       H = matid(n); swap(gel(H,1), gel(H,i));
     744                 :       7183 :       return mkvec2(qf_apply_tau(G,1,i),H);
     745                 :            :     }
     746                 :            :   /* case 2: G has a block +- [1,0;0,-1] on the diagonal */
     747         [ +  + ]:      33113 :   for (i = 2; i <= n; i++)
     748 [ +  + ][ +  + ]:      27810 :     if (!signe(gcoeff(G,i-1,i)) && both_pm1(gcoeff(G,i-1,i-1),gcoeff(G,i,i)))
     749                 :            :     {
     750                 :       2166 :       s = col_ei(n,i); gel(s,i-1) = gen_m1;
     751         [ +  + ]:       2166 :       if (!base) return s;
     752                 :       2105 :       H = matid(n); gel(H,i) = gel(H,1); gel(H,1) = s;
     753                 :       2105 :       return mkvec2(qf_apply_ZM(G,H),H);
     754                 :            :     }
     755         [ +  + ]:       5303 :   if (!R) return G; /* fail */
     756                 :            :   /* case 3: a principal minor is 0 */
     757                 :        173 :   s = keri(principal_minor(G, itos(R)));
     758                 :        173 :   s = vecextend(Q_primpart(gel(s,1)), n);
     759         [ +  + ]:        173 :   if (!base) return s;
     760                 :        147 :   H = completebasis(s, 0);
     761                 :        147 :   gel(H,n) = ZC_neg(gel(H,1)); gel(H,1) = s;
     762                 :      15090 :   return mkvec2(qf_apply_ZM(G,H),H);
     763                 :            : }
     764                 :            : 
     765                 :            : /* p a prime number, G 3x3 symetric. Finds X!=0 such that X^t G X = 0 mod p.
     766                 :            :  * Allow returning a shorter X: to be completed with 0s. */
     767                 :            : static GEN
     768                 :      12819 : qfsolvemodp(GEN G, GEN p)
     769                 :            : {
     770                 :            :   GEN a,b,c,d,e,f, v1,v2,v3,v4,v5, x1,x2,x3,N1,N2,N3,s,r;
     771                 :            : 
     772                 :            :   /* principal_minor(G,3) = [a,b,d; b,c,e; d,e,f] */
     773                 :      12819 :   a = modii(gcoeff(G,1,1), p);
     774         [ +  + ]:      12819 :   if (!signe(a)) return mkcol(gen_1);
     775                 :      10883 :   v1 = a;
     776                 :      10883 :   b = modii(gcoeff(G,1,2), p);
     777                 :      10883 :   c = modii(gcoeff(G,2,2), p);
     778                 :      10883 :   v2 = modii(subii(mulii(a,c), sqri(b)), p);
     779         [ +  + ]:      10883 :   if (!signe(v2)) return mkcol2(Fp_neg(b,p), a);
     780                 :       8695 :   d = modii(gcoeff(G,1,3), p);
     781                 :       8695 :   e = modii(gcoeff(G,2,3), p);
     782                 :       8695 :   f = modii(gcoeff(G,3,3), p);
     783                 :       8695 :   v4 = modii(subii(mulii(c,d), mulii(e,b)), p);
     784                 :       8695 :   v5 = modii(subii(mulii(a,e), mulii(d,b)), p);
     785                 :       8695 :   v3 = subii(mulii(v2,f), addii(mulii(v4,d), mulii(v5,e))); /* det(G) */
     786                 :       8695 :   v3 = modii(v3, p);
     787                 :       8695 :   N1 =  Fp_neg(v2,  p);
     788                 :       8695 :   x3 = mkcol3(v4, v5, N1);
     789         [ +  + ]:       8695 :   if (!signe(v3)) return x3;
     790                 :            : 
     791                 :            :   /* now, solve in dimension 3... reduction to the diagonal case: */
     792                 :       7544 :   x1 = mkcol3(gen_1, gen_0, gen_0);
     793                 :       7544 :   x2 = mkcol3(negi(b), a, gen_0);
     794         [ +  + ]:       7544 :   if (kronecker(N1,p) == 1) return ZC_lincomb(Fp_sqrt(N1,p),gen_1,x1,x2);
     795                 :       3002 :   N2 = Fp_div(Fp_neg(v3,p), v1, p);
     796         [ +  + ]:       3002 :   if (kronecker(N2,p) == 1) return ZC_lincomb(Fp_sqrt(N2,p),gen_1,x2,x3);
     797                 :       1499 :   N3 = Fp_mul(v2, N2, p);
     798         [ +  + ]:       1499 :   if (kronecker(N3,p) == 1) return ZC_lincomb(Fp_sqrt(N3,p),gen_1,x1,x3);
     799                 :        794 :   r = gen_1;
     800                 :            :   for(;;)
     801                 :            :   {
     802                 :       1378 :     s = Fp_sub(gen_1, Fp_mul(N1,Fp_sqr(r,p),p), p);
     803         [ +  + ]:       1378 :     if (kronecker(s, p) <= 0) break;
     804                 :        584 :     r = randomi(p);
     805                 :        584 :   }
     806                 :        794 :   s = Fp_sqrt(Fp_div(s,N3,p), p);
     807                 :      12819 :   return ZC_add(x1, ZC_lincomb(r,s,x2,x3));
     808                 :            : }
     809                 :            : 
     810                 :            : /* Given a square matrix G of dimension n >= 1, */
     811                 :            : /* solves over Z the quadratic equation X^tGX = 0. */
     812                 :            : /* G is assumed to have integral coprime coefficients. */
     813                 :            : /* The solution might be a vectorv or a matrix. */
     814                 :            : /* If no solution exists, returns an integer, that can */
     815                 :            : /* be a prime p such that there is no local solution at p, */
     816                 :            : /* or -1 if there is no real solution, */
     817                 :            : /* or 0 in some rare cases. */
     818                 :            : static  GEN
     819                 :       3045 : qfsolve_i(GEN G)
     820                 :            : {
     821                 :            :   GEN M, signG, Min, U, G1, M1, G2, M2, solG2, P, E;
     822                 :            :   GEN solG1, sol, Q, d, dQ, detG2, fam2detG;
     823                 :            :   long n, np, codim, dim;
     824                 :            : 
     825         [ -  + ]:       3045 :   if (typ(G) != t_MAT) pari_err_TYPE("qfsolve", G);
     826                 :       3045 :   n = lg(G)-1;
     827         [ -  + ]:       3045 :   if (n == 0) pari_err_DOMAIN("qfsolve", "dimension" , "=", gen_0, G);
     828         [ -  + ]:       3045 :   if (n != nbrows(G)) pari_err_DIM("qfsolve");
     829                 :       3045 :   G = Q_primpart(G); RgM_check_ZM(G, "qfsolve");
     830                 :            : 
     831                 :            :   /* Trivial case: det = 0 */
     832                 :       3045 :   d = ZM_det(G);
     833         [ +  + ]:       3045 :   if (!signe(d))
     834                 :            :   {
     835         [ +  - ]:          5 :     if (n == 1) return mkcol(gen_1);
     836                 :          0 :     sol = keri(G);
     837         [ #  # ]:          0 :     if (lg(sol) == 2) sol = gel(sol,1);
     838                 :          0 :     return sol;
     839                 :            :   }
     840                 :            : 
     841                 :            :   /* Small dimension: n <= 2 */
     842         [ +  + ]:       3040 :   if (n == 1) return gen_m1;
     843         [ +  + ]:       3035 :   if (n == 2)
     844                 :            :   {
     845                 :         15 :     GEN t, a =  gcoeff(G,1,1);
     846         [ +  + ]:         15 :     if (!signe(a)) return mkcol2(gen_1, gen_0);
     847         [ +  + ]:         10 :     if (!Z_issquareall(negi(d), &t)) return gen_m1;
     848                 :         15 :     return mkcol2(subii(t,gcoeff(G,1,2)), a);
     849                 :            :   }
     850                 :            : 
     851                 :            :   /* 1st reduction of the coefficients of G */
     852                 :       3020 :   M = qflllgram_indef(G,0);
     853         [ +  + ]:       3020 :   if (typ(M) == t_COL) return M;
     854                 :       3015 :   G = gel(M,1);
     855                 :       3015 :   M = gel(M,2);
     856                 :            : 
     857                 :            :   /* real solubility */
     858                 :       3015 :   signG = ZV_to_zv(qfsign(G));
     859                 :            :   {
     860                 :       3015 :     long r =  signG[1], s = signG[2];
     861 [ +  + ][ +  + ]:       3015 :     if (!r || !s) return gen_m1;
     862         [ +  + ]:       2965 :     if (r < s) { G = ZM_neg(G); signG = mkvecsmall2(s,r);  }
     863                 :            :   }
     864                 :            : 
     865                 :            :   /* factorization of the determinant */
     866                 :       2965 :   fam2detG = Z_factor( absi(d) );
     867                 :       2965 :   P = gel(fam2detG,1);
     868                 :       2965 :   E = ZV_to_zv(gel(fam2detG,2));
     869                 :            :   /* P,E = factor(|det(G)|) */
     870                 :            : 
     871                 :            :   /* Minimization and local solubility */
     872                 :       2965 :   Min = qfminimize(G, P, E);
     873         [ +  + ]:       2965 :   if (typ(Min) == t_INT) return Min;
     874                 :            : 
     875                 :       2625 :   M = RgM_mul(M, gel(Min,2));
     876                 :       2625 :   G = gel(Min,1);
     877                 :       2625 :   P = gel(Min,3);
     878                 :       2625 :   E = gel(Min,4);
     879                 :            :   /* P,E = factor(|det(G))| */
     880                 :            : 
     881                 :            :   /* Now, we know that local solutions exist (except maybe at 2 if n==4)
     882                 :            :    * if n==3, det(G) = +-1
     883                 :            :    * if n==4, or n is odd, det(G) is squarefree.
     884                 :            :    * if n>=6, det(G) has all its valuations <=2. */
     885                 :            : 
     886                 :            :   /* Reduction of G and search for trivial solutions. */
     887                 :            :   /* When |det G|=1, such trivial solutions always exist. */
     888                 :       2625 :   U = qflllgram_indef(G,0);
     889         [ +  + ]:       2625 :   if(typ(U) == t_COL) return Q_primpart(RgM_RgC_mul(M,U));
     890                 :       2105 :   G = gel(U,1);
     891                 :       2105 :   M = RgM_mul(M, gel(U,2));
     892                 :            :   /* P,E = factor(|det(G))| */
     893                 :            : 
     894                 :            :   /* If n >= 6 is even, need to increment the dimension by 1 to suppress all
     895                 :            :    * squares from det(G) */
     896                 :       2105 :   np = lg(P)-1;
     897 [ +  + ][ +  + ]:       2105 :   if (n < 6 || odd(n) || !np)
                 [ -  + ]
     898                 :            :   {
     899                 :       1115 :     codim = 0;
     900                 :       1115 :     G1 = G;
     901                 :       1115 :     M1 = NULL;
     902                 :            :   }
     903                 :            :   else
     904                 :            :   {
     905                 :            :     GEN aux;
     906                 :            :     long i;
     907                 :        990 :     codim = 1; n++;
     908                 :            :     /* largest square divisor of d */
     909                 :        990 :     aux = gen_1;
     910         [ +  + ]:       4660 :     for (i = 1; i <= np; i++)
     911         [ +  + ]:       3670 :       if (E[i] == 2) { aux = mulii(aux, gel(P,i)); E[i] = 3; }
     912                 :            :     /* Choose sign(aux) so as to balance the signature of G1 */
     913         [ +  + ]:        990 :     if (signG[1] > signG[2])
     914                 :            :     {
     915                 :        390 :       signG[2]++;
     916                 :        390 :       aux = negi(aux);
     917                 :            :     }
     918                 :            :     else
     919                 :        600 :       signG[1]++;
     920                 :        990 :     G1 = shallowmatconcat(diagonal_shallow(mkvec2(G,aux)));
     921                 :            :     /* P,E = factor(|det G1|) */
     922                 :        990 :     Min = qfminimize(G1, P, E);
     923                 :        990 :     G1 = gel(Min,1);
     924                 :        990 :     M1 = gel(Min,2);
     925                 :        990 :     P = gel(Min,3);
     926                 :        990 :     E = gel(Min,4);
     927                 :        990 :     np = lg(P)-1;
     928                 :            :   }
     929                 :            : 
     930                 :            :   /* now, d is squarefree */
     931         [ +  + ]:       2105 :   if (!np)
     932                 :            :   { /* |d| = 1 */
     933                 :        185 :      G2 = G1;
     934                 :        185 :      M2 = NULL;
     935                 :            :   }
     936                 :            :   else
     937                 :            :   { /* |d| > 1: increment dimension by 2 */
     938                 :            :     GEN factdP, factdE, W;
     939                 :            :     long i, lfactdP;
     940                 :       1920 :     codim += 2;
     941                 :       1920 :     d = ZV_prod(P); /* d = abs(matdet(G1)); */
     942         [ +  + ]:       1920 :     if (odd(signG[2])) togglesign_safe(&d); /* d = matdet(G1); */
     943                 :            :     /* solubility at 2 (this is the only remaining bad prime). */
     944 [ +  + ][ +  + ]:       1920 :     if (n == 4 && smodis(d,8) == 1 && qflocalinvariant(G,gen_2) == 1)
                 [ +  + ]
     945                 :         25 :       return gen_2;
     946                 :            : 
     947         [ +  + ]:       1895 :     P = shallowconcat(mpodd(d)? mkvec2(NULL,gen_2): mkvec(NULL), P);
     948                 :            :     /* build a binary quadratic form with given Witt invariants */
     949                 :       1895 :     W = const_vecsmall(lg(P)-1, 0);
     950                 :            :     /* choose signature of Q (real invariant and sign of the discriminant) */
     951                 :       1895 :     dQ = absi(d);
     952         [ +  + ]:       1895 :     if (signG[1] > signG[2]) togglesign_safe(&dQ); /* signQ = [2,0]; */
     953 [ +  + ][ +  + ]:       1895 :     if (n == 4 && smodis(dQ,4) != 1) dQ = shifti(dQ,2);
     954         [ +  + ]:       1895 :     if (n >= 5) dQ = shifti(dQ,3);
     955                 :            : 
     956                 :            :     /* p-adic invariants */
     957         [ +  + ]:       1895 :     if (n == 4)
     958                 :            :     {
     959                 :        475 :       GEN t = qflocalinvariants(ZM_neg(G1),P);
     960         [ +  + ]:       1905 :       for (i = 3; i < lg(P); i++) W[i] = ucoeff(t,i,1);
     961                 :            :     }
     962                 :            :     else
     963                 :            :     {
     964         [ +  + ]:       1420 :       long s = signe(dQ) == signe(d)? 1: -1;
     965                 :            :       GEN t;
     966         [ +  + ]:       1420 :       if (odd((n-3)/2)) s = -s;
     967         [ +  + ]:       1420 :       t = s > 0? utoipos(8): utoineg(8);
     968         [ +  + ]:       4283 :       for (i = 3; i < lg(P); i++)
     969                 :       2863 :         W[i] = hilbertii(t, gel(P,i), gel(P,i)) > 0;
     970                 :            :     }
     971                 :            :     /* for p = 2, the choice is fixed from the product formula */
     972                 :       1895 :     W[2] = Flv_sum(W, 2);
     973                 :            : 
     974                 :            :     /* Construction of the 2-class group of discriminant dQ until some product
     975                 :            :      * of the generators gives the desired invariants. */
     976                 :       1895 :     factdP = vecsplice(P, 1); lfactdP =  lg(factdP);
     977                 :       1895 :     factdE = cgetg(lfactdP, t_VECSMALL);
     978         [ +  + ]:       8083 :     for (i = 1; i < lfactdP; i++) factdE[i] = Z_pval(dQ, gel(factdP,i));
     979                 :       1895 :     factdE[1]++;
     980                 :            :     /* factdP,factdE = factor(2|dQ|), P = factor(-2|dQ|)[,1] */
     981                 :       1895 :     Q = quadclass2(dQ, factdP,factdE, P, W, n == 4);
     982                 :            :     /* Build a form of dim=n+2 potentially unimodular */
     983                 :       1895 :     G2 = shallowmatconcat(diagonal_shallow(mkvec2(G1,ZM_neg(Q))));
     984                 :            :     /* Minimization of G2 */
     985                 :       1895 :     detG2 = mulii(d, ZM_det(Q));
     986         [ +  + ]:       8083 :     for (i = 1; i < lfactdP; i++) factdE[i] = Z_pval(detG2, gel(factdP,i));
     987                 :            :     /* factdP,factdE = factor(|det G2|) */
     988                 :       1895 :     Min = qfminimize(G2, factdP,factdE);
     989                 :       1895 :     M2 = gel(Min,2);
     990                 :       1895 :     G2 = gel(Min,1);
     991                 :            :   }
     992                 :            :   /* |det(G2)| = 1, find a totally isotropic subspace for G2 */
     993                 :       2080 :   solG2 = qflllgram_indefgoon(G2);
     994                 :            :   /* G2 must have a subspace of solutions of dimension > codim */
     995                 :       2080 :   dim = codim+2;
     996         [ +  + ]:       3044 :   while(gequal0(principal_minor(gel(solG2,1), dim))) dim ++;
     997                 :       2080 :   solG2 = vecslice(gel(solG2,2), 1, dim-1);
     998                 :            : 
     999         [ +  + ]:       2080 :   if (!M2)
    1000                 :        185 :     solG1 = solG2;
    1001                 :            :   else
    1002                 :            :   { /* solution of G1 is simultaneously in solG2 and x[n+1] = x[n+2] = 0*/
    1003                 :            :     GEN K;
    1004                 :       1895 :     solG1 = RgM_mul(M2,solG2);
    1005                 :       1895 :     K = ker(rowslice(solG1,n+1,n+2));
    1006                 :       1895 :     solG1 = RgM_mul(rowslice(solG1,1,n), K);
    1007                 :            :   }
    1008         [ +  + ]:       2080 :   if (!M1)
    1009                 :       1090 :     sol = solG1;
    1010                 :            :   else
    1011                 :            :   { /* solution of G1 is simultaneously in solG2 and x[n] = 0 */
    1012                 :            :     GEN K;
    1013                 :        990 :     sol = RgM_mul(M1,solG1);
    1014                 :        990 :     K = ker(rowslice(sol,n,n));
    1015                 :        990 :     sol = RgM_mul(rowslice(sol,1,n-1), K);
    1016                 :            :   }
    1017                 :       2080 :   sol = Q_primpart(RgM_mul(M, sol));
    1018         [ +  + ]:       2080 :   if (lg(sol) == 2) sol = gel(sol,1);
    1019                 :       3045 :   return sol;
    1020                 :            : }
    1021                 :            : GEN
    1022                 :       3045 : qfsolve(GEN G)
    1023                 :            : {
    1024                 :       3045 :   pari_sp av = avma;
    1025                 :       3045 :   return gerepilecopy(av, qfsolve_i(G));
    1026                 :            : }
    1027                 :            : 
    1028                 :            : /* G is a symmetric 3x3 matrix, and sol a solution of sol~*G*sol=0.
    1029                 :            :  * Returns a parametrization of the solutions with the good invariants,
    1030                 :            :  * as a matrix 3x3, where each line contains
    1031                 :            :  * the coefficients of each of the 3 quadratic forms.
    1032                 :            :  * If fl!=0, the fl-th form is reduced. */
    1033                 :            : GEN
    1034                 :         20 : qfparam(GEN G, GEN sol, long fl)
    1035                 :            : {
    1036                 :         20 :   pari_sp av = avma;
    1037                 :            :   GEN U, G1, G2, a, b, c, d, e;
    1038                 :            :   long n;
    1039                 :            : 
    1040         [ -  + ]:         20 :   if (typ(G) != t_MAT) pari_err_TYPE("qfsolve", G);
    1041         [ +  + ]:         20 :   if (typ(sol) != t_COL) pari_err_TYPE("qfsolve", G);
    1042                 :         15 :   n = lg(G)-1;
    1043         [ -  + ]:         15 :   if (n == 0) pari_err_DOMAIN("qfsolve", "dimension" , "=", gen_0, G);
    1044 [ +  - ][ +  - ]:         15 :   if (n != nbrows(G) || n != 3 || lg(sol) != 4) pari_err_DIM("qfsolve");
                 [ -  + ]
    1045                 :         15 :   sol = Q_primpart(sol);
    1046                 :            :   /* build U such that U[,3] = sol, and |det(U)| = 1 */
    1047                 :         15 :   U = completebasis(sol,1);
    1048                 :         15 :   G1 = qf_apply_ZM(G,U); /* G1 has a 0 at the bottom right corner */
    1049                 :         15 :   a = shifti(gcoeff(G1,1,2),1);
    1050                 :         15 :   b = shifti(negi(gcoeff(G1,1,3)),1);
    1051                 :         15 :   c = shifti(negi(gcoeff(G1,2,3)),1);
    1052                 :         15 :   d = gcoeff(G1,1,1);
    1053                 :         15 :   e = gcoeff(G1,2,2);
    1054                 :         15 :   G2 = mkmat3(mkcol3(b,gen_0,d), mkcol3(c,b,a), mkcol3(gen_0,c,e));
    1055                 :         15 :   sol = ZM_mul(U,G2);
    1056         [ +  + ]:         15 :   if (fl)
    1057                 :            :   {
    1058                 :         10 :     GEN v = row(sol,fl);
    1059                 :         10 :     a = gel(v,1);
    1060                 :         10 :     b = gmul2n(gel(v,2),-1);
    1061                 :         10 :     c = gel(v,3);
    1062                 :         10 :     U = qflllgram_indef(mkmat2(mkcol2(a,b),mkcol2(b,c)), 1);
    1063                 :         10 :     U = gel(U,2);
    1064                 :         10 :     a = gcoeff(U,1,1); b = gcoeff(U,1,2);
    1065                 :         10 :     c = gcoeff(U,2,1); d = gcoeff(U,2,2);
    1066                 :         10 :     U = mkmat3(mkcol3(sqri(a),mulii(a,c),sqri(c)),
    1067                 :            :                mkcol3(shifti(mulii(a,b),1), addii(mulii(a,d),mulii(b,c)),
    1068                 :            :                       shifti(mulii(c,d),1)),
    1069                 :            :                mkcol3(sqri(b),mulii(b,d),sqri(d)));
    1070                 :         10 :     sol = ZM_mul(sol,U);
    1071                 :            :   }
    1072                 :         15 :   return gerepileupto(av, sol);
    1073                 :            : }

Generated by: LCOV version 1.9