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 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 - bibli1.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.10.0 lcov report (development 20777-d2a9243) Lines: 988 1036 95.4 %
Date: 2017-06-25 05:59:24 Functions: 59 65 90.8 %
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             : /**                 LLL Algorithm and close friends                **/
      17             : /**                                                                **/
      18             : /********************************************************************/
      19             : #include "pari.h"
      20             : #include "paripriv.h"
      21             : 
      22             : /********************************************************************/
      23             : /**             QR Factorization via Householder matrices          **/
      24             : /********************************************************************/
      25             : static int
      26     1037472 : no_prec_pb(GEN x)
      27             : {
      28     2834314 :   return (typ(x) != t_REAL || realprec(x) > LOWDEFAULTPREC
      29     1039898 :                            || expo(x) < BITS_IN_LONG/2);
      30             : }
      31             : /* Find a Householder transformation which, applied to x[k..#x], zeroes
      32             :  * x[k+1..#x]; fill L = (mu_{i,j}). Return 0 if precision problem [obtained
      33             :  * a 0 vector], 1 otherwise */
      34             : static int
      35     1037472 : FindApplyQ(GEN x, GEN L, GEN B, long k, GEN Q, long prec)
      36             : {
      37     1037472 :   long i, lx = lg(x)-1;
      38     1037472 :   GEN x2, x1, xd = x + (k-1);
      39             : 
      40     1037472 :   x1 = gel(xd,1);
      41     1037472 :   x2 = mpsqr(x1);
      42     1037472 :   if (k < lx)
      43             :   {
      44      690849 :     long lv = lx - (k-1) + 1;
      45      690849 :     GEN beta, Nx, v = cgetg(lv, t_VEC);
      46     2392978 :     for (i=2; i<lv; i++) {
      47     1702129 :       x2 = mpadd(x2, mpsqr(gel(xd,i)));
      48     1702129 :       gel(v,i) = gel(xd,i);
      49             :     }
      50      690849 :     if (!signe(x2)) return 0;
      51      690849 :     Nx = gsqrt(x2, prec); if (signe(x1) < 0) setsigne(Nx, -1);
      52      690849 :     gel(v,1) = mpadd(x1, Nx);
      53             : 
      54      690849 :     if (!signe(x1))
      55         120 :       beta = gtofp(x2, prec); /* make sure typ(beta) != t_INT */
      56             :     else
      57      690729 :       beta = mpadd(x2, mpmul(Nx,x1));
      58      690849 :     gel(Q,k) = mkvec2(invr(beta), v);
      59             : 
      60      690849 :     togglesign(Nx);
      61      690849 :     gcoeff(L,k,k) = Nx;
      62             :   }
      63             :   else
      64      346623 :     gcoeff(L,k,k) = gel(x,k);
      65     1037472 :   gel(B,k) = x2;
      66     1037472 :   for (i=1; i<k; i++) gcoeff(L,k,i) = gel(x,i);
      67     1037472 :   return no_prec_pb(x2);
      68             : }
      69             : 
      70             : /* apply Householder transformation Q = [beta,v] to r with t_INT/t_REAL
      71             :  * coefficients, in place: r -= ((0|v).r * beta) v */
      72             : static void
      73     1702111 : ApplyQ(GEN Q, GEN r)
      74             : {
      75     1702111 :   GEN s, rd, beta = gel(Q,1), v = gel(Q,2);
      76     1702111 :   long i, l = lg(v), lr = lg(r);
      77             : 
      78     1702111 :   rd = r + (lr - l);
      79     1702111 :   s = mpmul(gel(v,1), gel(rd,1));
      80     1702111 :   for (i=2; i<l; i++) s = mpadd(s, mpmul(gel(v,i) ,gel(rd,i)));
      81     1702111 :   s = mpmul(beta, s);
      82    12162681 :   for (i=1; i<l; i++)
      83    10460570 :     if (signe(gel(v,i))) gel(rd,i) = mpsub(gel(rd,i), mpmul(s, gel(v,i)));
      84     1702111 : }
      85             : /* apply Q[1], ..., Q[j-1] to r */
      86             : static GEN
      87      690843 : ApplyAllQ(GEN Q, GEN r, long j)
      88             : {
      89      690843 :   pari_sp av = avma;
      90             :   long i;
      91      690843 :   r = leafcopy(r);
      92      690843 :   for (i=1; i<j; i++) ApplyQ(gel(Q,i), r);
      93      690843 :   return gerepilecopy(av, r);
      94             : }
      95             : 
      96             : /* same, arbitrary coefficients [20% slower for t_REAL at DEFAULTPREC] */
      97             : static void
      98       22113 : RgC_ApplyQ(GEN Q, GEN r)
      99             : {
     100       22113 :   GEN s, rd, beta = gel(Q,1), v = gel(Q,2);
     101       22113 :   long i, l = lg(v), lr = lg(r);
     102             : 
     103       22113 :   rd = r + (lr - l);
     104       22113 :   s = gmul(gel(v,1), gel(rd,1));
     105       22113 :   for (i=2; i<l; i++) s = gadd(s, gmul(gel(v,i) ,gel(rd,i)));
     106       22113 :   s = gmul(beta, s);
     107      486486 :   for (i=1; i<l; i++)
     108      464373 :     if (signe(gel(v,i))) gel(rd,i) = gsub(gel(rd,i), gmul(s, gel(v,i)));
     109       22113 : }
     110             : static GEN
     111         567 : RgC_ApplyAllQ(GEN Q, GEN r, long j)
     112             : {
     113         567 :   pari_sp av = avma;
     114             :   long i;
     115         567 :   r = leafcopy(r);
     116         567 :   for (i=1; i<j; i++) RgC_ApplyQ(gel(Q,i), r);
     117         567 :   return gerepilecopy(av, r);
     118             : }
     119             : 
     120             : int
     121          14 : RgM_QR_init(GEN x, GEN *pB, GEN *pQ, GEN *pL, long prec)
     122             : {
     123          14 :   x = RgM_gtomp(x, prec);
     124          14 :   return QR_init(x, pB, pQ, pL, prec);
     125             : }
     126             : 
     127             : static void
     128          35 : check_householder(GEN Q)
     129             : {
     130          35 :   long i, l = lg(Q);
     131          35 :   if (typ(Q) != t_VEC) pari_err_TYPE("mathouseholder", Q);
     132        1120 :   for (i = 1; i < l; i++)
     133             :   {
     134        1092 :     GEN q = gel(Q,i), v;
     135        1092 :     if (typ(q) != t_VEC || lg(q) != 3) pari_err_TYPE("mathouseholder", Q);
     136        1092 :     v = gel(q,2);
     137        1092 :     if (typ(v) != t_VEC || lg(v)+i-2 != l) pari_err_TYPE("mathouseholder", Q);
     138             :   }
     139          28 : }
     140             : 
     141             : GEN
     142          35 : mathouseholder(GEN Q, GEN v)
     143             : {
     144          35 :   long l = lg(Q);
     145          35 :   check_householder(Q);
     146          28 :   switch(typ(v))
     147             :   {
     148             :     case t_MAT:
     149             :     {
     150             :       long lx, i;
     151          14 :       GEN M = cgetg_copy(v, &lx);
     152          14 :       if (lx == 1) return M;
     153          14 :       if (lgcols(v) != l+1) pari_err_TYPE("mathouseholder", v);
     154          14 :       for (i = 1; i < lx; i++) gel(M,i) = RgC_ApplyAllQ(Q, gel(v,i), l);
     155          14 :       return M;
     156             :     }
     157           7 :     case t_COL: if (lg(v) == l+1) break;
     158             :       /* fall through */
     159           7 :     default: pari_err_TYPE("mathouseholder", v);
     160             :   }
     161           7 :   return RgC_ApplyAllQ(Q, v, l);
     162             : }
     163             : 
     164             : GEN
     165          28 : matqr(GEN x, long flag, long prec)
     166             : {
     167          28 :   pari_sp av = avma;
     168             :   GEN B, Q, L;
     169          28 :   long n = lg(x)-1;
     170          28 :   if (typ(x) != t_MAT) pari_err_TYPE("matqr",x);
     171          28 :   if (!n)
     172             :   {
     173          14 :     if (!flag) retmkvec2(cgetg(1,t_MAT),cgetg(1,t_MAT));
     174           7 :     retmkvec2(cgetg(1,t_VEC),cgetg(1,t_MAT));
     175             :   }
     176          14 :   if (n != nbrows(x)) pari_err_DIM("matqr");
     177          14 :   if (!RgM_QR_init(x, &B,&Q,&L, prec)) pari_err_PREC("matqr");
     178          14 :   if (!flag) Q = shallowtrans(mathouseholder(Q, matid(n)));
     179          14 :   return gerepilecopy(av, mkvec2(Q, shallowtrans(L)));
     180             : }
     181             : 
     182             : /* compute B = | x[k] |^2, Q = Householder transforms and L = mu_{i,j} */
     183             : int
     184      346629 : QR_init(GEN x, GEN *pB, GEN *pQ, GEN *pL, long prec)
     185             : {
     186      346629 :   long j, k = lg(x)-1;
     187      346629 :   GEN B = cgetg(k+1, t_VEC), Q = cgetg(k, t_VEC), L = zeromatcopy(k,k);
     188     1384089 :   for (j=1; j<=k; j++)
     189             :   {
     190     1037472 :     GEN r = gel(x,j);
     191     1037472 :     if (j > 1) r = ApplyAllQ(Q, r, j);
     192     1037472 :     if (!FindApplyQ(r, L, B, j, Q, prec)) return 0;
     193             :   }
     194      346617 :   *pB = B; *pQ = Q; *pL = L; return 1;
     195             : }
     196             : /* x a square t_MAT with t_INT / t_REAL entries and maximal rank. Return
     197             :  * qfgaussred(x~*x) */
     198             : GEN
     199      346168 : gaussred_from_QR(GEN x, long prec)
     200             : {
     201      346168 :   long j, k = lg(x)-1;
     202             :   GEN B, Q, L;
     203      346168 :   if (!QR_init(x, &B,&Q,&L, prec)) return NULL;
     204     1034106 :   for (j=1; j<k; j++)
     205             :   {
     206      687938 :     GEN m = gel(L,j), invNx = invr(gel(m,j));
     207             :     long i;
     208      687938 :     gel(m,j) = gel(B,j);
     209      687938 :     for (i=j+1; i<=k; i++) gel(m,i) = mpmul(invNx, gel(m,i));
     210             :   }
     211      346168 :   gcoeff(L,j,j) = gel(B,j);
     212      346168 :   return shallowtrans(L);
     213             : }
     214             : GEN
     215         447 : R_from_QR(GEN x, long prec)
     216             : {
     217             :   GEN B, Q, L;
     218         447 :   if (!QR_init(x, &B,&Q,&L, prec)) return NULL;
     219         435 :   return shallowtrans(L);
     220             : }
     221             : 
     222             : /********************************************************************/
     223             : /**             QR Factorization via Gram-Schmidt                  **/
     224             : /********************************************************************/
     225             : 
     226             : /* return Gram-Schmidt orthogonal basis (f) attached to (e), B is the
     227             :  * vector of the (f_i . f_i)*/
     228             : GEN
     229         434 : RgM_gram_schmidt(GEN e, GEN *ptB)
     230             : {
     231         434 :   long i,j,lx = lg(e);
     232         434 :   GEN f = RgM_shallowcopy(e), B, iB;
     233             : 
     234         434 :   B = cgetg(lx, t_VEC);
     235         434 :   iB= cgetg(lx, t_VEC);
     236             : 
     237        1435 :   for (i=1;i<lx;i++)
     238             :   {
     239        1001 :     GEN p1 = NULL;
     240        1001 :     pari_sp av = avma;
     241        2205 :     for (j=1; j<i; j++)
     242             :     {
     243        1204 :       GEN mu = gmul(RgV_dotproduct(gel(e,i),gel(f,j)), gel(iB,j));
     244        1204 :       GEN p2 = gmul(mu, gel(f,j));
     245        1204 :       p1 = p1? gadd(p1,p2): p2;
     246             :     }
     247        1001 :     p1 = p1? gerepileupto(av, gsub(gel(e,i), p1)): gel(e,i);
     248        1001 :     gel(f,i) = p1;
     249        1001 :     gel(B,i) = RgV_dotsquare(gel(f,i));
     250        1001 :     gel(iB,i) = ginv(gel(B,i));
     251             :   }
     252         434 :   *ptB = B; return f;
     253             : }
     254             : 
     255             : /* Assume B an LLL-reduced basis, t a vector. Apply Babai's nearest plane
     256             :  * algorithm to (B,t) */
     257             : GEN
     258         434 : RgM_Babai(GEN B, GEN t)
     259             : {
     260         434 :   GEN C, N, G = RgM_gram_schmidt(B, &N), b = t;
     261         434 :   long j, n = lg(B)-1;
     262             : 
     263         434 :   C = cgetg(n+1,t_COL);
     264        1435 :   for (j = n; j > 0; j--)
     265             :   {
     266        1001 :     GEN c = gdiv( RgV_dotproduct(b, gel(G,j)), gel(N,j) );
     267             :     long e;
     268        1001 :     c = grndtoi(c,&e);
     269        1001 :     if (e >= 0) return NULL;
     270        1001 :     if (signe(c)) b = RgC_sub(b, RgC_Rg_mul(gel(G,j), c));
     271        1001 :     gel(C,j) = c;
     272             :   }
     273         434 :   return C;
     274             : }
     275             : 
     276             : /********************************************************************/
     277             : /**                                                                **/
     278             : /**                          LLL ALGORITHM                         **/
     279             : /**                                                                **/
     280             : /********************************************************************/
     281             : /* Def: a matrix M is said to be -partially reduced- if | m1 +- m2 | >= |m1|
     282             :  * for any two columns m1 != m2, in M.
     283             :  *
     284             :  * Input: an integer matrix mat whose columns are linearly independent. Find
     285             :  * another matrix T such that mat * T is partially reduced.
     286             :  *
     287             :  * Output: mat * T if flag = 0;  T if flag != 0,
     288             :  *
     289             :  * This routine is designed to quickly reduce lattices in which one row
     290             :  * is huge compared to the other rows.  For example, when searching for a
     291             :  * polynomial of degree 3 with root a mod N, the four input vectors might
     292             :  * be the coefficients of
     293             :  *     X^3 - (a^3 mod N), X^2 - (a^2 mod N), X - (a mod N), N.
     294             :  * All four constant coefficients are O(p) and the rest are O(1). By the
     295             :  * pigeon-hole principle, the coefficients of the smallest vector in the
     296             :  * lattice are O(p^(1/4)), hence significant reduction of vector lengths
     297             :  * can be anticipated.
     298             :  *
     299             :  * An improved algorithm would look only at the leading digits of dot*.  It
     300             :  * would use single-precision calculations as much as possible.
     301             :  *
     302             :  * Original code: Peter Montgomery (1994) */
     303             : static GEN
     304          35 : lllintpartialall(GEN m, long flag)
     305             : {
     306          35 :   const long ncol = lg(m)-1;
     307          35 :   const pari_sp av = avma;
     308             :   GEN tm1, tm2, mid;
     309             : 
     310          35 :   if (ncol <= 1) return flag? matid(ncol): gcopy(m);
     311             : 
     312          14 :   tm1 = flag? matid(ncol): NULL;
     313             :   {
     314          14 :     const pari_sp av2 = avma;
     315          14 :     GEN dot11 = ZV_dotsquare(gel(m,1));
     316          14 :     GEN dot22 = ZV_dotsquare(gel(m,2));
     317          14 :     GEN dot12 = ZV_dotproduct(gel(m,1), gel(m,2));
     318          14 :     GEN tm  = matid(2); /* For first two columns only */
     319             : 
     320          14 :     int progress = 0;
     321          14 :     long npass2 = 0;
     322             : 
     323             : /* Row reduce the first two columns of m. Our best result so far is
     324             :  * (first two columns of m)*tm.
     325             :  *
     326             :  * Initially tm = 2 x 2 identity matrix.
     327             :  * Inner products of the reduced matrix are in dot11, dot12, dot22. */
     328          63 :     while (npass2 < 2 || progress)
     329             :     {
     330          42 :       GEN dot12new, q = diviiround(dot12, dot22);
     331             : 
     332          35 :       npass2++; progress = signe(q);
     333          35 :       if (progress)
     334             :       {/* Conceptually replace (v1, v2) by (v1 - q*v2, v2), where v1 and v2
     335             :         * represent the reduced basis for the first two columns of the matrix.
     336             :         * We do this by updating tm and the inner products. */
     337          21 :         togglesign(q);
     338          21 :         dot12new = addii(dot12, mulii(q, dot22));
     339          21 :         dot11 = addii(dot11, mulii(q, addii(dot12, dot12new)));
     340          21 :         dot12 = dot12new;
     341          21 :         ZC_lincomb1_inplace(gel(tm,1), gel(tm,2), q);
     342             :       }
     343             : 
     344             :       /* Interchange the output vectors v1 and v2.  */
     345          35 :       swap(dot11,dot22);
     346          35 :       swap(gel(tm,1), gel(tm,2));
     347             : 
     348             :       /* Occasionally (including final pass) do garbage collection.  */
     349          35 :       if ((npass2 & 0xff) == 0 || !progress)
     350          14 :         gerepileall(av2, 4, &dot11,&dot12,&dot22,&tm);
     351             :     } /* while npass2 < 2 || progress */
     352             : 
     353             :     {
     354             :       long i;
     355           7 :       GEN det12 = subii(mulii(dot11, dot22), sqri(dot12));
     356             : 
     357           7 :       mid = cgetg(ncol+1, t_MAT);
     358          21 :       for (i = 1; i <= 2; i++)
     359             :       {
     360          14 :         GEN tmi = gel(tm,i);
     361          14 :         if (tm1)
     362             :         {
     363          14 :           GEN tm1i = gel(tm1,i);
     364          14 :           gel(tm1i,1) = gel(tmi,1);
     365          14 :           gel(tm1i,2) = gel(tmi,2);
     366             :         }
     367          14 :         gel(mid,i) = ZC_lincomb(gel(tmi,1),gel(tmi,2), gel(m,1),gel(m,2));
     368             :       }
     369          42 :       for (i = 3; i <= ncol; i++)
     370             :       {
     371          35 :         GEN c = gel(m,i);
     372          35 :         GEN dot1i = ZV_dotproduct(gel(mid,1), c);
     373          35 :         GEN dot2i = ZV_dotproduct(gel(mid,2), c);
     374             :        /* ( dot11  dot12 ) (q1)   ( dot1i )
     375             :         * ( dot12  dot22 ) (q2) = ( dot2i )
     376             :         *
     377             :         * Round -q1 and -q2 to nearest integer. Then compute
     378             :         *   c - q1*mid[1] - q2*mid[2].
     379             :         * This will be approximately orthogonal to the first two vectors in
     380             :         * the new basis. */
     381          35 :         GEN q1neg = subii(mulii(dot12, dot2i), mulii(dot22, dot1i));
     382          35 :         GEN q2neg = subii(mulii(dot12, dot1i), mulii(dot11, dot2i));
     383             : 
     384          35 :         q1neg = diviiround(q1neg, det12);
     385          35 :         q2neg = diviiround(q2neg, det12);
     386          35 :         if (tm1)
     387             :         {
     388          70 :           gcoeff(tm1,1,i) = addii(mulii(q1neg, gcoeff(tm,1,1)),
     389          35 :                                   mulii(q2neg, gcoeff(tm,1,2)));
     390          70 :           gcoeff(tm1,2,i) = addii(mulii(q1neg, gcoeff(tm,2,1)),
     391          35 :                                   mulii(q2neg, gcoeff(tm,2,2)));
     392             :         }
     393          35 :         gel(mid,i) = ZC_add(c, ZC_lincomb(q1neg,q2neg, gel(mid,1),gel(mid,2)));
     394             :       } /* for i */
     395             :     } /* local block */
     396             :   }
     397           7 :   if (DEBUGLEVEL>6)
     398             :   {
     399           0 :     if (tm1) err_printf("tm1 = %Ps",tm1);
     400           0 :     err_printf("mid = %Ps",mid); /* = m * tm1 */
     401             :   }
     402           7 :   gerepileall(av, tm1? 2: 1, &mid, &tm1);
     403             :   {
     404             :    /* For each pair of column vectors v and w in mid * tm2,
     405             :     * try to replace (v, w) by (v, v - q*w) for some q.
     406             :     * We compute all inner products and check them repeatedly. */
     407           7 :     const pari_sp av3 = avma;
     408           7 :     long i, j, npass = 0, e = LONG_MAX;
     409           7 :     GEN dot = cgetg(ncol+1, t_MAT); /* scalar products */
     410             : 
     411           7 :     tm2 = matid(ncol);
     412          56 :     for (i=1; i <= ncol; i++)
     413             :     {
     414          49 :       gel(dot,i) = cgetg(ncol+1,t_COL);
     415         245 :       for (j=1; j <= i; j++)
     416         196 :         gcoeff(dot,j,i) = gcoeff(dot,i,j) = ZV_dotproduct(gel(mid,i),gel(mid,j));
     417             :     }
     418             :     for(;;)
     419             :     {
     420          42 :       long reductions = 0, olde = e;
     421         336 :       for (i=1; i <= ncol; i++)
     422             :       {
     423             :         long ijdif;
     424        2058 :         for (ijdif=1; ijdif < ncol; ijdif++)
     425             :         {
     426             :           long d, k1, k2;
     427             :           GEN codi, q;
     428             : 
     429        1764 :           j = i + ijdif; if (j > ncol) j -= ncol;
     430             :           /* let k1, resp. k2,  index of larger, resp. smaller, column */
     431        1764 :           if (cmpii(gcoeff(dot,i,i), gcoeff(dot,j,j)) > 0) { k1 = i; k2 = j; }
     432        1022 :           else                                             { k1 = j; k2 = i; }
     433        1764 :           codi = gcoeff(dot,k2,k2);
     434        1764 :           q = signe(codi)? diviiround(gcoeff(dot,k1,k2), codi): gen_0;
     435        1764 :           if (!signe(q)) continue;
     436             : 
     437             :           /* Try to subtract a multiple of column k2 from column k1.  */
     438         700 :           reductions++; togglesign_safe(&q);
     439         700 :           ZC_lincomb1_inplace(gel(tm2,k1), gel(tm2,k2), q);
     440         700 :           ZC_lincomb1_inplace(gel(dot,k1), gel(dot,k2), q);
     441        1400 :           gcoeff(dot,k1,k1) = addii(gcoeff(dot,k1,k1),
     442         700 :                                     mulii(q, gcoeff(dot,k2,k1)));
     443         700 :           for (d = 1; d <= ncol; d++) gcoeff(dot,k1,d) = gcoeff(dot,d,k1);
     444             :         } /* for ijdif */
     445         294 :         if (gc_needed(av3,2))
     446             :         {
     447           0 :           if(DEBUGMEM>1) pari_warn(warnmem,"lllintpartialall");
     448           0 :           gerepileall(av3, 2, &dot,&tm2);
     449             :         }
     450             :       } /* for i */
     451          42 :       if (!reductions) break;
     452          35 :       e = 0;
     453          35 :       for (i = 1; i <= ncol; i++) e += expi( gcoeff(dot,i,i) );
     454          35 :       if (e == olde) break;
     455          35 :       if (DEBUGLEVEL>6)
     456             :       {
     457           0 :         npass++;
     458           0 :         err_printf("npass = %ld, red. last time = %ld, log_2(det) ~ %ld\n\n",
     459             :                     npass, reductions, e);
     460             :       }
     461          35 :     } /* for(;;)*/
     462             : 
     463             :    /* Sort columns so smallest comes first in m * tm1 * tm2.
     464             :     * Use insertion sort. */
     465          49 :     for (i = 1; i < ncol; i++)
     466             :     {
     467          42 :       long j, s = i;
     468             : 
     469         189 :       for (j = i+1; j <= ncol; j++)
     470         147 :         if (cmpii(gcoeff(dot,s,s),gcoeff(dot,j,j)) > 0) s = j;
     471          42 :       if (i != s)
     472             :       { /* Exchange with proper column; only the diagonal of dot is updated */
     473          28 :         swap(gel(tm2,i), gel(tm2,s));
     474          28 :         swap(gcoeff(dot,i,i), gcoeff(dot,s,s));
     475             :       }
     476             :     }
     477             :   } /* local block */
     478           7 :   return gerepileupto(av, ZM_mul(tm1? tm1: mid, tm2));
     479             : }
     480             : 
     481             : GEN
     482          35 : lllintpartial(GEN mat) { return lllintpartialall(mat,1); }
     483             : 
     484             : GEN
     485           0 : lllintpartial_inplace(GEN mat) { return lllintpartialall(mat,0); }
     486             : 
     487             : /********************************************************************/
     488             : /**                                                                **/
     489             : /**                    COPPERSMITH ALGORITHM                       **/
     490             : /**           Finding small roots of univariate equations.         **/
     491             : /**                                                                **/
     492             : /********************************************************************/
     493             : 
     494             : static int
     495         882 : check_condition(double beta, double tau, double rho, long d, long delta, long t)
     496             : {
     497         882 :   long dim = d*delta + t;
     498        1764 :   double cond = d*delta*(delta+1)/2 - beta*delta*dim
     499         882 :     + rho*delta*(delta - 1) / 2
     500         882 :     + rho * t * delta + tau*dim*(dim - 1)/2;
     501             : 
     502         882 :   if (DEBUGLEVEL >= 4)
     503           0 :     err_printf("delta = %d, t = %d, cond = %.1lf\n", delta, t, cond);
     504             : 
     505         882 :   return (cond <= 0);
     506             : }
     507             : 
     508             : static void
     509          14 : choose_params(GEN P, GEN N, GEN X, GEN B, long *pdelta, long *pt)
     510             : {
     511          14 :   long d = degpol(P);
     512          14 :   GEN P0 = leading_coeff(P);
     513          14 :   double logN = gtodouble(glog(N, DEFAULTPREC));
     514             :   double tau, beta, rho;
     515             :   long delta, t;
     516          14 :   tau = gtodouble(glog(X, DEFAULTPREC)) / logN;
     517          14 :   beta = B? gtodouble(glog(B, DEFAULTPREC)) / logN: 1.;
     518          14 :   if (tau >= beta * beta / d)
     519           0 :     pari_err_OVERFLOW("zncoppersmith [bound too large]");
     520             :   /* TODO : remove P0 completely ! */
     521          14 :   rho = gtodouble(glog(P0, DEFAULTPREC)) / logN;
     522             : 
     523             :   /* Enumerate (delta,t) by increasing dimension of resulting lattice.
     524             :    * Not subtle, but o(1) for computing time */
     525          14 :   t = d; delta = 0;
     526             :   for(;;)
     527             :   {
     528         175 :     t += d * delta + 1; delta = 0;
     529        1218 :     while (t >= 0) {
     530         882 :       if (check_condition(beta, tau, rho, d, delta, t)) {
     531          28 :         *pdelta = delta; *pt = t; return;
     532             :       }
     533         868 :       delta++; t -= d;
     534             :     }
     535         161 :   }
     536             : }
     537             : 
     538             : static int
     539       14021 : sol_OK(GEN x, GEN N, GEN B)
     540       14021 : { return B? (cmpii(gcdii(x,N),B) >= 0): !signe(remii(x,N)); }
     541             : /* deg(P) > 0, x >= 0. Find all j such that gcd(P(j), N) >= B, |j| <= x */
     542             : static GEN
     543           7 : do_exhaustive(GEN P, GEN N, long x, GEN B)
     544             : {
     545           7 :   GEN Pe, Po, sol = vecsmalltrunc_init(2*x + 2);
     546             :   pari_sp av;
     547             :   long j;
     548           7 :   RgX_even_odd(P, &Pe,&Po); av = avma;
     549           7 :   if (sol_OK(gel(P,2), N,B)) vecsmalltrunc_append(sol, 0);
     550        7007 :   for (j = 1; j <= x; j++, avma = av)
     551             :   {
     552        7000 :     GEN j2 = sqru(j), E = FpX_eval(Pe,j2,N), O = FpX_eval(Po,j2,N);
     553        7000 :     if (sol_OK(addmuliu(E,O,j), N,B)) vecsmalltrunc_append(sol, j);
     554        7000 :     if (sol_OK(submuliu(E,O,j), N,B)) vecsmalltrunc_append(sol,-j);
     555             :   }
     556           7 :   vecsmall_sort(sol); return zv_to_ZV(sol);
     557             : }
     558             : 
     559             : /* General Coppersmith, look for a root x0 <= p, p >= B, p | N, |x0| <= X.
     560             :  * B = N coded as NULL */
     561             : GEN
     562          21 : zncoppersmith(GEN P0, GEN N, GEN X, GEN B)
     563             : {
     564             :   GEN Q, R, N0, M, sh, short_pol, *Xpowers, sol, nsp, P, Z;
     565          21 :   long delta, i, j, row, d, l, dim, t, bnd = 10;
     566          21 :   const ulong X_SMALL = 1000;
     567             : 
     568          21 :   pari_sp av = avma;
     569             : 
     570          21 :   if (typ(P0) != t_POL) pari_err_TYPE("zncoppersmith",P0);
     571          21 :   if (typ(N) != t_INT) pari_err_TYPE("zncoppersmith",N);
     572          21 :   if (typ(X) != t_INT) {
     573           7 :     X = gfloor(X);
     574           7 :     if (typ(X) != t_INT) pari_err_TYPE("zncoppersmith",X);
     575             :   }
     576          21 :   if (signe(X) < 0) pari_err_DOMAIN("zncoppersmith", "X", "<", gen_0, X);
     577          21 :   d = degpol(P0);
     578          21 :   if (d == 0) { avma = av; return cgetg(1, t_VEC); }
     579          21 :   if (d < 0) pari_err_ROOTS0("zncoppersmith");
     580          21 :   if (B && typ(B) != t_INT) B = gceil(B);
     581             : 
     582          21 :   if (abscmpiu(X, X_SMALL) <= 0)
     583           7 :     return gerepileupto(av, do_exhaustive(P0, N, itos(X), B));
     584             : 
     585          14 :   if (B && equalii(B,N)) B = NULL;
     586          14 :   if (B) bnd = 1; /* bnd-hack is only for the case B = N */
     587          14 :   P = leafcopy(P0);
     588          14 :   if (!gequal1(gel(P,d+2)))
     589             :   {
     590             :     GEN r, z;
     591           7 :     gel(P,d+2) = bezout(gel(P,d+2), N, &z, &r);
     592           7 :     for (j = 0; j < d; j++) gel(P,j+2) = modii(mulii(gel(P,j+2), z), N);
     593             :   }
     594          14 :   if (DEBUGLEVEL >= 2) err_printf("Modified P: %Ps\n", P);
     595             : 
     596          14 :   choose_params(P, N, X, B, &delta, &t);
     597          14 :   if (DEBUGLEVEL >= 2)
     598           0 :     err_printf("Init: trying delta = %d, t = %d\n", delta, t);
     599             :   for(;;)
     600             :   {
     601          14 :     dim = d * delta + t;
     602             : 
     603             :     /* TODO: In case of failure do not recompute the full vector */
     604          14 :     Xpowers = (GEN*)new_chunk(dim + 1);
     605          14 :     Xpowers[0] = gen_1;
     606          14 :     for (j = 1; j <= dim; j++) Xpowers[j] = mulii(Xpowers[j-1], X);
     607             : 
     608             :     /* TODO: in case of failure, use the part of the matrix already computed */
     609          14 :     M = zeromatcopy(dim,dim);
     610             : 
     611             :     /* Rows of M correspond to the polynomials
     612             :      * N^delta, N^delta Xi, ... N^delta (Xi)^d-1,
     613             :      * N^(delta-1)P(Xi), N^(delta-1)XiP(Xi), ... N^(delta-1)P(Xi)(Xi)^d-1,
     614             :      * ...
     615             :      * P(Xi)^delta, XiP(Xi)^delta, ..., P(Xi)^delta(Xi)^t-1 */
     616          14 :     for (j = 1; j <= d;   j++) gcoeff(M, j, j) = gel(Xpowers,j-1);
     617             : 
     618             :     /* P-part */
     619          14 :     if (delta) row = d + 1; else row = 0;
     620             : 
     621          14 :     Q = P;
     622          70 :     for (i = 1; i < delta; i++)
     623             :     {
     624         182 :       for (j = 0; j < d; j++,row++)
     625        1239 :         for (l = j + 1; l <= row; l++)
     626        1113 :           gcoeff(M, l, row) = mulii(Xpowers[l-1], gel(Q,l-j+1));
     627          56 :       Q = ZX_mul(Q, P);
     628             :     }
     629          63 :     for (j = 0; j < t; row++, j++)
     630         490 :       for (l = j + 1; l <= row; l++)
     631         441 :         gcoeff(M, l, row) = mulii(Xpowers[l-1], gel(Q,l-j+1));
     632             : 
     633             :     /* N-part */
     634          14 :     row = dim - t; N0 = N;
     635          98 :     while (row >= 1)
     636             :     {
     637         224 :       for (j = 0; j < d; j++,row--)
     638        1421 :         for (l = 1; l <= row; l++)
     639        1267 :           gcoeff(M, l, row) = mulii(gmael(M, row, l), N0);
     640          70 :       if (row >= 1) N0 = mulii(N0, N);
     641             :     }
     642             :     /* Z is the upper bound for the L^1 norm of the polynomial,
     643             :        ie. N^delta if B = N, B^delta otherwise */
     644          14 :     if (B) Z = powiu(B, delta); else Z = N0;
     645             : 
     646          14 :     if (DEBUGLEVEL >= 2)
     647             :     {
     648           0 :       if (DEBUGLEVEL >= 6) err_printf("Matrix to be reduced:\n%Ps\n", M);
     649           0 :       err_printf("Entering LLL\nbitsize bound: %ld\n", expi(Z));
     650           0 :       err_printf("expected shvector bitsize: %ld\n", expi(ZM_det_triangular(M))/dim);
     651             :     }
     652             : 
     653          14 :     sh = ZM_lll(M, 0.75, LLL_INPLACE);
     654             :     /* Take the first vector if it is non constant */
     655          14 :     short_pol = gel(sh,1);
     656          14 :     if (ZV_isscalar(short_pol)) short_pol = gel(sh, 2);
     657             : 
     658          14 :     nsp = gen_0;
     659          14 :     for (j = 1; j <= dim; j++) nsp = addii(nsp, absi(gel(short_pol,j)));
     660             : 
     661          14 :     if (DEBUGLEVEL >= 2)
     662             :     {
     663           0 :       err_printf("Candidate: %Ps\n", short_pol);
     664           0 :       err_printf("bitsize Norm: %ld\n", expi(nsp));
     665           0 :       err_printf("bitsize bound: %ld\n", expi(mului(bnd, Z)));
     666             :     }
     667          14 :     if (cmpii(nsp, mului(bnd, Z)) < 0) break; /* SUCCESS */
     668             : 
     669             :     /* Failed with the precomputed or supplied value */
     670           0 :     t++; if (t == d) { delta++; t = 1; }
     671           0 :     if (DEBUGLEVEL >= 2)
     672           0 :       err_printf("Increasing dim, delta = %d t = %d\n", delta, t);
     673           0 :   }
     674          14 :   bnd = itos(divii(nsp, Z)) + 1;
     675             : 
     676          14 :   while (!signe(gel(short_pol,dim))) dim--;
     677             : 
     678          14 :   R = cgetg(dim + 2, t_POL); R[1] = P[1];
     679         217 :   for (j = 1; j <= dim; j++)
     680         203 :     gel(R,j+1) = diviiexact(gel(short_pol,j), Xpowers[j-1]);
     681          14 :   gel(R,2) = subii(gel(R,2), mului(bnd - 1, N0));
     682             : 
     683          14 :   sol = cgetg(1, t_VEC);
     684         112 :   for (i = -bnd + 1; i < bnd; i++)
     685             :   {
     686          98 :     GEN r = nfrootsQ(R);
     687          98 :     if (DEBUGLEVEL >= 2) err_printf("Roots: %Ps\n", r);
     688         119 :     for (j = 1; j < lg(r); j++)
     689             :     {
     690          21 :       GEN z = gel(r,j);
     691          21 :       if (typ(z) == t_INT && sol_OK(FpX_eval(P,z,N), N,B))
     692          14 :         sol = shallowconcat(sol, z);
     693             :     }
     694          98 :     if (i < bnd) gel(R,2) = addii(gel(R,2), Z);
     695             :   }
     696          14 :   return gerepileupto(av, ZV_sort_uniq(sol));
     697             : }
     698             : 
     699             : /********************************************************************/
     700             : /**                                                                **/
     701             : /**                   LINEAR & ALGEBRAIC DEPENDENCE                **/
     702             : /**                                                                **/
     703             : /********************************************************************/
     704             : 
     705             : static int
     706         735 : real_indep(GEN re, GEN im, long bit)
     707             : {
     708         735 :   GEN d = gsub(gmul(gel(re,1),gel(im,2)), gmul(gel(re,2),gel(im,1)));
     709         735 :   return (!gequal0(d) && gexpo(d) > - bit);
     710             : }
     711             : 
     712             : GEN
     713        2919 : lindep_bit(GEN x, long bit)
     714             : {
     715        2919 :   long lx = lg(x), ly, i, j;
     716        2919 :   pari_sp av = avma;
     717             :   GEN re, im, M;
     718             : 
     719        2919 :   if (! is_vec_t(typ(x))) pari_err_TYPE("lindep2",x);
     720        2919 :   if (lx <= 2)
     721             :   {
     722          21 :     if (lx == 2 && gequal0(x)) return mkcol(gen_1);
     723          14 :     return cgetg(1,t_COL);
     724             :   }
     725        2898 :   re = real_i(x);
     726        2898 :   im = imag_i(x);
     727             :   /* independent over R ? */
     728        2898 :   if (lx == 3 && real_indep(re,im,bit)) { avma = av; return cgetg(1, t_COL); }
     729        2891 :   if (gequal0(im)) im = NULL;
     730        2891 :   ly = im? lx+2: lx+1;
     731        2891 :   M = cgetg(lx,t_MAT);
     732       11613 :   for (i=1; i<lx; i++)
     733             :   {
     734        8722 :     GEN c = cgetg(ly,t_COL); gel(M,i) = c;
     735        8722 :     for (j=1; j<lx; j++) gel(c,j) = (i==j)? gen_1: gen_0;
     736        8722 :     gel(c,lx)           = gtrunc2n(gel(re,i), bit);
     737        8722 :     if (im) gel(c,lx+1) = gtrunc2n(gel(im,i), bit);
     738             :   }
     739        2891 :   M = ZM_lll(M, 0.99, LLL_INPLACE);
     740        2891 :   M = gel(M,1);
     741        2891 :   M[0] = evaltyp(t_COL) | evallg(lx);
     742        2891 :   return gerepilecopy(av, M);
     743             : }
     744             : /* deprecated */
     745             : GEN
     746         105 : lindep2(GEN x, long dig)
     747             : {
     748             :   long bit;
     749         105 :   if (dig < 0) pari_err_DOMAIN("lindep2", "accuracy", "<", gen_0, stoi(dig));
     750         105 :   if (dig) bit = (long) (dig/LOG10_2);
     751             :   else
     752             :   {
     753          91 :     bit = gprecision(x);
     754          91 :     if (!bit)
     755             :     {
     756          35 :       x = Q_primpart(x); /* left on stack */
     757          35 :       bit = 32 + gexpo(x);
     758             :     }
     759             :     else
     760          56 :       bit = (long)prec2nbits_mul(bit, 0.8);
     761             :   }
     762         105 :   return lindep_bit(x, bit);
     763             : }
     764             : 
     765             : /* x is a vector of elts of a p-adic field */
     766             : GEN
     767          14 : lindep_padic(GEN x)
     768             : {
     769          14 :   long i, j, prec = LONG_MAX, nx = lg(x)-1, v;
     770          14 :   pari_sp av = avma;
     771          14 :   GEN p = NULL, pn, m, a;
     772             : 
     773          14 :   if (nx < 2) return cgetg(1,t_COL);
     774          49 :   for (i=1; i<=nx; i++)
     775             :   {
     776          35 :     GEN c = gel(x,i), q;
     777          35 :     if (typ(c) != t_PADIC) continue;
     778             : 
     779          21 :     j = precp(c); if (j < prec) prec = j;
     780          21 :     q = gel(c,2);
     781          21 :     if (!p) p = q; else if (!equalii(p, q)) pari_err_MODULUS("lindep_padic", p, q);
     782             :   }
     783          14 :   if (!p) pari_err_TYPE("lindep_padic [not a p-adic vector]",x);
     784          14 :   v = gvaluation(x,p); pn = powiu(p,prec);
     785          14 :   if (v) x = gmul(x, powis(p, -v));
     786          14 :   x = RgV_to_FpV(x, pn);
     787             : 
     788          14 :   a = negi(gel(x,1));
     789          14 :   m = cgetg(nx,t_MAT);
     790          35 :   for (i=1; i<nx; i++)
     791             :   {
     792          21 :     GEN c = zerocol(nx);
     793          21 :     gel(c,1+i) = a;
     794          21 :     gel(c,1) = gel(x,i+1);
     795          21 :     gel(m,i) = c;
     796             :   }
     797          14 :   m = ZM_lll(ZM_hnfmodid(m, pn), 0.99, LLL_INPLACE);
     798          14 :   return gerepilecopy(av, gel(m,1));
     799             : }
     800             : /* x is a vector of t_POL/t_SER */
     801             : GEN
     802          42 : lindep_Xadic(GEN x)
     803             : {
     804          42 :   long i, prec = LONG_MAX, deg = 0, lx = lg(x), vx, v;
     805          42 :   pari_sp av = avma;
     806             :   GEN m;
     807             : 
     808          42 :   if (lx == 1) return cgetg(1,t_COL);
     809          42 :   vx = gvar(x);
     810          42 :   v = gvaluation(x, pol_x(vx));
     811          42 :   if (!v)         x = shallowcopy(x);
     812           0 :   else if (v > 0) x = gdiv(x, pol_xn(v, vx));
     813           0 :   else            x = gmul(x, pol_xn(-v, vx));
     814             :   /* all t_SER have valuation >= 0 */
     815         308 :   for (i=1; i<lx; i++)
     816             :   {
     817         266 :     GEN c = gel(x,i);
     818         266 :     if (gvar(c) != vx) { gel(x,i) = scalarpol_shallow(c, vx); continue; }
     819         259 :     switch(typ(c))
     820             :     {
     821         126 :       case t_POL: deg = maxss(deg, degpol(c)); break;
     822           0 :       case t_RFRAC: pari_err_TYPE("lindep_Xadic", c);
     823             :       case t_SER:
     824         133 :         prec = minss(prec, valp(c)+lg(c)-2);
     825         133 :         gel(x,i) = ser2rfrac_i(c);
     826             :     }
     827             :   }
     828          42 :   if (prec == LONG_MAX) prec = deg+1;
     829          42 :   m = RgXV_to_RgM(x, prec);
     830          42 :   return gerepileupto(av, deplin(m));
     831             : }
     832             : static GEN
     833          35 : vec_lindep(GEN x)
     834             : {
     835          35 :   pari_sp av = avma;
     836          35 :   long i, l = lg(x); /* > 1 */
     837          35 :   long t = typ(gel(x,1)), h = lg(gel(x,1));
     838          35 :   GEN m = cgetg(l, t_MAT);
     839         126 :   for (i = 1; i < l; i++)
     840             :   {
     841          98 :     GEN y = gel(x,i);
     842          98 :     if (lg(y) != h || typ(y) != t) pari_err_TYPE("lindep",x);
     843          91 :     if (t != t_COL) y = shallowtrans(y); /* Sigh */
     844          91 :     gel(m,i) = y;
     845             :   }
     846          28 :   return gerepileupto(av, deplin(m));
     847             : }
     848             : 
     849             : GEN
     850           0 : lindep(GEN x) { return lindep2(x, 0); }
     851             : 
     852             : GEN
     853         336 : lindep0(GEN x,long bit)
     854             : {
     855         336 :   long i, tx = typ(x);
     856         336 :   if (tx == t_MAT) return deplin(x);
     857         133 :   if (! is_vec_t(tx)) pari_err_TYPE("lindep",x);
     858         399 :   for (i = 1; i < lg(x); i++)
     859         322 :     switch(typ(gel(x,i)))
     860             :     {
     861           7 :       case t_PADIC: return lindep_padic(x);
     862             :       case t_POL:
     863             :       case t_RFRAC:
     864          14 :       case t_SER: return lindep_Xadic(x);
     865             :       case t_VEC:
     866          35 :       case t_COL: return vec_lindep(x);
     867             :     }
     868          77 :   return lindep2(x, bit);
     869             : }
     870             : 
     871             : GEN
     872          49 : algdep0(GEN x, long n, long bit)
     873             : {
     874          49 :   long tx = typ(x), i;
     875             :   pari_sp av;
     876             :   GEN y;
     877             : 
     878          49 :   if (! is_scalar_t(tx)) pari_err_TYPE("algdep0",x);
     879          49 :   if (tx==t_POLMOD) { y = RgX_copy(gel(x,1)); setvarn(y,0); return y; }
     880          49 :   if (gequal0(x)) return pol_x(0);
     881          49 :   if (n <= 0)
     882             :   {
     883          14 :     if (!n) return gen_1;
     884           7 :     pari_err_DOMAIN("algdep", "degree", "<", gen_0, stoi(n));
     885             :   }
     886             : 
     887          35 :   av = avma; y = cgetg(n+2,t_COL);
     888          35 :   gel(y,1) = gen_1;
     889          35 :   gel(y,2) = x; /* n >= 1 */
     890          35 :   for (i=3; i<=n+1; i++) gel(y,i) = gmul(gel(y,i-1),x);
     891          35 :   if (typ(x) == t_PADIC)
     892           7 :     y = lindep_padic(y);
     893             :   else
     894          28 :     y = lindep2(y, bit);
     895          35 :   if (lg(y) == 1) pari_err(e_DOMAIN,"algdep", "degree(x)",">", stoi(n), x);
     896          35 :   y = RgV_to_RgX(y, 0);
     897          35 :   if (signe(leading_coeff(y)) > 0) return gerepilecopy(av, y);
     898           0 :   return gerepileupto(av, ZX_neg(y));
     899             : }
     900             : 
     901             : GEN
     902           0 : algdep(GEN x, long n)
     903             : {
     904           0 :   return algdep0(x,n,0);
     905             : }
     906             : 
     907             : GEN
     908          28 : seralgdep(GEN s, long p, long r)
     909             : {
     910          28 :   pari_sp av = avma;
     911             :   long vy, i, m, n, prec;
     912             :   GEN S, v, D;
     913             : 
     914          28 :   if (typ(s) != t_SER) pari_err_TYPE("seralgdep",s);
     915          28 :   if (p <= 0) pari_err_DOMAIN("seralgdep", "p", "<=", gen_0, stoi(p));
     916          28 :   if (r < 0) pari_err_DOMAIN("seralgdep", "r", "<", gen_0, stoi(r));
     917          28 :   if (is_bigint(addiu(muluu(p, r), 1))) pari_err_OVERFLOW("seralgdep");
     918          28 :   vy = varn(s);
     919          28 :   if (!vy) pari_err_PRIORITY("seralgdep", s, ">", 0);
     920          28 :   r++; p++;
     921          28 :   prec = valp(s) + lg(s)-2;
     922          28 :   if (r > prec) r = prec;
     923          28 :   S = cgetg(p+1, t_VEC); gel(S, 1) = s;
     924          28 :   for (i = 2; i <= p; i++) gel(S,i) = gmul(gel(S,i-1), s);
     925          28 :   v = cgetg(r*p+1, t_VEC); /* v[r*n+m+1] = s^n * y^m */
     926             :   /* n = 0 */
     927          28 :   for (m = 0; m < r; m++) gel(v, m+1) = pol_xn(m, vy);
     928          70 :   for(n=1; n < p; n++)
     929         175 :     for (m = 0; m < r; m++)
     930             :     {
     931         133 :       GEN c = gel(S,n);
     932         133 :       if (m)
     933             :       {
     934          91 :         c = shallowcopy(c);
     935          91 :         setvalp(c, valp(c) + m);
     936             :       }
     937         133 :       gel(v, r*n + m + 1) = c;
     938             :     }
     939          28 :   D = lindep_Xadic(v);
     940          28 :   if (lg(D) == 1) { avma = av; return gen_0; }
     941          21 :   v = cgetg(p+1, t_VEC);
     942          70 :   for (n = 0; n < p; n++)
     943          49 :     gel(v, n+1) = RgV_to_RgX(vecslice(D, r*n+1, r*n+r), vy);
     944          21 :   return gerepilecopy(av, RgV_to_RgX(v, 0));
     945             : }
     946             : 
     947             : /********************************************************************/
     948             : /**                                                                **/
     949             : /**                              MINIM                             **/
     950             : /**                                                                **/
     951             : /********************************************************************/
     952             : void
     953       27776 : minim_alloc(long n, double ***q, GEN *x, double **y,  double **z, double **v)
     954             : {
     955             :   long i, s;
     956             : 
     957       27776 :   *x = cgetg(n, t_VECSMALL);
     958       27776 :   *q = (double**) new_chunk(n);
     959       27776 :   s = n * sizeof(double);
     960       27776 :   *y = (double*) stack_malloc_align(s, sizeof(double));
     961       27776 :   *z = (double*) stack_malloc_align(s, sizeof(double));
     962       27776 :   *v = (double*) stack_malloc_align(s, sizeof(double));
     963       27776 :   for (i=1; i<n; i++) (*q)[i] = (double*) stack_malloc_align(s, sizeof(double));
     964       27776 : }
     965             : 
     966             : static GEN
     967      202657 : ZC_canon(GEN V)
     968             : {
     969      202657 :   long l = lg(V), j;
     970      202657 :   for (j = 1; j < l  &&  signe(gel(V,j)) == 0; ++j);
     971      202657 :   return (j < l  &&  signe(gel(V,j)) < 0)? ZC_neg(V): V;
     972             : }
     973             : 
     974             : static GEN
     975      202657 : ZM_zc_mul_canon(GEN u, GEN x)
     976             : {
     977      202657 :   return ZC_canon(ZM_zc_mul(u,x));
     978             : }
     979             : 
     980             : struct qfvec
     981             : {
     982             :   GEN a, r, u;
     983             : };
     984             : 
     985             : static void
     986           0 : err_minim(GEN a)
     987             : {
     988           0 :   pari_err_DOMAIN("minim0","form","is not",
     989             :                   strtoGENstr("positive definite"),a);
     990           0 : }
     991             : 
     992             : static GEN
     993         580 : minim_lll(GEN a, GEN *u)
     994             : {
     995         580 :   *u = lllgramint(a);
     996         580 :   if (lg(*u) != lg(a)) err_minim(a);
     997         580 :   return qf_apply_ZM(a,*u);
     998             : }
     999             : 
    1000             : static void
    1001         580 : forqfvec_init_dolll(struct qfvec *qv, GEN a, long dolll)
    1002             : {
    1003             :   GEN r, u;
    1004         580 :   if (!dolll) u = NULL;
    1005             :   else
    1006             :   {
    1007         538 :     if (typ(a) != t_MAT || !RgM_is_ZM(a)) pari_err_TYPE("qfminim",a);
    1008         538 :     a = minim_lll(a, &u);
    1009             :   }
    1010         580 :   qv->a = RgM_gtofp(a, DEFAULTPREC);
    1011         580 :   r = qfgaussred_positive(qv->a);
    1012         580 :   if (!r)
    1013             :   {
    1014           0 :     r = qfgaussred_positive(a); /* exact computation */
    1015           0 :     if (!r) err_minim(a);
    1016           0 :     r = RgM_gtofp(r, DEFAULTPREC);
    1017             :   }
    1018         580 :   qv->r = r;
    1019         580 :   qv->u = u;
    1020         580 : }
    1021             : 
    1022             : static void
    1023          21 : forqfvec_init(struct qfvec *qv, GEN a)
    1024          21 : { forqfvec_init_dolll(qv, a, 1); }
    1025             : 
    1026             : static void
    1027          21 : forqfvec_i(void *E, long (*fun)(void *, GEN, GEN, double), struct qfvec *qv, GEN BORNE)
    1028             : {
    1029          21 :   GEN x, a = qv->a, r = qv->r, u = qv->u;
    1030          21 :   long n = lg(a), i, j, k;
    1031             :   double p,BOUND,*v,*y,*z,**q;
    1032          21 :   const double eps = 0.0001;
    1033          21 :   if (!BORNE) BORNE = gen_0;
    1034             :   else
    1035             :   {
    1036          14 :     BORNE = gfloor(BORNE);
    1037          14 :     if (typ(BORNE) != t_INT) pari_err_TYPE("minim0",BORNE);
    1038             :   }
    1039          42 :   if (n == 1) return;
    1040          14 :   minim_alloc(n, &q, &x, &y, &z, &v);
    1041          14 :   n--;
    1042          42 :   for (j=1; j<=n; j++)
    1043             :   {
    1044          28 :     v[j] = rtodbl(gcoeff(r,j,j));
    1045          28 :     for (i=1; i<j; i++) q[i][j] = rtodbl(gcoeff(r,i,j));
    1046             :   }
    1047             : 
    1048          14 :   if (gequal0(BORNE))
    1049             :   {
    1050             :     double c;
    1051           7 :     p = rtodbl(gcoeff(a,1,1));
    1052           7 :     for (i=2; i<=n; i++) { c = rtodbl(gcoeff(a,i,i)); if (c < p) p = c; }
    1053           7 :     BORNE = roundr(dbltor(p));
    1054             :   }
    1055             :   else
    1056           7 :     p = gtodouble(BORNE);
    1057          14 :   BOUND = p * (1 + eps);
    1058          14 :   if (BOUND == p) pari_err_PREC("minim0");
    1059             : 
    1060          14 :   k = n; y[n] = z[n] = 0;
    1061          14 :   x[n] = (long)sqrt(BOUND/v[n]);
    1062          28 :   for(;;x[1]--)
    1063             :   {
    1064             :     do
    1065             :     {
    1066          49 :       if (k>1)
    1067             :       {
    1068          21 :         long l = k-1;
    1069          21 :         z[l] = 0;
    1070          21 :         for (j=k; j<=n; j++) z[l] += q[l][j]*x[j];
    1071          21 :         p = (double)x[k] + z[k];
    1072          21 :         y[l] = y[k] + p*p*v[k];
    1073          21 :         x[l] = (long)floor(sqrt((BOUND-y[l])/v[l])-z[l]);
    1074          21 :         k = l;
    1075             :       }
    1076             :       for(;;)
    1077             :       {
    1078          56 :         p = (double)x[k] + z[k];
    1079          56 :         if (y[k] + p*p*v[k] <= BOUND) break;
    1080           7 :         k++; x[k]--;
    1081           7 :       }
    1082          49 :     } while (k > 1);
    1083          42 :     if (! x[1] && y[1]<=eps) break;
    1084             : 
    1085          28 :     p = (double)x[1] + z[1]; p = y[1] + p*p*v[1]; /* norm(x) */
    1086          28 :     if (fun(E, u, x, p)) break;
    1087          28 :   }
    1088             : }
    1089             : 
    1090             : void
    1091           0 : forqfvec(void *E, long (*fun)(void *, GEN, GEN, double), GEN a, GEN BORNE)
    1092             : {
    1093           0 :   pari_sp av = avma;
    1094             :   struct qfvec qv;
    1095           0 :   forqfvec_init(&qv, a);
    1096           0 :   forqfvec_i(E, fun, &qv, BORNE);
    1097           0 :   avma = av;
    1098           0 : }
    1099             : 
    1100             : static long
    1101          28 : _gp_forqf(void *E, GEN u, GEN x, double p/*unused*/)
    1102             : {
    1103          28 :   pari_sp av = avma;
    1104             :   (void)p;
    1105          28 :   set_lex(-1, ZM_zc_mul_canon(u, x));
    1106          28 :   closure_evalvoid((GEN)E);
    1107          28 :   avma = av;
    1108          28 :   return loop_break();
    1109             : }
    1110             : 
    1111             : void
    1112          21 : forqfvec0(GEN a, GEN BORNE, GEN code)
    1113             : {
    1114          21 :   pari_sp av = avma;
    1115             :   struct qfvec qv;
    1116          21 :   forqfvec_init(&qv, a);
    1117          21 :   push_lex(gen_0, code);
    1118          21 :   forqfvec_i((void*) code, &_gp_forqf, &qv, BORNE);
    1119          21 :   pop_lex(1);
    1120          21 :   avma = av;
    1121          21 : }
    1122             : 
    1123             : enum { min_ALL = 0, min_FIRST, min_VECSMALL, min_VECSMALL2 };
    1124             : 
    1125             : /* Minimal vectors for the integral definite quadratic form: a.
    1126             :  * Result u:
    1127             :  *   u[1]= Number of vectors of square norm <= BORNE
    1128             :  *   u[2]= maximum norm found
    1129             :  *   u[3]= list of vectors found (at most STOCKMAX, unless NULL)
    1130             :  *
    1131             :  *  If BORNE = NULL: Minimal non-zero vectors.
    1132             :  *  flag = min_ALL,   as above
    1133             :  *  flag = min_FIRST, exits when first suitable vector is found.
    1134             :  *  flag = min_VECSMALL, return a t_VECSMALL of (half) the number of vectors
    1135             :  *  for each norm
    1136             :  *  flag = min_VECSMALL2, same but count only vectors with even norm, and shift
    1137             :  *  the answer */
    1138             : static GEN
    1139         588 : minim0_dolll(GEN a, GEN BORNE, GEN STOCKMAX, long flag, long dolll)
    1140             : {
    1141             :   GEN x, u, r, L, gnorme;
    1142         588 :   long n = lg(a), i, j, k, s, maxrank, sBORNE;
    1143         588 :   pari_sp av = avma, av1;
    1144             :   double p,maxnorm,BOUND,*v,*y,*z,**q;
    1145         588 :   const double eps = 1e-10;
    1146         588 :   int stockall = 0;
    1147             :   struct qfvec qv;
    1148             : 
    1149         588 :   if (!BORNE)
    1150          56 :     sBORNE = 0;
    1151             :   else
    1152             :   {
    1153         532 :     BORNE = gfloor(BORNE);
    1154         532 :     if (typ(BORNE) != t_INT) pari_err_TYPE("minim0",BORNE);
    1155         532 :     if (is_bigint(BORNE)) pari_err_PREC( "qfminim");
    1156         531 :     sBORNE = itos(BORNE); avma = av;
    1157             :   }
    1158         587 :   if (!STOCKMAX)
    1159             :   {
    1160         244 :     stockall = 1;
    1161         244 :     maxrank = 200;
    1162             :   }
    1163             :   else
    1164             :   {
    1165         343 :     STOCKMAX = gfloor(STOCKMAX);
    1166         343 :     if (typ(STOCKMAX) != t_INT) pari_err_TYPE("minim0",STOCKMAX);
    1167         343 :     maxrank = itos(STOCKMAX);
    1168         343 :     if (maxrank < 0)
    1169           0 :       pari_err_TYPE("minim0 [negative number of vectors]",STOCKMAX);
    1170             :   }
    1171             : 
    1172         587 :   switch(flag)
    1173             :   {
    1174             :     case min_VECSMALL:
    1175             :     case min_VECSMALL2:
    1176         294 :       if (sBORNE <= 0) return cgetg(1, t_VECSMALL);
    1177         294 :       L = zero_zv(sBORNE);
    1178         294 :       if (flag == min_VECSMALL2) sBORNE <<= 1;
    1179         294 :       if (n == 1) return L;
    1180         287 :       break;
    1181             :     case min_FIRST:
    1182          35 :       if (n == 1 || (!sBORNE && BORNE)) return cgetg(1,t_VEC);
    1183          21 :       L = NULL; /* gcc -Wall */
    1184          21 :       break;
    1185             :     case min_ALL:
    1186         258 :       if (n == 1 || (!sBORNE && BORNE))
    1187           7 :         retmkvec3(gen_0, gen_0, cgetg(1, t_MAT));
    1188         251 :       L = new_chunk(1+maxrank);
    1189         251 :       break;
    1190             :     default:
    1191           0 :       return NULL;
    1192             :   }
    1193         559 :   minim_alloc(n, &q, &x, &y, &z, &v);
    1194             : 
    1195         559 :   forqfvec_init_dolll(&qv, a, dolll);
    1196         559 :   av1 = avma;
    1197         559 :   r = qv.r;
    1198         559 :   u = qv.u;
    1199         559 :   n--;
    1200        4834 :   for (j=1; j<=n; j++)
    1201             :   {
    1202        4275 :     v[j] = rtodbl(gcoeff(r,j,j));
    1203        4275 :     for (i=1; i<j; i++) q[i][j] = rtodbl(gcoeff(r,i,j)); /* |.| <= 1/2 */
    1204             :   }
    1205             : 
    1206         559 :   if (sBORNE) maxnorm = 0.;
    1207             :   else
    1208             :   {
    1209          56 :     GEN B = gcoeff(a,1,1);
    1210          56 :     long t = 1;
    1211         616 :     for (i=2; i<=n; i++)
    1212             :     {
    1213         560 :       GEN c = gcoeff(a,i,i);
    1214         560 :       if (cmpii(c, B) < 0) { B = c; t = i; }
    1215             :     }
    1216          56 :     if (flag == min_FIRST) return gerepilecopy(av, mkvec2(B, gel(u,t)));
    1217          49 :     maxnorm = -1.; /* don't update maxnorm */
    1218          49 :     if (is_bigint(B)) return NULL;
    1219          48 :     sBORNE = itos(B);
    1220             :   }
    1221         551 :   BOUND = sBORNE * (1 + eps);
    1222         551 :   if ((long)BOUND != sBORNE) return NULL;
    1223             : 
    1224         539 :   s = 0;
    1225         539 :   k = n; y[n] = z[n] = 0;
    1226         539 :   x[n] = (long)sqrt(BOUND/v[n]);
    1227     1034026 :   for(;;x[1]--)
    1228             :   {
    1229             :     do
    1230             :     {
    1231     1922550 :       if (k>1)
    1232             :       {
    1233      888517 :         long l = k-1;
    1234      888517 :         z[l] = 0;
    1235      888517 :         for (j=k; j<=n; j++) z[l] += q[l][j]*x[j];
    1236      888517 :         p = (double)x[k] + z[k];
    1237      888517 :         y[l] = y[k] + p*p*v[k];
    1238      888517 :         x[l] = (long)floor(sqrt((BOUND-y[l])/v[l])-z[l]);
    1239      888517 :         k = l;
    1240             :       }
    1241             :       for(;;)
    1242             :       {
    1243     2807539 :         p = (double)x[k] + z[k];
    1244     2807539 :         if (y[k] + p*p*v[k] <= BOUND) break;
    1245      884989 :         k++; x[k]--;
    1246      884989 :       }
    1247             :     }
    1248     1922550 :     while (k > 1);
    1249     1034565 :     if (! x[1] && y[1]<=eps) break;
    1250             : 
    1251     1034033 :     p = (double)x[1] + z[1]; p = y[1] + p*p*v[1]; /* norm(x) */
    1252     1034033 :     if (maxnorm >= 0)
    1253             :     {
    1254     1031485 :       if (p > maxnorm) maxnorm = p;
    1255             :     }
    1256             :     else
    1257             :     { /* maxnorm < 0 : only look for minimal vectors */
    1258        2548 :       pari_sp av2 = avma;
    1259        2548 :       gnorme = roundr(dbltor(p));
    1260        2548 :       if (cmpis(gnorme, sBORNE) >= 0) avma = av2;
    1261             :       else
    1262             :       {
    1263          14 :         sBORNE = itos(gnorme); avma = av1;
    1264          14 :         BOUND = sBORNE * (1+eps);
    1265          14 :         L = new_chunk(maxrank+1);
    1266          14 :         s = 0;
    1267             :       }
    1268             :     }
    1269     1034033 :     s++;
    1270             : 
    1271     1034033 :     switch(flag)
    1272             :     {
    1273             :       case min_FIRST:
    1274           7 :         if (dolll) x = ZM_zc_mul_canon(u,x);
    1275           7 :         return gerepilecopy(av, mkvec2(roundr(dbltor(p)), x));
    1276             : 
    1277             :       case min_ALL:
    1278      205058 :         if (s > maxrank && stockall) /* overflow */
    1279             :         {
    1280         406 :           long maxranknew = maxrank << 1;
    1281         406 :           GEN Lnew = new_chunk(1 + maxranknew);
    1282         406 :           for (i=1; i<=maxrank; i++) Lnew[i] = L[i];
    1283         406 :           L = Lnew; maxrank = maxranknew;
    1284             :         }
    1285      205058 :         if (s<=maxrank) gel(L,s) = leafcopy(x);
    1286      205058 :         break;
    1287             : 
    1288             :       case min_VECSMALL:
    1289       39200 :         { ulong norm = (ulong)(p + 0.5); L[norm]++; }
    1290       39200 :         break;
    1291             : 
    1292             :       case min_VECSMALL2:
    1293      789768 :         { ulong norm = (ulong)(p + 0.5); if (!odd(norm)) L[norm>>1]++; }
    1294      789768 :         break;
    1295             : 
    1296             :     }
    1297     1034026 :   }
    1298         532 :   switch(flag)
    1299             :   {
    1300             :     case min_FIRST:
    1301           7 :       avma = av; return cgetg(1,t_VEC);
    1302             :     case min_VECSMALL:
    1303             :     case min_VECSMALL2:
    1304         287 :       avma = (pari_sp)L; return L;
    1305             :   }
    1306         238 :   r = (maxnorm >= 0) ? roundr(dbltor(maxnorm)): stoi(sBORNE);
    1307         238 :   k = minss(s,maxrank);
    1308         238 :   L[0] = evaltyp(t_MAT) | evallg(k + 1);
    1309         238 :   if (dolll)
    1310         203 :     for (j=1; j<=k; j++) gel(L,j) = ZM_zc_mul_canon(u, gel(L,j));
    1311         238 :   return gerepilecopy(av, mkvec3(stoi(s<<1), r, L));
    1312             : }
    1313             : 
    1314             : static GEN
    1315         546 : minim0(GEN a, GEN BORNE, GEN STOCKMAX, long flag)
    1316             : {
    1317         546 :   GEN v = minim0_dolll(a, BORNE, STOCKMAX, flag, 1);
    1318         545 :   if (!v) pari_err_PREC("qfminim");
    1319         539 :   return v;
    1320             : }
    1321             : 
    1322             : GEN
    1323         294 : qfrep0(GEN a, GEN borne, long flag)
    1324         294 : { return minim0(a, borne, gen_0, (flag & 1)? min_VECSMALL2: min_VECSMALL); }
    1325             : 
    1326             : GEN
    1327         105 : qfminim0(GEN a, GEN borne, GEN stockmax, long flag, long prec)
    1328             : {
    1329         105 :   switch(flag)
    1330             :   {
    1331          42 :     case 0: return minim0(a,borne,stockmax,min_ALL);
    1332          35 :     case 1: return minim0(a,borne,gen_0   ,min_FIRST);
    1333             :     case 2:
    1334             :     {
    1335          28 :       long maxnum = -1;
    1336          28 :       if (typ(a) != t_MAT) pari_err_TYPE("qfminim",a);
    1337          28 :       if (stockmax) {
    1338          14 :         if (typ(stockmax) != t_INT) pari_err_TYPE("qfminim",stockmax);
    1339          14 :         maxnum = itos(stockmax);
    1340             :       }
    1341          28 :       a = fincke_pohst(a,borne,maxnum,prec,NULL);
    1342          28 :       if (!a) pari_err_PREC("qfminim");
    1343          28 :       return a;
    1344             :     }
    1345           0 :     default: pari_err_FLAG("qfminim");
    1346             :   }
    1347             :   return NULL; /* LCOV_EXCL_LINE */
    1348             : }
    1349             : 
    1350             : GEN
    1351         175 : minim(GEN a, GEN borne, GEN stockmax)
    1352         175 : { return minim0(a,borne,stockmax,min_ALL); }
    1353             : 
    1354             : GEN
    1355          42 : minim_raw(GEN a, GEN BORNE, GEN STOCKMAX)
    1356          42 : { return minim0_dolll(a, BORNE, STOCKMAX, min_ALL, 0); }
    1357             : 
    1358             : GEN
    1359           0 : minim2(GEN a, GEN borne, GEN stockmax)
    1360           0 : { return minim0(a,borne,stockmax,min_FIRST); }
    1361             : 
    1362             : /* If V depends linearly from the columns of the matrix, return 0.
    1363             :  * Otherwise, update INVP and L and return 1. No GC. */
    1364             : static int
    1365        1652 : addcolumntomatrix(GEN V, GEN invp, GEN L)
    1366             : {
    1367        1652 :   long i,j,k, n = lg(invp);
    1368        1652 :   GEN a = cgetg(n, t_COL), ak = NULL, mak;
    1369             : 
    1370       84231 :   for (k = 1; k < n; k++)
    1371       83706 :     if (!L[k])
    1372             :     {
    1373       27811 :       ak = RgMrow_zc_mul(invp, V, k);
    1374       27811 :       if (!gequal0(ak)) break;
    1375             :     }
    1376        1652 :   if (k == n) return 0;
    1377        1127 :   L[k] = 1;
    1378        1127 :   mak = gneg_i(ak);
    1379       43253 :   for (i=k+1; i<n; i++)
    1380       42126 :     gel(a,i) = gdiv(RgMrow_zc_mul(invp, V, i), mak);
    1381       43883 :   for (j=1; j<=k; j++)
    1382             :   {
    1383       42756 :     GEN c = gel(invp,j), ck = gel(c,k);
    1384       42756 :     if (gequal0(ck)) continue;
    1385        9471 :     gel(c,k) = gdiv(ck, ak);
    1386        9471 :     if (j==k)
    1387       43253 :       for (i=k+1; i<n; i++)
    1388       42126 :         gel(c,i) = gmul(gel(a,i), ck);
    1389             :     else
    1390      209979 :       for (i=k+1; i<n; i++)
    1391      201635 :         gel(c,i) = gadd(gel(c,i), gmul(gel(a,i), ck));
    1392             :   }
    1393        1127 :   return 1;
    1394             : }
    1395             : 
    1396             : GEN
    1397          42 : perf(GEN a)
    1398             : {
    1399          42 :   pari_sp av = avma;
    1400             :   GEN u, L;
    1401          42 :   long r, s, k, l, n = lg(a)-1;
    1402             : 
    1403          42 :   if (!n) return gen_0;
    1404          42 :   if (typ(a) != t_MAT || !RgM_is_ZM(a)) pari_err_TYPE("qfperfection",a);
    1405          42 :   a = minim_lll(a, &u);
    1406          42 :   L = minim_raw(a,NULL,NULL);
    1407          42 :   r = (n*(n+1)) >> 1;
    1408          42 :   if (L)
    1409             :   {
    1410             :     GEN D, V, invp;
    1411          35 :     L = gel(L, 3); l = lg(L);
    1412          35 :     if (l == 2) { avma = av; return gen_1; }
    1413             : 
    1414          21 :     D = zero_zv(r);
    1415          21 :     V = cgetg(r+1, t_VECSMALL);
    1416          21 :     invp = matid(r);
    1417          21 :     s = 0;
    1418        1659 :     for (k = 1; k < l; k++)
    1419             :     {
    1420        1652 :       pari_sp av2 = avma;
    1421        1652 :       GEN x = gel(L,k);
    1422             :       long i, j, I;
    1423       21098 :       for (i = I = 1; i<=n; i++)
    1424       19446 :         for (j=i; j<=n; j++,I++) V[I] = x[i]*x[j];
    1425        1652 :       if (!addcolumntomatrix(V,invp,D)) avma = av2;
    1426        1127 :       else if (++s == r) break;
    1427             :     }
    1428             :   }
    1429             :   else
    1430             :   {
    1431             :     GEN M;
    1432           7 :     L = fincke_pohst(a,NULL,-1, DEFAULTPREC, NULL);
    1433           7 :     if (!L) pari_err_PREC("qfminim");
    1434           7 :     L = gel(L, 3); l = lg(L);
    1435           7 :     if (l == 2) { avma = av; return gen_1; }
    1436           7 :     M = cgetg(l, t_MAT);
    1437         959 :     for (k = 1; k < l; k++)
    1438             :     {
    1439         952 :       GEN x = gel(L,k), c = cgetg(r+1, t_COL);
    1440             :       long i, I, j;
    1441       16184 :       for (i = I = 1; i<=n; i++)
    1442       15232 :         for (j=i; j<=n; j++,I++) gel(c,I) = mulii(gel(x,i), gel(x,j));
    1443         952 :       gel(M,k) = c;
    1444             :     }
    1445           7 :     s = ZM_rank(M);
    1446             :   }
    1447          28 :  avma = av; return utoipos(s);
    1448             : }
    1449             : 
    1450             : static GEN
    1451          64 : clonefill(GEN S, long s, long t)
    1452             : { /* initialize to dummy values */
    1453          64 :   GEN T = S, dummy = cgetg(1, t_STR);
    1454             :   long i;
    1455          64 :   for (i = s+1; i <= t; i++) gel(S,i) = dummy;
    1456          64 :   S = gclone(S); if (isclone(T)) gunclone(T);
    1457          64 :   return S;
    1458             : }
    1459             : 
    1460             : /* increment ZV x, by incrementing cell of index k. Initial value x0[k] was
    1461             :  * chosen to minimize qf(x) for given x0[1..k-1] and x0[k+1,..] = 0
    1462             :  * The last non-zero entry must be positive and goes through x0[k]+1,2,3,...
    1463             :  * Others entries go through: x0[k]+1,-1,2,-2,...*/
    1464             : INLINE void
    1465      916094 : step(GEN x, GEN y, GEN inc, long k)
    1466             : {
    1467      916094 :   if (!signe(gel(y,k))) /* x[k+1..] = 0 */
    1468        9471 :     gel(x,k) = addiu(gel(x,k), 1); /* leading coeff > 0 */
    1469             :   else
    1470             :   {
    1471      906623 :     long i = inc[k];
    1472      906623 :     gel(x,k) = addis(gel(x,k), i),
    1473      906623 :     inc[k] = (i > 0)? -1-i: 1-i;
    1474             :   }
    1475      916094 : }
    1476             : 
    1477             : /* 1 if we are "sure" that x < y, up to few rounding errors, i.e.
    1478             :  * x < y - epsilon. More precisely :
    1479             :  * - sign(x - y) < 0
    1480             :  * - lgprec(x-y) > 3 || expo(x - y) - expo(x) > -24 */
    1481             : static int
    1482      359614 : mplessthan(GEN x, GEN y)
    1483             : {
    1484      359614 :   pari_sp av = avma;
    1485      359614 :   GEN z = mpsub(x, y);
    1486      359614 :   avma = av;
    1487      359614 :   if (typ(z) == t_INT) return (signe(z) < 0);
    1488      359614 :   if (signe(z) >= 0) return 0;
    1489       58704 :   if (realprec(z) > LOWDEFAULTPREC) return 1;
    1490       58704 :   return ( expo(z) - mpexpo(x) > -24 );
    1491             : }
    1492             : 
    1493             : /* 1 if we are "sure" that x > y, up to few rounding errors, i.e.
    1494             :  * x > y + epsilon */
    1495             : static int
    1496     1450147 : mpgreaterthan(GEN x, GEN y)
    1497             : {
    1498     1450147 :   pari_sp av = avma;
    1499     1450147 :   GEN z = mpsub(x, y);
    1500     1450147 :   avma = av;
    1501     1450147 :   if (typ(z) == t_INT) return (signe(z) > 0);
    1502     1450147 :   if (signe(z) <= 0) return 0;
    1503      820635 :   if (realprec(z) > LOWDEFAULTPREC) return 1;
    1504      145158 :   return ( expo(z) - mpexpo(x) > -24 );
    1505             : }
    1506             : 
    1507             : /* x a t_INT, y  t_INT or t_REAL */
    1508             : INLINE GEN
    1509      358158 : mulimp(GEN x, GEN y)
    1510             : {
    1511      358158 :   if (typ(y) == t_INT) return mulii(x,y);
    1512      358158 :   return signe(x) ? mulir(x,y): gen_0;
    1513             : }
    1514             : /* x + y*z, x,z two mp's, y a t_INT */
    1515             : INLINE GEN
    1516     4507165 : addmulimp(GEN x, GEN y, GEN z)
    1517             : {
    1518     4507165 :   if (!signe(y)) return x;
    1519     1948991 :   if (typ(z) == t_INT) return mpadd(x, mulii(y, z));
    1520     1948991 :   return mpadd(x, mulir(y, z));
    1521             : }
    1522             : 
    1523             : /* yk + vk * (xk + zk)^2 */
    1524             : static GEN
    1525     1803279 : norm_aux(GEN xk, GEN yk, GEN zk, GEN vk)
    1526             : {
    1527     1803279 :   GEN t = mpadd(xk, zk);
    1528     1803279 :   if (typ(t) == t_INT) { /* probably gen_0, avoid loss of accuracy */
    1529       17649 :     yk = addmulimp(yk, sqri(t), vk);
    1530             :   } else {
    1531     1785630 :     yk = mpadd(yk, mpmul(sqrr(t), vk));
    1532             :   }
    1533     1803279 :   return yk;
    1534             : }
    1535             : /* yk + vk * (xk + zk)^2 < B + epsilon */
    1536             : static int
    1537     1273825 : check_bound(GEN B, GEN xk, GEN yk, GEN zk, GEN vk)
    1538             : {
    1539     1273825 :   pari_sp av = avma;
    1540     1273825 :   int f = mpgreaterthan(norm_aux(xk,yk,zk,vk), B);
    1541     1273825 :   avma = av; return !f;
    1542             : }
    1543             : 
    1544             : /* q(k-th canonical basis vector), where q is given in Cholesky form
    1545             :  * q(x) = sum_{i = 1}^n q[i,i] (x[i] + sum_{j > i} q[i,j] x[j])^2.
    1546             :  * Namely q(e_k) = q[k,k] + sum_{i < k} q[i,i] q[i,k]^2
    1547             :  * Assume 1 <= k <= n. */
    1548             : static GEN
    1549         182 : cholesky_norm_ek(GEN q, long k)
    1550             : {
    1551         182 :   GEN t = gcoeff(q,k,k);
    1552             :   long i;
    1553         182 :   for (i = 1; i < k; i++) t = norm_aux(gen_0, t, gcoeff(q,i,k), gcoeff(q,i,i));
    1554         182 :   return t;
    1555             : }
    1556             : 
    1557             : /* q is the Cholesky decomposition of a quadratic form
    1558             :  * Enumerate vectors whose norm is less than BORNE (Algo 2.5.7),
    1559             :  * minimal vectors if BORNE = NULL (implies check = NULL).
    1560             :  * If (check != NULL) consider only vectors passing the check, and assumes
    1561             :  *   we only want the smallest possible vectors */
    1562             : static GEN
    1563         867 : smallvectors(GEN q, GEN BORNE, long maxnum, FP_chk_fun *CHECK)
    1564             : {
    1565         867 :   long N = lg(q), n = N-1, i, j, k, s, stockmax, checkcnt = 1;
    1566             :   pari_sp av, av1;
    1567             :   GEN inc, S, x, y, z, v, p1, alpha, norms;
    1568             :   GEN norme1, normax1, borne1, borne2;
    1569         867 :   GEN (*check)(void *,GEN) = CHECK? CHECK->f: NULL;
    1570         867 :   void *data = CHECK? CHECK->data: NULL;
    1571         867 :   const long skipfirst = CHECK? CHECK->skipfirst: 0;
    1572         867 :   const int stockall = (maxnum == -1);
    1573             : 
    1574         867 :   alpha = dbltor(0.95);
    1575         867 :   normax1 = gen_0;
    1576             : 
    1577         867 :   v = cgetg(N,t_VEC);
    1578         867 :   inc = const_vecsmall(n, 1);
    1579             : 
    1580         867 :   av = avma;
    1581         867 :   stockmax = stockall? 2000: maxnum;
    1582         867 :   norms = cgetg(check?(stockmax+1): 1,t_VEC); /* unused if (!check) */
    1583         867 :   S = cgetg(stockmax+1,t_VEC);
    1584         867 :   x = cgetg(N,t_COL);
    1585         867 :   y = cgetg(N,t_COL);
    1586         867 :   z = cgetg(N,t_COL);
    1587        5162 :   for (i=1; i<N; i++) {
    1588        4295 :     gel(v,i) = gcoeff(q,i,i);
    1589        4295 :     gel(x,i) = gel(y,i) = gel(z,i) = gen_0;
    1590             :   }
    1591         867 :   if (BORNE)
    1592             :   {
    1593         853 :     borne1 = BORNE;
    1594         853 :     if (typ(borne1) != t_REAL)
    1595             :     {
    1596             :       long prec;
    1597         426 :       if (gequal0(borne1)) retmkvec3(gen_0, gen_0, cgetg(1,t_MAT));
    1598         419 :       prec = nbits2prec(gexpo(borne1) + 10);
    1599         419 :       borne1 = gtofp(borne1, maxss(prec, DEFAULTPREC));
    1600             :     }
    1601             :   }
    1602             :   else
    1603             :   {
    1604          14 :     borne1 = gcoeff(q,1,1);
    1605         196 :     for (i=2; i<N; i++)
    1606             :     {
    1607         182 :       GEN b = cholesky_norm_ek(q, i);
    1608         182 :       if (gcmp(b, borne1) < 0) borne1 = b;
    1609             :     }
    1610             :     /* borne1 = norm of smallest basis vector */
    1611             :   }
    1612         860 :   borne2 = mulrr(borne1,alpha);
    1613         860 :   if (DEBUGLEVEL>2)
    1614           0 :     err_printf("smallvectors looking for norm < %P.4G\n",borne1);
    1615         860 :   s = 0; k = n;
    1616      170427 :   for(;; step(x,y,inc,k)) /* main */
    1617             :   { /* x (supposedly) small vector, ZV.
    1618             :      * For all t >= k, we have
    1619             :      *   z[t] = sum_{j > t} q[t,j] * x[j]
    1620             :      *   y[t] = sum_{i > t} q[i,i] * (x[i] + z[i])^2
    1621             :      *        = 0 <=> x[i]=0 for all i>t */
    1622             :     do
    1623             :     {
    1624      528585 :       int skip = 0;
    1625      528585 :       if (k > 1)
    1626             :       {
    1627      358158 :         long l = k-1;
    1628      358158 :         av1 = avma;
    1629      358158 :         p1 = mulimp(gel(x,k), gcoeff(q,l,k));
    1630      358158 :         for (j=k+1; j<N; j++) p1 = addmulimp(p1, gel(x,j), gcoeff(q,l,j));
    1631      358158 :         gel(z,l) = gerepileuptoleaf(av1,p1);
    1632             : 
    1633      358158 :         av1 = avma;
    1634      358158 :         p1 = norm_aux(gel(x,k), gel(y,k), gel(z,k), gel(v,k));
    1635      358158 :         gel(y,l) = gerepileuptoleaf(av1, p1);
    1636             :         /* skip the [x_1,...,x_skipfirst,0,...,0] */
    1637      358158 :         if ((l <= skipfirst && !signe(gel(y,skipfirst)))
    1638      357731 :          || mplessthan(borne1, gel(y,l))) skip = 1;
    1639             :         else /* initial value, minimizing (x[l] + z[l])^2, hence qf(x) for
    1640             :                 the given x[1..l-1] */
    1641      357731 :           gel(x,l) = mpround( mpneg(gel(z,l)) );
    1642      358158 :         k = l;
    1643             :       }
    1644      358158 :       for(;; step(x,y,inc,k))
    1645             :       { /* at most 2n loops */
    1646      886743 :         if (!skip)
    1647             :         {
    1648      886316 :           if (check_bound(borne1, gel(x,k),gel(y,k),gel(z,k),gel(v,k))) break;
    1649      387509 :           step(x,y,inc,k);
    1650      387509 :           if (check_bound(borne1, gel(x,k),gel(y,k),gel(z,k),gel(v,k))) break;
    1651             :         }
    1652      359018 :         skip = 0; inc[k] = 1;
    1653      359018 :         if (++k > n) goto END;
    1654      358158 :       }
    1655             : 
    1656      527725 :       if (gc_needed(av,2))
    1657             :       {
    1658           8 :         if(DEBUGMEM>1) pari_warn(warnmem,"smallvectors");
    1659           8 :         if (stockmax) S = clonefill(S, s, stockmax);
    1660           8 :         if (check) {
    1661           8 :           GEN dummy = cgetg(1, t_STR);
    1662           8 :           for (i=s+1; i<=stockmax; i++) gel(norms,i) = dummy;
    1663             :         }
    1664           8 :         gerepileall(av,7,&x,&y,&z,&normax1,&borne1,&borne2,&norms);
    1665             :       }
    1666             :     }
    1667      527725 :     while (k > 1);
    1668      170427 :     if (!signe(gel(x,1)) && !signe(gel(y,1))) continue; /* exclude 0 */
    1669             : 
    1670      169994 :     av1 = avma;
    1671      169994 :     norme1 = norm_aux(gel(x,1),gel(y,1),gel(z,1),gel(v,1));
    1672      169994 :     if (mpgreaterthan(norme1,borne1)) { avma = av1; continue; /* main */ }
    1673             : 
    1674      169994 :     norme1 = gerepileuptoleaf(av1,norme1);
    1675      169994 :     if (check)
    1676             :     {
    1677      101415 :       if (checkcnt < 5 && mpcmp(norme1, borne2) < 0)
    1678             :       {
    1679         426 :         if (!check(data,x)) { checkcnt++ ; continue; /* main */}
    1680         174 :         if (DEBUGLEVEL>4) err_printf("New bound: %Ps", norme1);
    1681         174 :         borne1 = norme1;
    1682         174 :         borne2 = mulrr(borne1, alpha);
    1683         174 :         s = 0; checkcnt = 0;
    1684             :       }
    1685             :     }
    1686             :     else
    1687             :     {
    1688       68579 :       if (!BORNE) /* find minimal vectors */
    1689             :       {
    1690        1883 :         if (mplessthan(norme1, borne1))
    1691             :         { /* strictly smaller vector than previously known */
    1692           0 :           borne1 = norme1; /* + epsilon */
    1693           0 :           s = 0;
    1694             :         }
    1695             :       }
    1696             :       else
    1697       66696 :         if (mpcmp(norme1,normax1) > 0) normax1 = norme1;
    1698             :     }
    1699      169742 :     if (++s > stockmax) continue; /* too many vectors: no longer remember */
    1700      168811 :     if (check) gel(norms,s) = norme1;
    1701      168811 :     gel(S,s) = leafcopy(x);
    1702      168811 :     if (s != stockmax) continue; /* still room, get next vector */
    1703             : 
    1704             :     /* overflow, eliminate vectors failing "check" */
    1705          56 :     if (check)
    1706             :     {
    1707          35 :       pari_sp av2 = avma;
    1708             :       long imin, imax;
    1709          35 :       GEN per = indexsort(norms);
    1710          35 :       if (DEBUGLEVEL>2) err_printf("sorting... [%ld elts]\n",s);
    1711             :       /* let N be the minimal norm so far for x satisfying 'check'. Keep
    1712             :        * all elements of norm N */
    1713        3675 :       for (i = 1; i <= s; i++)
    1714             :       {
    1715        3675 :         long k = per[i];
    1716        3675 :         if (check(data,gel(S,k))) { borne1 = gel(norms,k); break; }
    1717             :       }
    1718          35 :       imin = i;
    1719         273 :       for (; i <= s; i++)
    1720         273 :         if (mpgreaterthan(gel(norms,per[i]), borne1)) break;
    1721          35 :       imax = i;
    1722          35 :       for (i=imin, s=0; i < imax; i++) gel(S,++s) = gel(S,per[i]);
    1723          35 :       avma = av2;
    1724          35 :       if (s)
    1725             :       {
    1726          35 :         borne2 = mulrr(borne1, alpha);
    1727          35 :         checkcnt = 0;
    1728             :       }
    1729          35 :       if (!stockall) continue;
    1730          35 :       if (s > stockmax/2) stockmax <<= 1;
    1731          35 :       norms = cgetg(stockmax+1, t_VEC);
    1732          35 :       for (i = 1; i <= s; i++) gel(norms,i) = borne1;
    1733             :     }
    1734             :     else
    1735             :     {
    1736          21 :       if (!stockall && BORNE) goto END;
    1737          21 :       if (!stockall) continue;
    1738          21 :       stockmax <<= 1;
    1739             :     }
    1740             : 
    1741             :     {
    1742          56 :       GEN Snew = cgetg(stockmax + 1, t_VEC);
    1743          56 :       for (i = 1; i <= s; i++) gel(Snew,i) = gel(S,i);
    1744          56 :       Snew = clonefill(Snew, s, stockmax);
    1745          56 :       if (isclone(S)) gunclone(S);
    1746          56 :       S = Snew;
    1747             :     }
    1748      170427 :   }
    1749             : END:
    1750         860 :   if (s < stockmax) stockmax = s;
    1751         860 :   if (check)
    1752             :   {
    1753             :     GEN per, alph, pols, p;
    1754         839 :     if (DEBUGLEVEL>2) err_printf("final sort & check...\n");
    1755         839 :     setlg(norms,stockmax+1); per = indexsort(norms);
    1756         839 :     alph = cgetg(stockmax+1,t_VEC);
    1757         839 :     pols = cgetg(stockmax+1,t_VEC);
    1758        7812 :     for (j=0,i=1; i<=stockmax; i++)
    1759             :     {
    1760        6993 :       long t = per[i];
    1761        6993 :       GEN N = gel(norms,t);
    1762        6993 :       if (j && mpgreaterthan(N, borne1)) break;
    1763        6973 :       if ((p = check(data,gel(S,t))))
    1764             :       {
    1765        6308 :         if (!j) borne1 = N;
    1766        6308 :         j++;
    1767        6308 :         gel(pols,j) = p;
    1768        6308 :         gel(alph,j) = gel(S,t);
    1769             :       }
    1770             :     }
    1771         839 :     setlg(pols,j+1);
    1772         839 :     setlg(alph,j+1);
    1773         839 :     if (stockmax && isclone(S)) { alph = gcopy(alph); gunclone(S); }
    1774         839 :     return mkvec2(pols, alph);
    1775             :   }
    1776          21 :   if (stockmax)
    1777             :   {
    1778          14 :     setlg(S,stockmax+1);
    1779          14 :     settyp(S,t_MAT);
    1780          14 :     if (isclone(S)) { p1 = S; S = gcopy(S); gunclone(p1); }
    1781             :   }
    1782             :   else
    1783           7 :     S = cgetg(1,t_MAT);
    1784          21 :   return mkvec3(utoi(s<<1), borne1, S);
    1785             : }
    1786             : 
    1787             : /* solve q(x) = x~.a.x <= bound, a > 0.
    1788             :  * If check is non-NULL keep x only if check(x).
    1789             :  * If a is a vector, assume a[1] is the LLL-reduced Cholesky form of q */
    1790             : GEN
    1791         882 : fincke_pohst(GEN a, GEN B0, long stockmax, long PREC, FP_chk_fun *CHECK)
    1792             : {
    1793         882 :   pari_sp av = avma;
    1794             :   VOLATILE long i,j,l;
    1795         882 :   VOLATILE GEN r,rinv,rinvtrans,u,v,res,z,vnorm,rperm,perm,uperm, bound = B0;
    1796             : 
    1797         882 :   if (typ(a) == t_VEC)
    1798             :   {
    1799         435 :     r = gel(a,1);
    1800         435 :     u = NULL;
    1801             :   }
    1802             :   else
    1803             :   {
    1804         447 :     long prec = PREC;
    1805         447 :     l = lg(a);
    1806         447 :     if (l == 1)
    1807             :     {
    1808           7 :       if (CHECK) pari_err_TYPE("fincke_pohst [dimension 0]", a);
    1809           7 :       retmkvec3(gen_0, gen_0, cgetg(1,t_MAT));
    1810             :     }
    1811         440 :     u = lllfp(a, 0.75, LLL_GRAM);
    1812         440 :     if (lg(u) != lg(a)) return NULL;
    1813         440 :     r = qf_apply_RgM(a,u);
    1814         440 :     i = gprecision(r);
    1815         440 :     if (i)
    1816         412 :       prec = i;
    1817             :     else {
    1818          28 :       prec = DEFAULTPREC + nbits2extraprec(gexpo(r));
    1819          28 :       if (prec < PREC) prec = PREC;
    1820             :     }
    1821         440 :     if (DEBUGLEVEL>2) err_printf("first LLL: prec = %ld\n", prec);
    1822         440 :     r = qfgaussred_positive(r);
    1823         440 :     if (!r) return NULL;
    1824        1991 :     for (i=1; i<l; i++)
    1825             :     {
    1826        1551 :       GEN s = gsqrt(gcoeff(r,i,i), prec);
    1827        1551 :       gcoeff(r,i,i) = s;
    1828        1551 :       for (j=i+1; j<l; j++) gcoeff(r,i,j) = gmul(s, gcoeff(r,i,j));
    1829             :     }
    1830             :   }
    1831             :   /* now r~ * r = a in LLL basis */
    1832         875 :   rinv = RgM_inv_upper(r);
    1833         875 :   if (!rinv) return NULL;
    1834         875 :   rinvtrans = shallowtrans(rinv);
    1835         875 :   if (DEBUGLEVEL>2)
    1836           0 :     err_printf("Fincke-Pohst, final LLL: prec = %ld\n", gprecision(rinvtrans));
    1837         875 :   v = lll(rinvtrans);
    1838         875 :   if (lg(v) != lg(rinvtrans)) return NULL;
    1839             : 
    1840         875 :   rinvtrans = RgM_mul(rinvtrans, v);
    1841         875 :   v = ZM_inv(shallowtrans(v),gen_1);
    1842         875 :   r = RgM_mul(r,v);
    1843         875 :   u = u? ZM_mul(u,v): v;
    1844             : 
    1845         875 :   l = lg(r);
    1846         875 :   vnorm = cgetg(l,t_VEC);
    1847         875 :   for (j=1; j<l; j++) gel(vnorm,j) = gnorml2(gel(rinvtrans,j));
    1848         875 :   rperm = cgetg(l,t_MAT);
    1849         875 :   uperm = cgetg(l,t_MAT); perm = indexsort(vnorm);
    1850         875 :   for (i=1; i<l; i++) { uperm[l-i] = u[perm[i]]; rperm[l-i] = r[perm[i]]; }
    1851         875 :   u = uperm;
    1852         875 :   r = rperm; res = NULL;
    1853         875 :   pari_CATCH(e_PREC) { }
    1854             :   pari_TRY {
    1855             :     GEN q;
    1856         875 :     if (CHECK && CHECK->f_init) bound = CHECK->f_init(CHECK, r, u);
    1857         867 :     q = gaussred_from_QR(r, gprecision(vnorm));
    1858         867 :     if (!q) pari_err_PREC("fincke_pohst");
    1859         867 :     res = smallvectors(q, bound, stockmax, CHECK);
    1860         867 :   } pari_ENDCATCH;
    1861         875 :   if (CHECK)
    1862             :   {
    1863         847 :     if (CHECK->f_post) res = CHECK->f_post(CHECK, res, u);
    1864         847 :     return res;
    1865             :   }
    1866          28 :   if (!res) pari_err_PREC("fincke_pohst");
    1867             : 
    1868          28 :   z = cgetg(4,t_VEC);
    1869          28 :   gel(z,1) = gcopy(gel(res,1));
    1870          28 :   gel(z,2) = gcopy(gel(res,2));
    1871          28 :   gel(z,3) = ZM_mul(u, gel(res,3)); return gerepileupto(av,z);
    1872             : }

Generated by: LCOV version 1.11