Testing latest pari + WASM + node.js... and it works?! Wow.
License: GPL3
ubuntu2004
/* Copyright (C) 2000 The PARI group.12This file is part of the PARI/GP package.34PARI/GP is free software; you can redistribute it and/or modify it under the5terms of the GNU General Public License as published by the Free Software6Foundation; either version 2 of the License, or (at your option) any later7version. It is distributed in the hope that it will be useful, but WITHOUT8ANY WARRANTY WHATSOEVER.910Check the License for details. You should have received a copy of it, along11with the package; see the file 'COPYING'. If not, write to the Free Software12Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */1314/********************************************************************/15/** **/16/** LLL Algorithm and close friends **/17/** **/18/********************************************************************/19#include "pari.h"20#include "paripriv.h"2122#define DEBUGLEVEL DEBUGLEVEL_qf2324/********************************************************************/25/** QR Factorization via Householder matrices **/26/********************************************************************/27static int28no_prec_pb(GEN x)29{30return (typ(x) != t_REAL || realprec(x) > LOWDEFAULTPREC31|| expo(x) < BITS_IN_LONG/2);32}33/* Find a Householder transformation which, applied to x[k..#x], zeroes34* x[k+1..#x]; fill L = (mu_{i,j}). Return 0 if precision problem [obtained35* a 0 vector], 1 otherwise */36static int37FindApplyQ(GEN x, GEN L, GEN B, long k, GEN Q, long prec)38{39long i, lx = lg(x)-1;40GEN x2, x1, xd = x + (k-1);4142x1 = gel(xd,1);43x2 = mpsqr(x1);44if (k < lx)45{46long lv = lx - (k-1) + 1;47GEN beta, Nx, v = cgetg(lv, t_VEC);48for (i=2; i<lv; i++) {49x2 = mpadd(x2, mpsqr(gel(xd,i)));50gel(v,i) = gel(xd,i);51}52if (!signe(x2)) return 0;53Nx = gsqrt(x2, prec); if (signe(x1) < 0) setsigne(Nx, -1);54gel(v,1) = mpadd(x1, Nx);5556if (!signe(x1))57beta = gtofp(x2, prec); /* make sure typ(beta) != t_INT */58else59beta = mpadd(x2, mpmul(Nx,x1));60gel(Q,k) = mkvec2(invr(beta), v);6162togglesign(Nx);63gcoeff(L,k,k) = Nx;64}65else66gcoeff(L,k,k) = gel(x,k);67gel(B,k) = x2;68for (i=1; i<k; i++) gcoeff(L,k,i) = gel(x,i);69return no_prec_pb(x2);70}7172/* apply Householder transformation Q = [beta,v] to r with t_INT/t_REAL73* coefficients, in place: r -= ((0|v).r * beta) v */74static void75ApplyQ(GEN Q, GEN r)76{77GEN s, rd, beta = gel(Q,1), v = gel(Q,2);78long i, l = lg(v), lr = lg(r);7980rd = r + (lr - l);81s = mpmul(gel(v,1), gel(rd,1));82for (i=2; i<l; i++) s = mpadd(s, mpmul(gel(v,i) ,gel(rd,i)));83s = mpmul(beta, s);84for (i=1; i<l; i++)85if (signe(gel(v,i))) gel(rd,i) = mpsub(gel(rd,i), mpmul(s, gel(v,i)));86}87/* apply Q[1], ..., Q[j-1] to r */88static GEN89ApplyAllQ(GEN Q, GEN r, long j)90{91pari_sp av = avma;92long i;93r = leafcopy(r);94for (i=1; i<j; i++) ApplyQ(gel(Q,i), r);95return gerepilecopy(av, r);96}9798/* same, arbitrary coefficients [20% slower for t_REAL at DEFAULTPREC] */99static void100RgC_ApplyQ(GEN Q, GEN r)101{102GEN s, rd, beta = gel(Q,1), v = gel(Q,2);103long i, l = lg(v), lr = lg(r);104105rd = r + (lr - l);106s = gmul(gel(v,1), gel(rd,1));107for (i=2; i<l; i++) s = gadd(s, gmul(gel(v,i) ,gel(rd,i)));108s = gmul(beta, s);109for (i=1; i<l; i++)110if (signe(gel(v,i))) gel(rd,i) = gsub(gel(rd,i), gmul(s, gel(v,i)));111}112static GEN113RgC_ApplyAllQ(GEN Q, GEN r, long j)114{115pari_sp av = avma;116long i;117r = leafcopy(r);118for (i=1; i<j; i++) RgC_ApplyQ(gel(Q,i), r);119return gerepilecopy(av, r);120}121122int123RgM_QR_init(GEN x, GEN *pB, GEN *pQ, GEN *pL, long prec)124{125x = RgM_gtomp(x, prec);126return QR_init(x, pB, pQ, pL, prec);127}128129static void130check_householder(GEN Q)131{132long i, l = lg(Q);133if (typ(Q) != t_VEC) pari_err_TYPE("mathouseholder", Q);134for (i = 1; i < l; i++)135{136GEN q = gel(Q,i), v;137if (typ(q) != t_VEC || lg(q) != 3) pari_err_TYPE("mathouseholder", Q);138v = gel(q,2);139if (typ(v) != t_VEC || lg(v)+i-2 != l) pari_err_TYPE("mathouseholder", Q);140}141}142143GEN144mathouseholder(GEN Q, GEN v)145{146long l = lg(Q);147check_householder(Q);148switch(typ(v))149{150case t_MAT:151{152long lx, i;153GEN M = cgetg_copy(v, &lx);154if (lx == 1) return M;155if (lgcols(v) != l+1) pari_err_TYPE("mathouseholder", v);156for (i = 1; i < lx; i++) gel(M,i) = RgC_ApplyAllQ(Q, gel(v,i), l);157return M;158}159case t_COL: if (lg(v) == l+1) break;160/* fall through */161default: pari_err_TYPE("mathouseholder", v);162}163return RgC_ApplyAllQ(Q, v, l);164}165166GEN167matqr(GEN x, long flag, long prec)168{169pari_sp av = avma;170GEN B, Q, L;171long n = lg(x)-1;172if (typ(x) != t_MAT) pari_err_TYPE("matqr",x);173if (!n)174{175if (!flag) retmkvec2(cgetg(1,t_MAT),cgetg(1,t_MAT));176retmkvec2(cgetg(1,t_VEC),cgetg(1,t_MAT));177}178if (n != nbrows(x)) pari_err_DIM("matqr");179if (!RgM_QR_init(x, &B,&Q,&L, prec)) pari_err_PREC("matqr");180if (!flag) Q = shallowtrans(mathouseholder(Q, matid(n)));181return gerepilecopy(av, mkvec2(Q, shallowtrans(L)));182}183184/* compute B = | x[k] |^2, Q = Householder transforms and L = mu_{i,j} */185int186QR_init(GEN x, GEN *pB, GEN *pQ, GEN *pL, long prec)187{188long j, k = lg(x)-1;189GEN B = cgetg(k+1, t_VEC), Q = cgetg(k, t_VEC), L = zeromatcopy(k,k);190for (j=1; j<=k; j++)191{192GEN r = gel(x,j);193if (j > 1) r = ApplyAllQ(Q, r, j);194if (!FindApplyQ(r, L, B, j, Q, prec)) return 0;195}196*pB = B; *pQ = Q; *pL = L; return 1;197}198/* x a square t_MAT with t_INT / t_REAL entries and maximal rank. Return199* qfgaussred(x~*x) */200GEN201gaussred_from_QR(GEN x, long prec)202{203long j, k = lg(x)-1;204GEN B, Q, L;205if (!QR_init(x, &B,&Q,&L, prec)) return NULL;206for (j=1; j<k; j++)207{208GEN m = gel(L,j), invNx = invr(gel(m,j));209long i;210gel(m,j) = gel(B,j);211for (i=j+1; i<=k; i++) gel(m,i) = mpmul(invNx, gel(m,i));212}213gcoeff(L,j,j) = gel(B,j);214return shallowtrans(L);215}216GEN217R_from_QR(GEN x, long prec)218{219GEN B, Q, L;220if (!QR_init(x, &B,&Q,&L, prec)) return NULL;221return shallowtrans(L);222}223224/********************************************************************/225/** QR Factorization via Gram-Schmidt **/226/********************************************************************/227228/* return Gram-Schmidt orthogonal basis (f) attached to (e), B is the229* vector of the (f_i . f_i)*/230GEN231RgM_gram_schmidt(GEN e, GEN *ptB)232{233long i,j,lx = lg(e);234GEN f = RgM_shallowcopy(e), B, iB;235236B = cgetg(lx, t_VEC);237iB= cgetg(lx, t_VEC);238239for (i=1;i<lx;i++)240{241GEN p1 = NULL;242pari_sp av = avma;243for (j=1; j<i; j++)244{245GEN mu = gmul(RgV_dotproduct(gel(e,i),gel(f,j)), gel(iB,j));246GEN p2 = gmul(mu, gel(f,j));247p1 = p1? gadd(p1,p2): p2;248}249p1 = p1? gerepileupto(av, gsub(gel(e,i), p1)): gel(e,i);250gel(f,i) = p1;251gel(B,i) = RgV_dotsquare(gel(f,i));252gel(iB,i) = ginv(gel(B,i));253}254*ptB = B; return f;255}256257/* B a Z-basis (which the caller should LLL-reduce for efficiency), t a vector.258* Apply Babai's nearest plane algorithm to (B,t) */259GEN260RgM_Babai(GEN B, GEN t)261{262GEN C, N, G = RgM_gram_schmidt(B, &N), b = t;263long j, n = lg(B)-1;264265C = cgetg(n+1,t_COL);266for (j = n; j > 0; j--)267{268GEN c = gdiv( RgV_dotproduct(b, gel(G,j)), gel(N,j) );269long e;270c = grndtoi(c,&e);271if (e >= 0) return NULL;272if (signe(c)) b = RgC_sub(b, RgC_Rg_mul(gel(B,j), c));273gel(C,j) = c;274}275return C;276}277278/********************************************************************/279/** **/280/** LLL ALGORITHM **/281/** **/282/********************************************************************/283/* Def: a matrix M is said to be -partially reduced- if | m1 +- m2 | >= |m1|284* for any two columns m1 != m2, in M.285*286* Input: an integer matrix mat whose columns are linearly independent. Find287* another matrix T such that mat * T is partially reduced.288*289* Output: mat * T if flag = 0; T if flag != 0,290*291* This routine is designed to quickly reduce lattices in which one row292* is huge compared to the other rows. For example, when searching for a293* polynomial of degree 3 with root a mod N, the four input vectors might294* be the coefficients of295* X^3 - (a^3 mod N), X^2 - (a^2 mod N), X - (a mod N), N.296* All four constant coefficients are O(p) and the rest are O(1). By the297* pigeon-hole principle, the coefficients of the smallest vector in the298* lattice are O(p^(1/4)), hence significant reduction of vector lengths299* can be anticipated.300*301* An improved algorithm would look only at the leading digits of dot*. It302* would use single-precision calculations as much as possible.303*304* Original code: Peter Montgomery (1994) */305static GEN306lllintpartialall(GEN m, long flag)307{308const long ncol = lg(m)-1;309const pari_sp av = avma;310GEN tm1, tm2, mid;311312if (ncol <= 1) return flag? matid(ncol): gcopy(m);313314tm1 = flag? matid(ncol): NULL;315{316const pari_sp av2 = avma;317GEN dot11 = ZV_dotsquare(gel(m,1));318GEN dot22 = ZV_dotsquare(gel(m,2));319GEN dot12 = ZV_dotproduct(gel(m,1), gel(m,2));320GEN tm = matid(2); /* For first two columns only */321322int progress = 0;323long npass2 = 0;324325/* Row reduce the first two columns of m. Our best result so far is326* (first two columns of m)*tm.327*328* Initially tm = 2 x 2 identity matrix.329* Inner products of the reduced matrix are in dot11, dot12, dot22. */330while (npass2 < 2 || progress)331{332GEN dot12new, q = diviiround(dot12, dot22);333334npass2++; progress = signe(q);335if (progress)336{/* Conceptually replace (v1, v2) by (v1 - q*v2, v2), where v1 and v2337* represent the reduced basis for the first two columns of the matrix.338* We do this by updating tm and the inner products. */339togglesign(q);340dot12new = addii(dot12, mulii(q, dot22));341dot11 = addii(dot11, mulii(q, addii(dot12, dot12new)));342dot12 = dot12new;343ZC_lincomb1_inplace(gel(tm,1), gel(tm,2), q);344}345346/* Interchange the output vectors v1 and v2. */347swap(dot11,dot22);348swap(gel(tm,1), gel(tm,2));349350/* Occasionally (including final pass) do garbage collection. */351if ((npass2 & 0xff) == 0 || !progress)352gerepileall(av2, 4, &dot11,&dot12,&dot22,&tm);353} /* while npass2 < 2 || progress */354355{356long i;357GEN det12 = subii(mulii(dot11, dot22), sqri(dot12));358359mid = cgetg(ncol+1, t_MAT);360for (i = 1; i <= 2; i++)361{362GEN tmi = gel(tm,i);363if (tm1)364{365GEN tm1i = gel(tm1,i);366gel(tm1i,1) = gel(tmi,1);367gel(tm1i,2) = gel(tmi,2);368}369gel(mid,i) = ZC_lincomb(gel(tmi,1),gel(tmi,2), gel(m,1),gel(m,2));370}371for (i = 3; i <= ncol; i++)372{373GEN c = gel(m,i);374GEN dot1i = ZV_dotproduct(gel(mid,1), c);375GEN dot2i = ZV_dotproduct(gel(mid,2), c);376/* ( dot11 dot12 ) (q1) ( dot1i )377* ( dot12 dot22 ) (q2) = ( dot2i )378*379* Round -q1 and -q2 to nearest integer. Then compute380* c - q1*mid[1] - q2*mid[2].381* This will be approximately orthogonal to the first two vectors in382* the new basis. */383GEN q1neg = subii(mulii(dot12, dot2i), mulii(dot22, dot1i));384GEN q2neg = subii(mulii(dot12, dot1i), mulii(dot11, dot2i));385386q1neg = diviiround(q1neg, det12);387q2neg = diviiround(q2neg, det12);388if (tm1)389{390gcoeff(tm1,1,i) = addii(mulii(q1neg, gcoeff(tm,1,1)),391mulii(q2neg, gcoeff(tm,1,2)));392gcoeff(tm1,2,i) = addii(mulii(q1neg, gcoeff(tm,2,1)),393mulii(q2neg, gcoeff(tm,2,2)));394}395gel(mid,i) = ZC_add(c, ZC_lincomb(q1neg,q2neg, gel(mid,1),gel(mid,2)));396} /* for i */397} /* local block */398}399if (DEBUGLEVEL>6)400{401if (tm1) err_printf("tm1 = %Ps",tm1);402err_printf("mid = %Ps",mid); /* = m * tm1 */403}404gerepileall(av, tm1? 2: 1, &mid, &tm1);405{406/* For each pair of column vectors v and w in mid * tm2,407* try to replace (v, w) by (v, v - q*w) for some q.408* We compute all inner products and check them repeatedly. */409const pari_sp av3 = avma;410long i, j, npass = 0, e = LONG_MAX;411GEN dot = cgetg(ncol+1, t_MAT); /* scalar products */412413tm2 = matid(ncol);414for (i=1; i <= ncol; i++)415{416gel(dot,i) = cgetg(ncol+1,t_COL);417for (j=1; j <= i; j++)418gcoeff(dot,j,i) = gcoeff(dot,i,j) = ZV_dotproduct(gel(mid,i),gel(mid,j));419}420for(;;)421{422long reductions = 0, olde = e;423for (i=1; i <= ncol; i++)424{425long ijdif;426for (ijdif=1; ijdif < ncol; ijdif++)427{428long d, k1, k2;429GEN codi, q;430431j = i + ijdif; if (j > ncol) j -= ncol;432/* let k1, resp. k2, index of larger, resp. smaller, column */433if (cmpii(gcoeff(dot,i,i), gcoeff(dot,j,j)) > 0) { k1 = i; k2 = j; }434else { k1 = j; k2 = i; }435codi = gcoeff(dot,k2,k2);436q = signe(codi)? diviiround(gcoeff(dot,k1,k2), codi): gen_0;437if (!signe(q)) continue;438439/* Try to subtract a multiple of column k2 from column k1. */440reductions++; togglesign_safe(&q);441ZC_lincomb1_inplace(gel(tm2,k1), gel(tm2,k2), q);442ZC_lincomb1_inplace(gel(dot,k1), gel(dot,k2), q);443gcoeff(dot,k1,k1) = addii(gcoeff(dot,k1,k1),444mulii(q, gcoeff(dot,k2,k1)));445for (d = 1; d <= ncol; d++) gcoeff(dot,k1,d) = gcoeff(dot,d,k1);446} /* for ijdif */447if (gc_needed(av3,2))448{449if(DEBUGMEM>1) pari_warn(warnmem,"lllintpartialall");450gerepileall(av3, 2, &dot,&tm2);451}452} /* for i */453if (!reductions) break;454e = 0;455for (i = 1; i <= ncol; i++) e += expi( gcoeff(dot,i,i) );456if (e == olde) break;457if (DEBUGLEVEL>6)458{459npass++;460err_printf("npass = %ld, red. last time = %ld, log_2(det) ~ %ld\n\n",461npass, reductions, e);462}463} /* for(;;)*/464465/* Sort columns so smallest comes first in m * tm1 * tm2.466* Use insertion sort. */467for (i = 1; i < ncol; i++)468{469long j, s = i;470471for (j = i+1; j <= ncol; j++)472if (cmpii(gcoeff(dot,s,s),gcoeff(dot,j,j)) > 0) s = j;473if (i != s)474{ /* Exchange with proper column; only the diagonal of dot is updated */475swap(gel(tm2,i), gel(tm2,s));476swap(gcoeff(dot,i,i), gcoeff(dot,s,s));477}478}479} /* local block */480return gerepileupto(av, ZM_mul(tm1? tm1: mid, tm2));481}482483GEN484lllintpartial(GEN mat) { return lllintpartialall(mat,1); }485486GEN487lllintpartial_inplace(GEN mat) { return lllintpartialall(mat,0); }488489/********************************************************************/490/** **/491/** COPPERSMITH ALGORITHM **/492/** Finding small roots of univariate equations. **/493/** **/494/********************************************************************/495496static int497check(double b, double x, double rho, long d, long dim, long delta, long t)498{499double cond = delta * (d * (delta+1) - 2*b*dim + rho * (delta-1 + 2*t))500+ x*dim*(dim - 1);501if (DEBUGLEVEL >= 4)502err_printf("delta = %d, t = %d (%.1lf)\n", delta, t, cond);503return (cond <= 0);504}505506static void507choose_params(GEN P, GEN N, GEN X, GEN B, long *pdelta, long *pt)508{509long d = degpol(P), dim;510GEN P0 = leading_coeff(P);511double logN = gtodouble(glog(N, DEFAULTPREC)), x, b, rho;512x = gtodouble(glog(X, DEFAULTPREC)) / logN;513b = B? gtodouble(glog(B, DEFAULTPREC)) / logN: 1.;514if (x * d >= b * b) pari_err_OVERFLOW("zncoppersmith [bound too large]");515/* TODO : remove P0 completely */516rho = is_pm1(P0)? 0: gtodouble(glog(P0, DEFAULTPREC)) / logN;517518/* Enumerate (delta,t) by increasing lattice dimension */519for(dim = d + 1;; dim++)520{521long delta, t; /* dim = d*delta + t in the loop */522for (delta = 0, t = dim; t >= 0; delta++, t -= d)523if (check(b,x,rho,d,dim,delta,t)) { *pdelta = delta; *pt = t; return; }524}525}526527static int528sol_OK(GEN x, GEN N, GEN B)529{ return B? (cmpii(gcdii(x,N),B) >= 0): dvdii(x,N); }530/* deg(P) > 0, x >= 0. Find all j such that gcd(P(j), N) >= B, |j| <= x */531static GEN532do_exhaustive(GEN P, GEN N, long x, GEN B)533{534GEN Pe, Po, sol = vecsmalltrunc_init(2*x + 2);535pari_sp av;536long j;537RgX_even_odd(P, &Pe,&Po); av = avma;538if (sol_OK(gel(P,2), N,B)) vecsmalltrunc_append(sol, 0);539for (j = 1; j <= x; j++, set_avma(av))540{541GEN j2 = sqru(j), E = FpX_eval(Pe,j2,N), O = FpX_eval(Po,j2,N);542if (sol_OK(addmuliu(E,O,j), N,B)) vecsmalltrunc_append(sol, j);543if (sol_OK(submuliu(E,O,j), N,B)) vecsmalltrunc_append(sol,-j);544}545vecsmall_sort(sol); return zv_to_ZV(sol);546}547548/* General Coppersmith, look for a root x0 <= p, p >= B, p | N, |x0| <= X.549* B = N coded as NULL */550GEN551zncoppersmith(GEN P, GEN N, GEN X, GEN B)552{553GEN Q, R, N0, M, sh, short_pol, *Xpowers, sol, nsp, cP, Z;554long delta, i, j, row, d, l, t, dim, bnd = 10;555const ulong X_SMALL = 1000;556pari_sp av = avma;557558if (typ(P) != t_POL || !RgX_is_ZX(P)) pari_err_TYPE("zncoppersmith",P);559if (typ(N) != t_INT) pari_err_TYPE("zncoppersmith",N);560if (typ(X) != t_INT) {561X = gfloor(X);562if (typ(X) != t_INT) pari_err_TYPE("zncoppersmith",X);563}564if (signe(X) < 0) pari_err_DOMAIN("zncoppersmith", "X", "<", gen_0, X);565P = FpX_red(P, N); d = degpol(P);566if (d == 0) { set_avma(av); return cgetg(1, t_VEC); }567if (d < 0) pari_err_ROOTS0("zncoppersmith");568if (B && typ(B) != t_INT) B = gceil(B);569if (abscmpiu(X, X_SMALL) <= 0)570return gerepileupto(av, do_exhaustive(P, N, itos(X), B));571572if (B && equalii(B,N)) B = NULL;573if (B) bnd = 1; /* bnd-hack is only for the case B = N */574cP = gel(P,d+2);575if (!gequal1(cP))576{577GEN r, z;578gel(P,d+2) = cP = bezout(cP, N, &z, &r);579for (j = 0; j < d; j++) gel(P,j+2) = Fp_mul(gel(P,j+2), z, N);580if (!is_pm1(cP))581{582P = Q_primitive_part(P, &cP);583if (cP) { N = diviiexact(N,cP); B = gceil(gdiv(B, cP)); }584}585}586if (DEBUGLEVEL >= 2) err_printf("Modified P: %Ps\n", P);587588choose_params(P, N, X, B, &delta, &t);589if (DEBUGLEVEL >= 2)590err_printf("Init: trying delta = %d, t = %d\n", delta, t);591for(;;)592{593dim = d * delta + t;594/* TODO: In case of failure do not recompute the full vector */595Xpowers = (GEN*)new_chunk(dim + 1);596Xpowers[0] = gen_1;597for (j = 1; j <= dim; j++) Xpowers[j] = mulii(Xpowers[j-1], X);598599/* TODO: in case of failure, use the part of the matrix already computed */600M = zeromatcopy(dim,dim);601602/* Rows of M correspond to the polynomials603* N^delta, N^delta Xi, ... N^delta (Xi)^d-1,604* N^(delta-1)P(Xi), N^(delta-1)XiP(Xi), ... N^(delta-1)P(Xi)(Xi)^d-1,605* ...606* P(Xi)^delta, XiP(Xi)^delta, ..., P(Xi)^delta(Xi)^t-1 */607for (j = 1; j <= d; j++) gcoeff(M, j, j) = gel(Xpowers,j-1);608609/* P-part */610if (delta) row = d + 1; else row = 0;611612Q = P;613for (i = 1; i < delta; i++)614{615for (j = 0; j < d; j++,row++)616for (l = j + 1; l <= row; l++)617gcoeff(M, l, row) = mulii(Xpowers[l-1], gel(Q,l-j+1));618Q = ZX_mul(Q, P);619}620for (j = 0; j < t; row++, j++)621for (l = j + 1; l <= row; l++)622gcoeff(M, l, row) = mulii(Xpowers[l-1], gel(Q,l-j+1));623624/* N-part */625row = dim - t; N0 = N;626while (row >= 1)627{628for (j = 0; j < d; j++,row--)629for (l = 1; l <= row; l++)630gcoeff(M, l, row) = mulii(gmael(M, row, l), N0);631if (row >= 1) N0 = mulii(N0, N);632}633/* Z is the upper bound for the L^1 norm of the polynomial,634ie. N^delta if B = N, B^delta otherwise */635if (B) Z = powiu(B, delta); else Z = N0;636637if (DEBUGLEVEL >= 2)638{639if (DEBUGLEVEL >= 6) err_printf("Matrix to be reduced:\n%Ps\n", M);640err_printf("Entering LLL\nbitsize bound: %ld\n", expi(Z));641err_printf("expected shvector bitsize: %ld\n", expi(ZM_det_triangular(M))/dim);642}643644sh = ZM_lll(M, 0.75, LLL_INPLACE);645/* Take the first vector if it is non constant */646short_pol = gel(sh,1);647if (ZV_isscalar(short_pol)) short_pol = gel(sh, 2);648649nsp = gen_0;650for (j = 1; j <= dim; j++) nsp = addii(nsp, absi_shallow(gel(short_pol,j)));651652if (DEBUGLEVEL >= 2)653{654err_printf("Candidate: %Ps\n", short_pol);655err_printf("bitsize Norm: %ld\n", expi(nsp));656err_printf("bitsize bound: %ld\n", expi(mului(bnd, Z)));657}658if (cmpii(nsp, mului(bnd, Z)) < 0) break; /* SUCCESS */659660/* Failed with the precomputed or supplied value */661if (++t == d) { delta++; t = 1; }662if (DEBUGLEVEL >= 2)663err_printf("Increasing dim, delta = %d t = %d\n", delta, t);664}665bnd = itos(divii(nsp, Z)) + 1;666667while (!signe(gel(short_pol,dim))) dim--;668669R = cgetg(dim + 2, t_POL); R[1] = P[1];670for (j = 1; j <= dim; j++)671gel(R,j+1) = diviiexact(gel(short_pol,j), Xpowers[j-1]);672gel(R,2) = subii(gel(R,2), mului(bnd - 1, N0));673674sol = cgetg(1, t_VEC);675for (i = -bnd + 1; i < bnd; i++)676{677GEN r = nfrootsQ(R);678if (DEBUGLEVEL >= 2) err_printf("Roots: %Ps\n", r);679for (j = 1; j < lg(r); j++)680{681GEN z = gel(r,j);682if (typ(z) == t_INT && sol_OK(FpX_eval(P,z,N), N,B))683sol = shallowconcat(sol, z);684}685if (i < bnd) gel(R,2) = addii(gel(R,2), Z);686}687return gerepileupto(av, ZV_sort_uniq(sol));688}689690/********************************************************************/691/** **/692/** LINEAR & ALGEBRAIC DEPENDENCE **/693/** **/694/********************************************************************/695696static int697real_indep(GEN re, GEN im, long bit)698{699GEN d = gsub(gmul(gel(re,1),gel(im,2)), gmul(gel(re,2),gel(im,1)));700return (!gequal0(d) && gexpo(d) > - bit);701}702703GEN704lindepfull_bit(GEN x, long bit)705{706long lx = lg(x), ly, i, j;707GEN re, im, M;708709if (! is_vec_t(typ(x))) pari_err_TYPE("lindep2",x);710if (lx <= 2)711{712if (lx == 2 && gequal0(x)) return matid(1);713return NULL;714}715re = real_i(x);716im = imag_i(x);717/* independent over R ? */718if (lx == 3 && real_indep(re,im,bit)) return NULL;719if (gequal0(im)) im = NULL;720ly = im? lx+2: lx+1;721M = cgetg(lx,t_MAT);722for (i=1; i<lx; i++)723{724GEN c = cgetg(ly,t_COL); gel(M,i) = c;725for (j=1; j<lx; j++) gel(c,j) = gen_0;726gel(c,i) = gen_1;727gel(c,lx) = gtrunc2n(gel(re,i), bit);728if (im) gel(c,lx+1) = gtrunc2n(gel(im,i), bit);729}730return ZM_lll(M, 0.99, LLL_INPLACE);731}732GEN733lindep_bit(GEN x, long bit)734{735pari_sp av = avma;736GEN v, M = lindepfull_bit(x,bit);737if (!M) { set_avma(av); return cgetg(1, t_COL); }738v = gel(M,1); setlg(v, lg(M));739return gerepilecopy(av, v);740}741/* deprecated */742GEN743lindep2(GEN x, long dig)744{745long bit;746if (dig < 0) pari_err_DOMAIN("lindep2", "accuracy", "<", gen_0, stoi(dig));747if (dig) bit = (long) (dig/LOG10_2);748else749{750bit = gprecision(x);751if (!bit)752{753x = Q_primpart(x); /* left on stack */754bit = 32 + gexpo(x);755}756else757bit = (long)prec2nbits_mul(bit, 0.8);758}759return lindep_bit(x, bit);760}761762/* x is a vector of elts of a p-adic field */763GEN764lindep_padic(GEN x)765{766long i, j, prec = LONG_MAX, nx = lg(x)-1, v;767pari_sp av = avma;768GEN p = NULL, pn, m, a;769770if (nx < 2) return cgetg(1,t_COL);771for (i=1; i<=nx; i++)772{773GEN c = gel(x,i), q;774if (typ(c) != t_PADIC) continue;775776j = precp(c); if (j < prec) prec = j;777q = gel(c,2);778if (!p) p = q; else if (!equalii(p, q)) pari_err_MODULUS("lindep_padic", p, q);779}780if (!p) pari_err_TYPE("lindep_padic [not a p-adic vector]",x);781v = gvaluation(x,p); pn = powiu(p,prec);782if (v) x = gmul(x, powis(p, -v));783x = RgV_to_FpV(x, pn);784785a = negi(gel(x,1));786m = cgetg(nx,t_MAT);787for (i=1; i<nx; i++)788{789GEN c = zerocol(nx);790gel(c,1+i) = a;791gel(c,1) = gel(x,i+1);792gel(m,i) = c;793}794m = ZM_lll(ZM_hnfmodid(m, pn), 0.99, LLL_INPLACE);795return gerepilecopy(av, gel(m,1));796}797/* x is a vector of t_POL/t_SER */798GEN799lindep_Xadic(GEN x)800{801long i, prec = LONG_MAX, deg = 0, lx = lg(x), vx, v;802pari_sp av = avma;803GEN m;804805if (lx == 1) return cgetg(1,t_COL);806vx = gvar(x);807v = gvaluation(x, pol_x(vx));808if (!v) x = shallowcopy(x);809else if (v > 0) x = gdiv(x, pol_xn(v, vx));810else x = gmul(x, pol_xn(-v, vx));811/* all t_SER have valuation >= 0 */812for (i=1; i<lx; i++)813{814GEN c = gel(x,i);815if (gvar(c) != vx) { gel(x,i) = scalarpol_shallow(c, vx); continue; }816switch(typ(c))817{818case t_POL: deg = maxss(deg, degpol(c)); break;819case t_RFRAC: pari_err_TYPE("lindep_Xadic", c);820case t_SER:821prec = minss(prec, valp(c)+lg(c)-2);822gel(x,i) = ser2rfrac_i(c);823}824}825if (prec == LONG_MAX) prec = deg+1;826m = RgXV_to_RgM(x, prec);827return gerepileupto(av, deplin(m));828}829static GEN830vec_lindep(GEN x)831{832pari_sp av = avma;833long i, l = lg(x); /* > 1 */834long t = typ(gel(x,1)), h = lg(gel(x,1));835GEN m = cgetg(l, t_MAT);836for (i = 1; i < l; i++)837{838GEN y = gel(x,i);839if (lg(y) != h || typ(y) != t) pari_err_TYPE("lindep",x);840if (t != t_COL) y = shallowtrans(y); /* Sigh */841gel(m,i) = y;842}843return gerepileupto(av, deplin(m));844}845846GEN847lindep(GEN x) { return lindep2(x, 0); }848849GEN850lindep0(GEN x,long bit)851{852long i, tx = typ(x);853if (tx == t_MAT) return deplin(x);854if (! is_vec_t(tx)) pari_err_TYPE("lindep",x);855for (i = 1; i < lg(x); i++)856switch(typ(gel(x,i)))857{858case t_PADIC: return lindep_padic(x);859case t_POL:860case t_RFRAC:861case t_SER: return lindep_Xadic(x);862case t_VEC:863case t_COL: return vec_lindep(x);864}865return lindep2(x, bit);866}867868GEN869algdep0(GEN x, long n, long bit)870{871long tx = typ(x), i;872pari_sp av;873GEN y;874875if (! is_scalar_t(tx)) pari_err_TYPE("algdep0",x);876if (tx==t_POLMOD) { y = RgX_copy(gel(x,1)); setvarn(y,0); return y; }877if (gequal0(x)) return pol_x(0);878if (n <= 0)879{880if (!n) return gen_1;881pari_err_DOMAIN("algdep", "degree", "<", gen_0, stoi(n));882}883884av = avma; y = cgetg(n+2,t_COL);885gel(y,1) = gen_1;886gel(y,2) = x; /* n >= 1 */887for (i=3; i<=n+1; i++) gel(y,i) = gmul(gel(y,i-1),x);888if (typ(x) == t_PADIC)889y = lindep_padic(y);890else891y = lindep2(y, bit);892if (lg(y) == 1) pari_err(e_DOMAIN,"algdep", "degree(x)",">", stoi(n), x);893y = RgV_to_RgX(y, 0);894if (signe(leading_coeff(y)) > 0) return gerepilecopy(av, y);895return gerepileupto(av, ZX_neg(y));896}897898GEN899algdep(GEN x, long n)900{901return algdep0(x,n,0);902}903904GEN905seralgdep(GEN s, long p, long r)906{907pari_sp av = avma;908long vy, i, m, n, prec;909GEN S, v, D;910911if (typ(s) != t_SER) pari_err_TYPE("seralgdep",s);912if (p <= 0) pari_err_DOMAIN("seralgdep", "p", "<=", gen_0, stoi(p));913if (r < 0) pari_err_DOMAIN("seralgdep", "r", "<", gen_0, stoi(r));914if (is_bigint(addiu(muluu(p, r), 1))) pari_err_OVERFLOW("seralgdep");915vy = varn(s);916if (!vy) pari_err_PRIORITY("seralgdep", s, ">", 0);917r++; p++;918prec = valp(s) + lg(s)-2;919if (r > prec) r = prec;920S = cgetg(p+1, t_VEC); gel(S, 1) = s;921for (i = 2; i <= p; i++) gel(S,i) = gmul(gel(S,i-1), s);922v = cgetg(r*p+1, t_VEC); /* v[r*n+m+1] = s^n * y^m */923/* n = 0 */924for (m = 0; m < r; m++) gel(v, m+1) = pol_xn(m, vy);925for(n=1; n < p; n++)926for (m = 0; m < r; m++)927{928GEN c = gel(S,n);929if (m)930{931c = shallowcopy(c);932setvalp(c, valp(c) + m);933}934gel(v, r*n + m + 1) = c;935}936D = lindep_Xadic(v);937if (lg(D) == 1) { set_avma(av); return gen_0; }938v = cgetg(p+1, t_VEC);939for (n = 0; n < p; n++)940gel(v, n+1) = RgV_to_RgX(vecslice(D, r*n+1, r*n+r), vy);941return gerepilecopy(av, RgV_to_RgX(v, 0));942}943944/* FIXME: could precompute ZM_lll attached to V[2..] */945static GEN946lindepcx(GEN V, long bit)947{948GEN Vr = real_i(V), Vi = imag_i(V);949if (gexpo(Vr) < -bit) V = Vi;950else if (gexpo(Vi) < -bit) V = Vr;951return lindepfull_bit(V, bit);952}953/* c floating point t_REAL or t_COMPLEX, T ZX, recognize in Q[x]/(T).954* V helper vector (containing complex roots of T), MODIFIED */955static GEN956cx_bestapprnf(GEN c, GEN T, GEN V, long bit)957{958GEN M, a, v = NULL;959long i, l;960gel(V,1) = gneg(c); M = lindepcx(V, bit);961if (!M) pari_err(e_MISC, "cannot rationalize coeff in bestapprnf");962l = lg(M); a = NULL;963for (i = 1; i < l; i ++) { v = gel(M,i); a = gel(v,1); if (signe(a)) break; }964v = RgC_Rg_div(vecslice(v, 2, lg(M)-1), a);965if (!T) return gel(v,1);966v = RgV_to_RgX(v, varn(T)); l = lg(v);967if (l == 2) return gen_0;968if (l == 3) return gel(v,2);969return mkpolmod(v, T);970}971static GEN972bestapprnf_i(GEN x, GEN T, GEN V, long bit)973{974long i, l, tx = typ(x);975GEN z;976switch (tx)977{978case t_INT: case t_FRAC: return x;979case t_REAL: case t_COMPLEX: return cx_bestapprnf(x, T, V, bit);980case t_POLMOD: if (RgX_equal(gel(x,1),T)) return x;981break;982case t_POL: case t_SER: case t_VEC: case t_COL: case t_MAT:983l = lg(x); z = cgetg(l, tx);984for (i = 1; i < lontyp[tx]; i++) z[i] = x[i];985for (; i < l; i++) gel(z,i) = bestapprnf_i(gel(x,i), T, V, bit);986return z;987}988pari_err_TYPE("mfcxtoQ", x);989return NULL;/*LCOV_EXCL_LINE*/990}991992GEN993bestapprnf(GEN x, GEN T, GEN roT, long prec)994{995pari_sp av = avma;996long tx = typ(x), dT = 1, bit;997GEN V;998999if (T)1000{1001if (typ(T) != t_POL) T = nf_get_pol(checknf(T));1002else if (!RgX_is_ZX(T)) pari_err_TYPE("bestapprnf", T);1003dT = degpol(T);1004}1005if (is_rational_t(tx)) return gcopy(x);1006if (tx == t_POLMOD)1007{1008if (!T || !RgX_equal(T, gel(x,1))) pari_err_TYPE("bestapprnf",x);1009return gcopy(x);1010}10111012if (roT)1013{1014long l = gprecision(roT);1015switch(typ(roT))1016{1017case t_INT: case t_FRAC: case t_REAL: case t_COMPLEX: break;1018default: pari_err_TYPE("bestapprnf", roT);1019}1020if (prec < l) prec = l;1021}1022else if (!T)1023roT = gen_1;1024else1025{1026long n = poliscyclo(T); /* cyclotomic is an important special case */1027roT = n? rootsof1u_cx(n,prec): gel(QX_complex_roots(T,prec), 1);1028}1029V = vec_prepend(gpowers(roT, dT-1), NULL);1030bit = prec2nbits_mul(prec, 0.8);1031return gerepilecopy(av, bestapprnf_i(x, T, V, bit));1032}10331034/********************************************************************/1035/** **/1036/** MINIM **/1037/** **/1038/********************************************************************/1039void1040minim_alloc(long n, double ***q, GEN *x, double **y, double **z, double **v)1041{1042long i, s;10431044*x = cgetg(n, t_VECSMALL);1045*q = (double**) new_chunk(n);1046s = n * sizeof(double);1047*y = (double*) stack_malloc_align(s, sizeof(double));1048*z = (double*) stack_malloc_align(s, sizeof(double));1049*v = (double*) stack_malloc_align(s, sizeof(double));1050for (i=1; i<n; i++) (*q)[i] = (double*) stack_malloc_align(s, sizeof(double));1051}10521053static GEN1054ZC_canon(GEN V)1055{1056long l = lg(V), j;1057for (j = 1; j < l && signe(gel(V,j)) == 0; ++j);1058return (j < l && signe(gel(V,j)) < 0)? ZC_neg(V): V;1059}10601061static GEN1062ZM_zc_mul_canon(GEN u, GEN x)1063{1064return ZC_canon(ZM_zc_mul(u,x));1065}10661067static GEN1068ZM_zc_mul_canon_zm(GEN u, GEN x)1069{1070pari_sp av = avma;1071GEN M = ZV_to_zv(ZC_canon(ZM_zc_mul(u,x)));1072return gerepileupto(av, M);1073}10741075struct qfvec1076{1077GEN a, r, u;1078};10791080static void1081err_minim(GEN a)1082{1083pari_err_DOMAIN("minim0","form","is not",1084strtoGENstr("positive definite"),a);1085}10861087static GEN1088minim_lll(GEN a, GEN *u)1089{1090*u = lllgramint(a);1091if (lg(*u) != lg(a)) err_minim(a);1092return qf_apply_ZM(a,*u);1093}10941095static void1096forqfvec_init_dolll(struct qfvec *qv, GEN *pa, long dolll)1097{1098GEN r, u, a = *pa;1099if (!dolll) u = NULL;1100else1101{1102if (typ(a) != t_MAT || !RgM_is_ZM(a)) pari_err_TYPE("qfminim",a);1103a = *pa = minim_lll(a, &u);1104}1105qv->a = RgM_gtofp(a, DEFAULTPREC);1106r = qfgaussred_positive(qv->a);1107if (!r)1108{1109r = qfgaussred_positive(a); /* exact computation */1110if (!r) err_minim(a);1111r = RgM_gtofp(r, DEFAULTPREC);1112}1113qv->r = r;1114qv->u = u;1115}11161117static void1118forqfvec_init(struct qfvec *qv, GEN a)1119{ forqfvec_init_dolll(qv, &a, 1); }11201121static void1122forqfvec_i(void *E, long (*fun)(void *, GEN, GEN, double), struct qfvec *qv, GEN BORNE)1123{1124GEN x, a = qv->a, r = qv->r, u = qv->u;1125long n = lg(a), i, j, k;1126double p,BOUND,*v,*y,*z,**q;1127const double eps = 0.0001;1128if (!BORNE) BORNE = gen_0;1129else1130{1131BORNE = gfloor(BORNE);1132if (typ(BORNE) != t_INT) pari_err_TYPE("minim0",BORNE);1133if (signe(BORNE) <= 0) return;1134}1135if (n == 1) return;1136minim_alloc(n, &q, &x, &y, &z, &v);1137n--;1138for (j=1; j<=n; j++)1139{1140v[j] = rtodbl(gcoeff(r,j,j));1141for (i=1; i<j; i++) q[i][j] = rtodbl(gcoeff(r,i,j));1142}11431144if (gequal0(BORNE))1145{1146double c;1147p = rtodbl(gcoeff(a,1,1));1148for (i=2; i<=n; i++) { c = rtodbl(gcoeff(a,i,i)); if (c < p) p = c; }1149BORNE = roundr(dbltor(p));1150}1151else1152p = gtodouble(BORNE);1153BOUND = p * (1 + eps);1154if (BOUND == p) pari_err_PREC("minim0");11551156k = n; y[n] = z[n] = 0;1157x[n] = (long)sqrt(BOUND/v[n]);1158for(;;x[1]--)1159{1160do1161{1162if (k>1)1163{1164long l = k-1;1165z[l] = 0;1166for (j=k; j<=n; j++) z[l] += q[l][j]*x[j];1167p = (double)x[k] + z[k];1168y[l] = y[k] + p*p*v[k];1169x[l] = (long)floor(sqrt((BOUND-y[l])/v[l])-z[l]);1170k = l;1171}1172for(;;)1173{1174p = (double)x[k] + z[k];1175if (y[k] + p*p*v[k] <= BOUND) break;1176k++; x[k]--;1177}1178} while (k > 1);1179if (! x[1] && y[1]<=eps) break;11801181p = (double)x[1] + z[1]; p = y[1] + p*p*v[1]; /* norm(x) */1182if (fun(E, u, x, p)) break;1183}1184}11851186void1187forqfvec(void *E, long (*fun)(void *, GEN, GEN, double), GEN a, GEN BORNE)1188{1189pari_sp av = avma;1190struct qfvec qv;1191forqfvec_init(&qv, a);1192forqfvec_i(E, fun, &qv, BORNE);1193set_avma(av);1194}11951196struct qfvecwrap1197{1198void *E;1199long (*fun)(void *, GEN);1200};12011202static long1203forqfvec_wrap(void *E, GEN u, GEN x, double d)1204{1205pari_sp av = avma;1206struct qfvecwrap *W = (struct qfvecwrap *) E;1207(void) d;1208return gc_long(av, W->fun(W->E, ZM_zc_mul_canon(u, x)));1209}12101211void1212forqfvec1(void *E, long (*fun)(void *, GEN), GEN a, GEN BORNE)1213{1214pari_sp av = avma;1215struct qfvecwrap wr;1216struct qfvec qv;1217wr.E = E; wr.fun = fun;1218forqfvec_init(&qv, a);1219forqfvec_i((void*) &wr, forqfvec_wrap, &qv, BORNE);1220set_avma(av);1221}12221223void1224forqfvec0(GEN a, GEN BORNE, GEN code)1225{ EXPRVOID_WRAP(code, forqfvec1(EXPR_ARGVOID, a, BORNE)) }12261227enum { min_ALL = 0, min_FIRST, min_VECSMALL, min_VECSMALL2 };12281229/* Minimal vectors for the integral definite quadratic form: a.1230* Result u:1231* u[1]= Number of vectors of square norm <= BORNE1232* u[2]= maximum norm found1233* u[3]= list of vectors found (at most STOCKMAX, unless NULL)1234*1235* If BORNE = NULL: Minimal nonzero vectors.1236* flag = min_ALL, as above1237* flag = min_FIRST, exits when first suitable vector is found.1238* flag = min_VECSMALL, return a t_VECSMALL of (half) the number of vectors1239* for each norm1240* flag = min_VECSMALL2, same but count only vectors with even norm, and shift1241* the answer */1242static GEN1243minim0_dolll(GEN a, GEN BORNE, GEN STOCKMAX, long flag, long dolll)1244{1245GEN x, u, r, L, gnorme;1246long n = lg(a), i, j, k, s, maxrank, sBORNE;1247pari_sp av = avma, av1;1248double p,maxnorm,BOUND,*v,*y,*z,**q;1249const double eps = 1e-10;1250int stockall = 0;1251struct qfvec qv;12521253if (!BORNE)1254sBORNE = 0;1255else1256{1257BORNE = gfloor(BORNE);1258if (typ(BORNE) != t_INT) pari_err_TYPE("minim0",BORNE);1259if (is_bigint(BORNE)) pari_err_PREC( "qfminim");1260sBORNE = itos(BORNE); set_avma(av);1261if (sBORNE < 0) sBORNE = 0;1262}1263if (!STOCKMAX)1264{1265stockall = 1;1266maxrank = 200;1267}1268else1269{1270STOCKMAX = gfloor(STOCKMAX);1271if (typ(STOCKMAX) != t_INT) pari_err_TYPE("minim0",STOCKMAX);1272maxrank = itos(STOCKMAX);1273if (maxrank < 0)1274pari_err_TYPE("minim0 [negative number of vectors]",STOCKMAX);1275}12761277switch(flag)1278{1279case min_VECSMALL:1280case min_VECSMALL2:1281if (sBORNE <= 0) return cgetg(1, t_VECSMALL);1282L = zero_zv(sBORNE);1283if (flag == min_VECSMALL2) sBORNE <<= 1;1284if (n == 1) return L;1285break;1286case min_FIRST:1287if (n == 1 || (!sBORNE && BORNE)) return cgetg(1,t_VEC);1288L = NULL; /* gcc -Wall */1289break;1290case min_ALL:1291if (n == 1 || (!sBORNE && BORNE))1292retmkvec3(gen_0, gen_0, cgetg(1, t_MAT));1293L = new_chunk(1+maxrank);1294break;1295default:1296return NULL;1297}1298minim_alloc(n, &q, &x, &y, &z, &v);12991300forqfvec_init_dolll(&qv, &a, dolll);1301av1 = avma;1302r = qv.r;1303u = qv.u;1304n--;1305for (j=1; j<=n; j++)1306{1307v[j] = rtodbl(gcoeff(r,j,j));1308for (i=1; i<j; i++) q[i][j] = rtodbl(gcoeff(r,i,j)); /* |.| <= 1/2 */1309}13101311if (sBORNE) maxnorm = 0.;1312else1313{1314GEN B = gcoeff(a,1,1);1315long t = 1;1316for (i=2; i<=n; i++)1317{1318GEN c = gcoeff(a,i,i);1319if (cmpii(c, B) < 0) { B = c; t = i; }1320}1321if (flag == min_FIRST) return gerepilecopy(av, mkvec2(B, gel(u,t)));1322maxnorm = -1.; /* don't update maxnorm */1323if (is_bigint(B)) return NULL;1324sBORNE = itos(B);1325}1326BOUND = sBORNE * (1 + eps);1327if ((long)BOUND != sBORNE) return NULL;13281329s = 0;1330k = n; y[n] = z[n] = 0;1331x[n] = (long)sqrt(BOUND/v[n]);1332for(;;x[1]--)1333{1334do1335{1336if (k>1)1337{1338long l = k-1;1339z[l] = 0;1340for (j=k; j<=n; j++) z[l] += q[l][j]*x[j];1341p = (double)x[k] + z[k];1342y[l] = y[k] + p*p*v[k];1343x[l] = (long)floor(sqrt((BOUND-y[l])/v[l])-z[l]);1344k = l;1345}1346for(;;)1347{1348p = (double)x[k] + z[k];1349if (y[k] + p*p*v[k] <= BOUND) break;1350k++; x[k]--;1351}1352}1353while (k > 1);1354if (! x[1] && y[1]<=eps) break;13551356p = (double)x[1] + z[1]; p = y[1] + p*p*v[1]; /* norm(x) */1357if (maxnorm >= 0)1358{1359if (p > maxnorm) maxnorm = p;1360}1361else1362{ /* maxnorm < 0 : only look for minimal vectors */1363pari_sp av2 = avma;1364gnorme = roundr(dbltor(p));1365if (cmpis(gnorme, sBORNE) >= 0) set_avma(av2);1366else1367{1368sBORNE = itos(gnorme); set_avma(av1);1369BOUND = sBORNE * (1+eps);1370L = new_chunk(maxrank+1);1371s = 0;1372}1373}1374s++;13751376switch(flag)1377{1378case min_FIRST:1379if (dolll) x = ZM_zc_mul_canon(u,x);1380return gerepilecopy(av, mkvec2(roundr(dbltor(p)), x));13811382case min_ALL:1383if (s > maxrank && stockall) /* overflow */1384{1385long maxranknew = maxrank << 1;1386GEN Lnew = new_chunk(1 + maxranknew);1387for (i=1; i<=maxrank; i++) Lnew[i] = L[i];1388L = Lnew; maxrank = maxranknew;1389}1390if (s<=maxrank) gel(L,s) = leafcopy(x);1391break;13921393case min_VECSMALL:1394{ ulong norm = (ulong)(p + 0.5); L[norm]++; }1395break;13961397case min_VECSMALL2:1398{ ulong norm = (ulong)(p + 0.5); if (!odd(norm)) L[norm>>1]++; }1399break;14001401}1402}1403switch(flag)1404{1405case min_FIRST:1406set_avma(av); return cgetg(1,t_VEC);1407case min_VECSMALL:1408case min_VECSMALL2:1409set_avma((pari_sp)L); return L;1410}1411r = (maxnorm >= 0) ? roundr(dbltor(maxnorm)): stoi(sBORNE);1412k = minss(s,maxrank);1413L[0] = evaltyp(t_MAT) | evallg(k + 1);1414if (dolll)1415for (j=1; j<=k; j++)1416gel(L,j) = dolll==1 ? ZM_zc_mul_canon(u, gel(L,j))1417: ZM_zc_mul_canon_zm(u, gel(L,j));1418return gerepilecopy(av, mkvec3(stoi(s<<1), r, L));1419}14201421static GEN1422minim0(GEN a, GEN BORNE, GEN STOCKMAX, long flag)1423{1424GEN v = minim0_dolll(a, BORNE, STOCKMAX, flag, 1);1425if (!v) pari_err_PREC("qfminim");1426return v;1427}14281429static GEN1430minim0_zm(GEN a, GEN BORNE, GEN STOCKMAX, long flag)1431{1432GEN v = minim0_dolll(a, BORNE, STOCKMAX, flag, 2);1433if (!v) pari_err_PREC("qfminim");1434return v;1435}14361437GEN1438qfrep0(GEN a, GEN borne, long flag)1439{ return minim0(a, borne, gen_0, (flag & 1)? min_VECSMALL2: min_VECSMALL); }14401441GEN1442qfminim0(GEN a, GEN borne, GEN stockmax, long flag, long prec)1443{1444switch(flag)1445{1446case 0: return minim0(a,borne,stockmax,min_ALL);1447case 1: return minim0(a,borne,gen_0 ,min_FIRST);1448case 2:1449{1450long maxnum = -1;1451if (typ(a) != t_MAT) pari_err_TYPE("qfminim",a);1452if (stockmax) {1453if (typ(stockmax) != t_INT) pari_err_TYPE("qfminim",stockmax);1454maxnum = itos(stockmax);1455}1456a = fincke_pohst(a,borne,maxnum,prec,NULL);1457if (!a) pari_err_PREC("qfminim");1458return a;1459}1460default: pari_err_FLAG("qfminim");1461}1462return NULL; /* LCOV_EXCL_LINE */1463}14641465GEN1466minim(GEN a, GEN borne, GEN stockmax)1467{ return minim0(a,borne,stockmax,min_ALL); }14681469GEN1470minim_zm(GEN a, GEN borne, GEN stockmax)1471{ return minim0_zm(a,borne,stockmax,min_ALL); }14721473GEN1474minim_raw(GEN a, GEN BORNE, GEN STOCKMAX)1475{ return minim0_dolll(a, BORNE, STOCKMAX, min_ALL, 0); }14761477GEN1478minim2(GEN a, GEN borne, GEN stockmax)1479{ return minim0(a,borne,stockmax,min_FIRST); }14801481/* If V depends linearly from the columns of the matrix, return 0.1482* Otherwise, update INVP and L and return 1. No GC. */1483static int1484addcolumntomatrix(GEN V, GEN invp, GEN L)1485{1486long i,j,k, n = lg(invp);1487GEN a = cgetg(n, t_COL), ak = NULL, mak;14881489for (k = 1; k < n; k++)1490if (!L[k])1491{1492ak = RgMrow_zc_mul(invp, V, k);1493if (!gequal0(ak)) break;1494}1495if (k == n) return 0;1496L[k] = 1;1497mak = gneg_i(ak);1498for (i=k+1; i<n; i++)1499gel(a,i) = gdiv(RgMrow_zc_mul(invp, V, i), mak);1500for (j=1; j<=k; j++)1501{1502GEN c = gel(invp,j), ck = gel(c,k);1503if (gequal0(ck)) continue;1504gel(c,k) = gdiv(ck, ak);1505if (j==k)1506for (i=k+1; i<n; i++)1507gel(c,i) = gmul(gel(a,i), ck);1508else1509for (i=k+1; i<n; i++)1510gel(c,i) = gadd(gel(c,i), gmul(gel(a,i), ck));1511}1512return 1;1513}15141515GEN1516qfperfection(GEN a)1517{1518pari_sp av = avma;1519GEN u, L;1520long r, s, k, l, n = lg(a)-1;15211522if (!n) return gen_0;1523if (typ(a) != t_MAT || !RgM_is_ZM(a)) pari_err_TYPE("qfperfection",a);1524a = minim_lll(a, &u);1525L = minim_raw(a,NULL,NULL);1526r = (n*(n+1)) >> 1;1527if (L)1528{1529GEN D, V, invp;1530L = gel(L, 3); l = lg(L);1531if (l == 2) { set_avma(av); return gen_1; }1532/* |L[i]|^2 fits into a long for all i */1533D = zero_zv(r);1534V = cgetg(r+1, t_VECSMALL);1535invp = matid(r);1536s = 0;1537for (k = 1; k < l; k++)1538{1539pari_sp av2 = avma;1540GEN x = gel(L,k);1541long i, j, I;1542for (i = I = 1; i<=n; i++)1543for (j=i; j<=n; j++,I++) V[I] = x[i]*x[j];1544if (!addcolumntomatrix(V,invp,D)) set_avma(av2);1545else if (++s == r) break;1546}1547}1548else1549{1550GEN M;1551L = fincke_pohst(a,NULL,-1, DEFAULTPREC, NULL);1552if (!L) pari_err_PREC("qfminim");1553L = gel(L, 3); l = lg(L);1554if (l == 2) { set_avma(av); return gen_1; }1555M = cgetg(l, t_MAT);1556for (k = 1; k < l; k++)1557{1558GEN x = gel(L,k), c = cgetg(r+1, t_COL);1559long i, I, j;1560for (i = I = 1; i<=n; i++)1561for (j=i; j<=n; j++,I++) gel(c,I) = mulii(gel(x,i), gel(x,j));1562gel(M,k) = c;1563}1564s = ZM_rank(M);1565}1566set_avma(av); return utoipos(s);1567}15681569static GEN1570clonefill(GEN S, long s, long t)1571{ /* initialize to dummy values */1572GEN T = S, dummy = cgetg(1, t_STR);1573long i;1574for (i = s+1; i <= t; i++) gel(S,i) = dummy;1575S = gclone(S); if (isclone(T)) gunclone(T);1576return S;1577}15781579/* increment ZV x, by incrementing cell of index k. Initial value x0[k] was1580* chosen to minimize qf(x) for given x0[1..k-1] and x0[k+1,..] = 01581* The last nonzero entry must be positive and goes through x0[k]+1,2,3,...1582* Others entries go through: x0[k]+1,-1,2,-2,...*/1583INLINE void1584step(GEN x, GEN y, GEN inc, long k)1585{1586if (!signe(gel(y,k))) /* x[k+1..] = 0 */1587gel(x,k) = addiu(gel(x,k), 1); /* leading coeff > 0 */1588else1589{1590long i = inc[k];1591gel(x,k) = addis(gel(x,k), i),1592inc[k] = (i > 0)? -1-i: 1-i;1593}1594}15951596/* 1 if we are "sure" that x < y, up to few rounding errors, i.e.1597* x < y - epsilon. More precisely :1598* - sign(x - y) < 01599* - lgprec(x-y) > 3 || expo(x - y) - expo(x) > -24 */1600static int1601mplessthan(GEN x, GEN y)1602{1603pari_sp av = avma;1604GEN z = mpsub(x, y);1605set_avma(av);1606if (typ(z) == t_INT) return (signe(z) < 0);1607if (signe(z) >= 0) return 0;1608if (realprec(z) > LOWDEFAULTPREC) return 1;1609return ( expo(z) - mpexpo(x) > -24 );1610}16111612/* 1 if we are "sure" that x > y, up to few rounding errors, i.e.1613* x > y + epsilon */1614static int1615mpgreaterthan(GEN x, GEN y)1616{1617pari_sp av = avma;1618GEN z = mpsub(x, y);1619set_avma(av);1620if (typ(z) == t_INT) return (signe(z) > 0);1621if (signe(z) <= 0) return 0;1622if (realprec(z) > LOWDEFAULTPREC) return 1;1623return ( expo(z) - mpexpo(x) > -24 );1624}16251626/* x a t_INT, y t_INT or t_REAL */1627INLINE GEN1628mulimp(GEN x, GEN y)1629{1630if (typ(y) == t_INT) return mulii(x,y);1631return signe(x) ? mulir(x,y): gen_0;1632}1633/* x + y*z, x,z two mp's, y a t_INT */1634INLINE GEN1635addmulimp(GEN x, GEN y, GEN z)1636{1637if (!signe(y)) return x;1638if (typ(z) == t_INT) return mpadd(x, mulii(y, z));1639return mpadd(x, mulir(y, z));1640}16411642/* yk + vk * (xk + zk)^2 */1643static GEN1644norm_aux(GEN xk, GEN yk, GEN zk, GEN vk)1645{1646GEN t = mpadd(xk, zk);1647if (typ(t) == t_INT) { /* probably gen_0, avoid loss of accuracy */1648yk = addmulimp(yk, sqri(t), vk);1649} else {1650yk = mpadd(yk, mpmul(sqrr(t), vk));1651}1652return yk;1653}1654/* yk + vk * (xk + zk)^2 < B + epsilon */1655static int1656check_bound(GEN B, GEN xk, GEN yk, GEN zk, GEN vk)1657{1658pari_sp av = avma;1659int f = mpgreaterthan(norm_aux(xk,yk,zk,vk), B);1660return gc_bool(av, !f);1661}16621663/* q(k-th canonical basis vector), where q is given in Cholesky form1664* q(x) = sum_{i = 1}^n q[i,i] (x[i] + sum_{j > i} q[i,j] x[j])^2.1665* Namely q(e_k) = q[k,k] + sum_{i < k} q[i,i] q[i,k]^21666* Assume 1 <= k <= n. */1667static GEN1668cholesky_norm_ek(GEN q, long k)1669{1670GEN t = gcoeff(q,k,k);1671long i;1672for (i = 1; i < k; i++) t = norm_aux(gen_0, t, gcoeff(q,i,k), gcoeff(q,i,i));1673return t;1674}16751676/* q is the Cholesky decomposition of a quadratic form1677* Enumerate vectors whose norm is less than BORNE (Algo 2.5.7),1678* minimal vectors if BORNE = NULL (implies check = NULL).1679* If (check != NULL) consider only vectors passing the check, and assumes1680* we only want the smallest possible vectors */1681static GEN1682smallvectors(GEN q, GEN BORNE, long maxnum, FP_chk_fun *CHECK)1683{1684long N = lg(q), n = N-1, i, j, k, s, stockmax, checkcnt = 1;1685pari_sp av, av1;1686GEN inc, S, x, y, z, v, p1, alpha, norms;1687GEN norme1, normax1, borne1, borne2;1688GEN (*check)(void *,GEN) = CHECK? CHECK->f: NULL;1689void *data = CHECK? CHECK->data: NULL;1690const long skipfirst = CHECK? CHECK->skipfirst: 0;1691const int stockall = (maxnum == -1);16921693alpha = dbltor(0.95);1694normax1 = gen_0;16951696v = cgetg(N,t_VEC);1697inc = const_vecsmall(n, 1);16981699av = avma;1700stockmax = stockall? 2000: maxnum;1701norms = cgetg(check?(stockmax+1): 1,t_VEC); /* unused if (!check) */1702S = cgetg(stockmax+1,t_VEC);1703x = cgetg(N,t_COL);1704y = cgetg(N,t_COL);1705z = cgetg(N,t_COL);1706for (i=1; i<N; i++) {1707gel(v,i) = gcoeff(q,i,i);1708gel(x,i) = gel(y,i) = gel(z,i) = gen_0;1709}1710if (BORNE)1711{1712borne1 = BORNE;1713if (gsigne(borne1) <= 0) retmkvec3(gen_0, gen_0, cgetg(1,t_MAT));1714if (typ(borne1) != t_REAL)1715{1716long prec;1717prec = nbits2prec(gexpo(borne1) + 10);1718borne1 = gtofp(borne1, maxss(prec, DEFAULTPREC));1719}1720}1721else1722{1723borne1 = gcoeff(q,1,1);1724for (i=2; i<N; i++)1725{1726GEN b = cholesky_norm_ek(q, i);1727if (gcmp(b, borne1) < 0) borne1 = b;1728}1729/* borne1 = norm of smallest basis vector */1730}1731borne2 = mulrr(borne1,alpha);1732if (DEBUGLEVEL>2)1733err_printf("smallvectors looking for norm < %P.4G\n",borne1);1734s = 0; k = n;1735for(;; step(x,y,inc,k)) /* main */1736{ /* x (supposedly) small vector, ZV.1737* For all t >= k, we have1738* z[t] = sum_{j > t} q[t,j] * x[j]1739* y[t] = sum_{i > t} q[i,i] * (x[i] + z[i])^21740* = 0 <=> x[i]=0 for all i>t */1741do1742{1743int skip = 0;1744if (k > 1)1745{1746long l = k-1;1747av1 = avma;1748p1 = mulimp(gel(x,k), gcoeff(q,l,k));1749for (j=k+1; j<N; j++) p1 = addmulimp(p1, gel(x,j), gcoeff(q,l,j));1750gel(z,l) = gerepileuptoleaf(av1,p1);17511752av1 = avma;1753p1 = norm_aux(gel(x,k), gel(y,k), gel(z,k), gel(v,k));1754gel(y,l) = gerepileuptoleaf(av1, p1);1755/* skip the [x_1,...,x_skipfirst,0,...,0] */1756if ((l <= skipfirst && !signe(gel(y,skipfirst)))1757|| mplessthan(borne1, gel(y,l))) skip = 1;1758else /* initial value, minimizing (x[l] + z[l])^2, hence qf(x) for1759the given x[1..l-1] */1760gel(x,l) = mpround( mpneg(gel(z,l)) );1761k = l;1762}1763for(;; step(x,y,inc,k))1764{ /* at most 2n loops */1765if (!skip)1766{1767if (check_bound(borne1, gel(x,k),gel(y,k),gel(z,k),gel(v,k))) break;1768step(x,y,inc,k);1769if (check_bound(borne1, gel(x,k),gel(y,k),gel(z,k),gel(v,k))) break;1770}1771skip = 0; inc[k] = 1;1772if (++k > n) goto END;1773}17741775if (gc_needed(av,2))1776{1777if(DEBUGMEM>1) pari_warn(warnmem,"smallvectors");1778if (stockmax) S = clonefill(S, s, stockmax);1779if (check) {1780GEN dummy = cgetg(1, t_STR);1781for (i=s+1; i<=stockmax; i++) gel(norms,i) = dummy;1782}1783gerepileall(av,7,&x,&y,&z,&normax1,&borne1,&borne2,&norms);1784}1785}1786while (k > 1);1787if (!signe(gel(x,1)) && !signe(gel(y,1))) continue; /* exclude 0 */17881789av1 = avma;1790norme1 = norm_aux(gel(x,1),gel(y,1),gel(z,1),gel(v,1));1791if (mpgreaterthan(norme1,borne1)) { set_avma(av1); continue; /* main */ }17921793norme1 = gerepileuptoleaf(av1,norme1);1794if (check)1795{1796if (checkcnt < 5 && mpcmp(norme1, borne2) < 0)1797{1798if (!check(data,x)) { checkcnt++ ; continue; /* main */}1799if (DEBUGLEVEL>4) err_printf("New bound: %Ps", norme1);1800borne1 = norme1;1801borne2 = mulrr(borne1, alpha);1802s = 0; checkcnt = 0;1803}1804}1805else1806{1807if (!BORNE) /* find minimal vectors */1808{1809if (mplessthan(norme1, borne1))1810{ /* strictly smaller vector than previously known */1811borne1 = norme1; /* + epsilon */1812s = 0;1813}1814}1815else1816if (mpcmp(norme1,normax1) > 0) normax1 = norme1;1817}1818if (++s > stockmax) continue; /* too many vectors: no longer remember */1819if (check) gel(norms,s) = norme1;1820gel(S,s) = leafcopy(x);1821if (s != stockmax) continue; /* still room, get next vector */18221823if (check)1824{ /* overflow, eliminate vectors failing "check" */1825pari_sp av2 = avma;1826long imin, imax;1827GEN per = indexsort(norms), S2 = cgetg(stockmax+1, t_VEC);1828if (DEBUGLEVEL>2) err_printf("sorting... [%ld elts]\n",s);1829/* let N be the minimal norm so far for x satisfying 'check'. Keep1830* all elements of norm N */1831for (i = 1; i <= s; i++)1832{1833long k = per[i];1834if (check(data,gel(S,k))) { borne1 = gel(norms,k); break; }1835}1836imin = i;1837for (; i <= s; i++)1838if (mpgreaterthan(gel(norms,per[i]), borne1)) break;1839imax = i;1840for (i=imin, s=0; i < imax; i++) gel(S2,++s) = gel(S,per[i]);1841for (i = 1; i <= s; i++) gel(S,i) = gel(S2,i);1842set_avma(av2);1843if (s) { borne2 = mulrr(borne1, alpha); checkcnt = 0; }1844if (!stockall) continue;1845if (s > stockmax/2) stockmax <<= 1;1846norms = cgetg(stockmax+1, t_VEC);1847for (i = 1; i <= s; i++) gel(norms,i) = borne1;1848}1849else1850{1851if (!stockall && BORNE) goto END;1852if (!stockall) continue;1853stockmax <<= 1;1854}18551856{1857GEN Snew = clonefill(vec_lengthen(S,stockmax), s, stockmax);1858if (isclone(S)) gunclone(S);1859S = Snew;1860}1861}1862END:1863if (s < stockmax) stockmax = s;1864if (check)1865{1866GEN per, alph, pols, p;1867if (DEBUGLEVEL>2) err_printf("final sort & check...\n");1868setlg(norms,stockmax+1); per = indexsort(norms);1869alph = cgetg(stockmax+1,t_VEC);1870pols = cgetg(stockmax+1,t_VEC);1871for (j=0,i=1; i<=stockmax; i++)1872{1873long t = per[i];1874GEN N = gel(norms,t);1875if (j && mpgreaterthan(N, borne1)) break;1876if ((p = check(data,gel(S,t))))1877{1878if (!j) borne1 = N;1879j++;1880gel(pols,j) = p;1881gel(alph,j) = gel(S,t);1882}1883}1884setlg(pols,j+1);1885setlg(alph,j+1);1886if (stockmax && isclone(S)) { alph = gcopy(alph); gunclone(S); }1887return mkvec2(pols, alph);1888}1889if (stockmax)1890{1891setlg(S,stockmax+1);1892settyp(S,t_MAT);1893if (isclone(S)) { p1 = S; S = gcopy(S); gunclone(p1); }1894}1895else1896S = cgetg(1,t_MAT);1897return mkvec3(utoi(s<<1), borne1, S);1898}18991900/* solve q(x) = x~.a.x <= bound, a > 0.1901* If check is non-NULL keep x only if check(x).1902* If a is a vector, assume a[1] is the LLL-reduced Cholesky form of q */1903GEN1904fincke_pohst(GEN a, GEN B0, long stockmax, long PREC, FP_chk_fun *CHECK)1905{1906pari_sp av = avma;1907VOLATILE long i,j,l;1908VOLATILE GEN r,rinv,rinvtrans,u,v,res,z,vnorm,rperm,perm,uperm, bound = B0;19091910if (typ(a) == t_VEC)1911{1912r = gel(a,1);1913u = NULL;1914}1915else1916{1917long prec = PREC;1918l = lg(a);1919if (l == 1)1920{1921if (CHECK) pari_err_TYPE("fincke_pohst [dimension 0]", a);1922retmkvec3(gen_0, gen_0, cgetg(1,t_MAT));1923}1924u = lllfp(a, 0.75, LLL_GRAM | LLL_IM);1925if (lg(u) != lg(a)) return NULL;1926r = qf_apply_RgM(a,u);1927i = gprecision(r);1928if (i)1929prec = i;1930else {1931prec = DEFAULTPREC + nbits2extraprec(gexpo(r));1932if (prec < PREC) prec = PREC;1933}1934if (DEBUGLEVEL>2) err_printf("first LLL: prec = %ld\n", prec);1935r = qfgaussred_positive(r);1936if (!r) return NULL;1937for (i=1; i<l; i++)1938{1939GEN s = gsqrt(gcoeff(r,i,i), prec);1940gcoeff(r,i,i) = s;1941for (j=i+1; j<l; j++) gcoeff(r,i,j) = gmul(s, gcoeff(r,i,j));1942}1943}1944/* now r~ * r = a in LLL basis */1945rinv = RgM_inv_upper(r);1946if (!rinv) return NULL;1947rinvtrans = shallowtrans(rinv);1948if (DEBUGLEVEL>2)1949err_printf("Fincke-Pohst, final LLL: prec = %ld\n", gprecision(rinvtrans));1950v = lll(rinvtrans);1951if (lg(v) != lg(rinvtrans)) return NULL;19521953rinvtrans = RgM_mul(rinvtrans, v);1954v = ZM_inv(shallowtrans(v),NULL);1955r = RgM_mul(r,v);1956u = u? ZM_mul(u,v): v;19571958l = lg(r);1959vnorm = cgetg(l,t_VEC);1960for (j=1; j<l; j++) gel(vnorm,j) = gnorml2(gel(rinvtrans,j));1961rperm = cgetg(l,t_MAT);1962uperm = cgetg(l,t_MAT); perm = indexsort(vnorm);1963for (i=1; i<l; i++) { uperm[l-i] = u[perm[i]]; rperm[l-i] = r[perm[i]]; }1964u = uperm;1965r = rperm; res = NULL;1966pari_CATCH(e_PREC) { }1967pari_TRY {1968GEN q;1969if (CHECK && CHECK->f_init) bound = CHECK->f_init(CHECK, r, u);1970q = gaussred_from_QR(r, gprecision(vnorm));1971if (!q) pari_err_PREC("fincke_pohst");1972res = smallvectors(q, bound, stockmax, CHECK);1973} pari_ENDCATCH;1974if (CHECK)1975{1976if (CHECK->f_post) res = CHECK->f_post(CHECK, res, u);1977return res;1978}1979if (!res) pari_err_PREC("fincke_pohst");19801981z = cgetg(4,t_VEC);1982gel(z,1) = gcopy(gel(res,1));1983gel(z,2) = gcopy(gel(res,2));1984gel(z,3) = ZM_mul(u, gel(res,3)); return gerepileupto(av,z);1985}198619871988