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; either version 2 of the License, or (at your option) any later
8 : version. It is distributed in the hope that it will be useful, but WITHOUT
9 : ANY WARRANTY WHATSOEVER.
10 :
11 : Check the License for details. You should have received a copy of it, along
12 : with the package; see the file 'COPYING'. If not, write to the Free Software
13 : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
14 :
15 : /********************************************************************/
16 : /** **/
17 : /** LLL Algorithm and close friends **/
18 : /** **/
19 : /********************************************************************/
20 : #include "pari.h"
21 : #include "paripriv.h"
22 :
23 : #define DEBUGLEVEL DEBUGLEVEL_qf
24 :
25 : /********************************************************************/
26 : /** QR Factorization via Householder matrices **/
27 : /********************************************************************/
28 : static int
29 24725719 : no_prec_pb(GEN x)
30 : {
31 24650520 : return (typ(x) != t_REAL || realprec(x) > DEFAULTPREC
32 49376239 : || expo(x) < DEFAULTPREC>>1);
33 : }
34 : /* Find a Householder transformation which, applied to x[k..#x], zeroes
35 : * x[k+1..#x]; fill L = (mu_{i,j}). Return 0 if precision problem [obtained
36 : * a 0 vector], 1 otherwise */
37 : static int
38 24734539 : FindApplyQ(GEN x, GEN L, GEN B, long k, GEN Q, long prec)
39 : {
40 24734539 : long i, nx = lg(x)-1;
41 24734539 : GEN x2, x1, xd = x + (k-1);
42 :
43 24734539 : x1 = gel(xd,1);
44 24734539 : x2 = mpsqr(x1);
45 24733577 : if (k < nx)
46 : {
47 19493533 : long lv = nx - (k-1) + 1;
48 19493533 : GEN beta, Nx, v = cgetg(lv, t_VEC);
49 76999217 : for (i=2; i<lv; i++) {
50 57505988 : x2 = mpadd(x2, mpsqr(gel(xd,i)));
51 57505368 : gel(v,i) = gel(xd,i);
52 : }
53 19493229 : if (!signe(x2)) return 0;
54 19485009 : Nx = gsqrt(x2, prec); if (signe(x1) < 0) setsigne(Nx, -1);
55 19485969 : gel(v,1) = mpadd(x1, Nx);
56 :
57 19485244 : if (!signe(x1))
58 732924 : beta = gtofp(x2, prec); /* make sure typ(beta) != t_INT */
59 : else
60 18752320 : beta = mpadd(x2, mpmul(Nx,x1));
61 19485464 : gel(Q,k) = mkvec2(invr(beta), v); /* [t_REAL, vector of t_INT/t_REALs] */
62 :
63 19485731 : togglesign(Nx);
64 19485501 : gcoeff(L,k,k) = Nx; /* t_REAL */
65 : }
66 : else /* k = nx */
67 5240044 : gcoeff(L,k,k) = gel(x,k); /* t_INT or t_REAL */
68 24725545 : gel(B,k) = x2;
69 70820039 : for (i=1; i<k; i++) gcoeff(L,k,i) = gel(x,i); /* t_INT or t_REAL */
70 24725545 : return no_prec_pb(x2);
71 : }
72 :
73 : /* apply Householder transformation Q = [beta,v] to r with t_INT/t_REAL
74 : * coefficients, in place: r -= ((0|v).r * beta) v */
75 : static void
76 46104272 : ApplyQ(GEN Q, GEN r)
77 : {
78 46104272 : GEN s, rd, beta = gel(Q,1), v = gel(Q,2);
79 46104272 : long i, l = lg(v), lr = lg(r);
80 :
81 46104272 : rd = r + (lr - l);
82 46104272 : s = mpmul(gel(v,1), gel(rd,1));
83 478873108 : for (i=2; i<l; i++) s = mpadd(s, mpmul(gel(v,i), gel(rd,i)));
84 46100666 : s = mpmul(beta, s);
85 525159505 : for (i=1; i<l; i++)
86 479053288 : if (signe(gel(v,i))) gel(rd,i) = mpsub(gel(rd,i), mpmul(s, gel(v,i)));
87 46106217 : }
88 : /* apply Q[1], ..., Q[j-1] to r */
89 : static GEN
90 16980584 : ApplyAllQ(GEN Q, GEN r, long j)
91 : {
92 16980584 : pari_sp av = avma;
93 : long i;
94 16980584 : r = leafcopy(r);
95 63082727 : for (i=1; i<j; i++) ApplyQ(gel(Q,i), r);
96 16978649 : return gc_GEN(av, r);
97 : }
98 :
99 : /* same, arbitrary coefficients [20% slower for t_REAL at DEFAULTPREC] */
100 : static void
101 22113 : RgC_ApplyQ(GEN Q, GEN r)
102 : {
103 22113 : GEN s, rd, beta = gel(Q,1), v = gel(Q,2);
104 22113 : long i, l = lg(v), lr = lg(r);
105 :
106 22113 : rd = r + (lr - l);
107 22113 : s = gmul(gel(v,1), gel(rd,1));
108 464373 : for (i=2; i<l; i++) s = gadd(s, gmul(gel(v,i), gel(rd,i)));
109 22113 : s = gmul(beta, s);
110 486486 : for (i=1; i<l; i++)
111 464373 : if (signe(gel(v,i))) gel(rd,i) = gsub(gel(rd,i), gmul(s, gel(v,i)));
112 22113 : }
113 : static GEN
114 567 : RgC_ApplyAllQ(GEN Q, GEN r, long j)
115 : {
116 567 : pari_sp av = avma;
117 : long i;
118 567 : r = leafcopy(r);
119 22680 : for (i=1; i<j; i++) RgC_ApplyQ(gel(Q,i), r);
120 567 : return gc_GEN(av, r);
121 : }
122 :
123 : int
124 21 : RgM_QR_init(GEN x, GEN *pB, GEN *pQ, GEN *pL, long prec)
125 : {
126 21 : x = RgM_gtomp(x, prec);
127 21 : return QR_init(x, pB, pQ, pL, prec);
128 : }
129 :
130 : static void
131 35 : check_householder(GEN Q)
132 : {
133 35 : long i, l = lg(Q);
134 35 : if (typ(Q) != t_VEC) pari_err_TYPE("mathouseholder", Q);
135 854 : for (i = 1; i < l; i++)
136 : {
137 826 : GEN q = gel(Q,i), v;
138 826 : if (typ(q) != t_VEC || lg(q) != 3) pari_err_TYPE("mathouseholder", Q);
139 826 : v = gel(q,2);
140 826 : if (typ(v) != t_VEC || lg(v)+i-2 != l) pari_err_TYPE("mathouseholder", Q);
141 : }
142 28 : }
143 :
144 : GEN
145 35 : mathouseholder(GEN Q, GEN x)
146 : {
147 35 : long l = lg(Q);
148 35 : check_householder(Q);
149 28 : switch(typ(x))
150 : {
151 14 : case t_MAT:
152 14 : if (lg(x) == 1) return cgetg(1, t_MAT);
153 14 : if (lgcols(x) != l+1) pari_err_TYPE("mathouseholder", x);
154 574 : pari_APPLY_same(RgC_ApplyAllQ(Q, gel(x,i), l));
155 7 : case t_COL:
156 7 : if (lg(x) == l+1) return RgC_ApplyAllQ(Q, x, l);
157 : }
158 7 : pari_err_TYPE("mathouseholder", x);
159 : return NULL; /* LCOV_EXCL_LINE */
160 : }
161 :
162 : GEN
163 35 : matqr(GEN x, long flag, long prec)
164 : {
165 35 : pari_sp av = avma;
166 : GEN B, Q, L;
167 35 : long n = lg(x)-1;
168 35 : if (typ(x) != t_MAT) pari_err_TYPE("matqr",x);
169 35 : if (!n)
170 : {
171 14 : if (!flag) retmkvec2(cgetg(1,t_MAT),cgetg(1,t_MAT));
172 7 : retmkvec2(cgetg(1,t_VEC),cgetg(1,t_MAT));
173 : }
174 21 : if (n != nbrows(x)) pari_err_DIM("matqr");
175 21 : if (!RgM_QR_init(x, &B,&Q,&L, prec)) pari_err_PREC("matqr");
176 21 : if (!flag) Q = shallowtrans(mathouseholder(Q, matid(n)));
177 21 : return gc_GEN(av, mkvec2(Q, shallowtrans(L)));
178 : }
179 :
180 : /* compute B = squared length of orthogonalized vectors x[k]^*,
181 : * Q = Householder transforms and L = mu_{i,j}. B[k] a t_REAL for k > 1,
182 : * t_INT/t_REAL for k = 1; L[j,j] a t_REAL for j < #x */
183 : int
184 7753689 : QR_init(GEN x, GEN *pB, GEN *pQ, GEN *pL, long prec)
185 : {
186 7753689 : long j, k = lg(x)-1;
187 7753689 : GEN B = cgetg(k+1, t_VEC), Q = cgetg(k, t_VEC), L = zeromatcopy(k,k);
188 30244634 : for (j=1; j<=k; j++)
189 : {
190 24734248 : GEN r = gel(x,j);
191 24734248 : if (j > 1) r = ApplyAllQ(Q, r, j);
192 24734524 : if (!FindApplyQ(r, L, B, j, Q, prec)) return 0;
193 : }
194 5510386 : *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 300857 : gaussred_from_QR(GEN x, long prec)
200 : {
201 300857 : long j, k = lg(x)-1;
202 : GEN B, Q, L;
203 300857 : if (!QR_init(x, &B,&Q,&L, prec)) return NULL;
204 1073035 : for (j=1; j<k; j++)
205 : {
206 772181 : GEN m = gel(L,j), invNx = invr(gel(m,j));
207 : long i;
208 772168 : gel(m,j) = gel(B,j);
209 2980609 : for (i=j+1; i<=k; i++) gel(m,i) = mpmul(invNx, gel(m,i));
210 : }
211 300854 : gcoeff(L,j,j) = gel(B,j); /* t_REAL for j > 1, t_INT or t_REAL for j = 1 */
212 300854 : return shallowtrans(L);
213 : }
214 : GEN
215 14280 : R_from_QR(GEN x, long prec)
216 : {
217 : GEN B, Q, L;
218 14280 : if (!QR_init(x, &B,&Q,&L, prec)) return NULL;
219 14266 : 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 56780 : RgM_gram_schmidt(GEN e, GEN *ptB)
230 : {
231 56780 : long i,j,lx = lg(e);
232 56780 : GEN f = RgM_shallowcopy(e), B, iB;
233 :
234 56780 : B = cgetg(lx, t_VEC);
235 56780 : iB= cgetg(lx, t_VEC);
236 :
237 120645 : for (i=1;i<lx;i++)
238 : {
239 63865 : GEN p1 = NULL;
240 63865 : pari_sp av = avma;
241 126602 : for (j=1; j<i; j++)
242 : {
243 62737 : GEN mu = gmul(RgV_dotproduct(gel(e,i),gel(f,j)), gel(iB,j));
244 62737 : GEN p2 = gmul(mu, gel(f,j));
245 62737 : p1 = p1? gadd(p1,p2): p2;
246 : }
247 63865 : p1 = p1? gc_upto(av, gsub(gel(e,i), p1)): gel(e,i);
248 63865 : gel(f,i) = p1;
249 63865 : gel(B,i) = RgV_dotsquare(gel(f,i));
250 63865 : gel(iB,i) = ginv(gel(B,i));
251 : }
252 56780 : *ptB = B; return f;
253 : }
254 :
255 : /* B a Z-basis (which the caller should LLL-reduce for efficiency), t a vector.
256 : * Apply Babai's nearest plane algorithm to (B,t) */
257 : GEN
258 56780 : RgM_Babai(GEN B, GEN t)
259 : {
260 56780 : GEN C, N, G = RgM_gram_schmidt(B, &N), b = t;
261 56780 : long j, n = lg(B)-1;
262 :
263 56780 : C = cgetg(n+1,t_COL);
264 120645 : for (j = n; j > 0; j--)
265 : {
266 63865 : GEN c = gdiv( RgV_dotproduct(b, gel(G,j)), gel(N,j) );
267 : long e;
268 63865 : c = grndtoi(c,&e);
269 63865 : if (e >= 0) return NULL;
270 63865 : if (signe(c)) b = RgC_sub(b, RgC_Rg_mul(gel(B,j), c));
271 63865 : gel(C,j) = c;
272 : }
273 56780 : 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 49 : 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 : (void)gc_all(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 35 : gcoeff(tm1,1,i) = addii(mulii(q1neg, gcoeff(tm,1,1)),
389 35 : mulii(q2neg, gcoeff(tm,1,2)));
390 35 : 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 : (void)gc_all(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 35 : {
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 700 : gcoeff(dot,k1,k1) = addii(gcoeff(dot,k1,k1),
442 700 : mulii(q, gcoeff(dot,k2,k1)));
443 5600 : 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 : (void)gc_all(av3, 2, &dot,&tm2);
449 : }
450 : } /* for i */
451 42 : if (!reductions) break;
452 35 : e = 0;
453 280 : 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 : } /* 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 gc_upto(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(double b, double x, double rho, long d, long dim, long delta, long t)
496 : {
497 882 : double cond = delta * (d * (delta+1) - 2*b*dim + rho * (delta-1 + 2*t))
498 882 : + x*dim*(dim - 1);
499 882 : if (DEBUGLEVEL >= 4)
500 0 : err_printf("delta = %d, t = %d (%.1lf)\n", delta, t, cond);
501 882 : return (cond <= 0);
502 : }
503 :
504 : static void
505 21 : choose_params(GEN P, GEN N, GEN X, GEN B, long *pdelta, long *pt)
506 : {
507 21 : long d = degpol(P), dim;
508 21 : GEN P0 = leading_coeff(P);
509 21 : double logN = dbllog2(N), x, b, rho;
510 21 : x = dbllog2(X) / logN;
511 21 : b = B? dbllog2(B) / logN: 1.;
512 21 : if (x * d >= b * b) pari_err_OVERFLOW("zncoppersmith [bound too large]");
513 : /* TODO : remove P0 completely */
514 14 : rho = is_pm1(P0)? 0: dbllog2(P0) / logN;
515 :
516 : /* Enumerate (delta,t) by increasing lattice dimension */
517 14 : for(dim = d + 1;; dim++)
518 161 : {
519 : long delta, t; /* dim = d*delta + t in the loop */
520 1043 : for (delta = 0, t = dim; t >= 0; delta++, t -= d)
521 882 : if (check(b,x,rho,d,dim,delta,t)) { *pdelta = delta; *pt = t; return; }
522 : }
523 : }
524 :
525 : static int
526 14021 : sol_OK(GEN x, GEN N, GEN B)
527 14021 : { return B? (cmpii(gcdii(x,N),B) >= 0): dvdii(x,N); }
528 : /* deg(P) > 0, x >= 0. Find all j such that gcd(P(j), N) >= B, |j| <= x */
529 : static GEN
530 7 : do_exhaustive(GEN P, GEN N, long x, GEN B)
531 : {
532 7 : GEN Pe, Po, sol = vecsmalltrunc_init(2*x + 2);
533 : pari_sp av;
534 : long j;
535 7 : RgX_even_odd(P, &Pe,&Po); av = avma;
536 7 : if (sol_OK(gel(P,2), N,B)) vecsmalltrunc_append(sol, 0);
537 7007 : for (j = 1; j <= x; j++, set_avma(av))
538 : {
539 7000 : GEN j2 = sqru(j), E = FpX_eval(Pe,j2,N), O = FpX_eval(Po,j2,N);
540 7000 : if (sol_OK(addmuliu(E,O,j), N,B)) vecsmalltrunc_append(sol, j);
541 7000 : if (sol_OK(submuliu(E,O,j), N,B)) vecsmalltrunc_append(sol,-j);
542 : }
543 7 : vecsmall_sort(sol); return zv_to_ZV(sol);
544 : }
545 :
546 : /* General Coppersmith, look for a root x0 <= p, p >= B, p | N, |x0| <= X.
547 : * B = N coded as NULL */
548 : GEN
549 35 : zncoppersmith(GEN P, GEN N, GEN X, GEN B)
550 : {
551 : GEN Q, R, N0, M, sh, short_pol, *Xpowers, sol, nsp, cP, Z;
552 35 : long delta, i, j, row, d, l, t, dim, bnd = 10;
553 35 : const ulong X_SMALL = 1000;
554 35 : pari_sp av = avma;
555 :
556 35 : if (typ(P) != t_POL || !RgX_is_ZX(P)) pari_err_TYPE("zncoppersmith",P);
557 28 : if (typ(N) != t_INT) pari_err_TYPE("zncoppersmith",N);
558 28 : if (typ(X) != t_INT) {
559 7 : X = gfloor(X);
560 7 : if (typ(X) != t_INT) pari_err_TYPE("zncoppersmith",X);
561 : }
562 28 : if (signe(X) < 0) pari_err_DOMAIN("zncoppersmith", "X", "<", gen_0, X);
563 28 : P = FpX_red(P, N); d = degpol(P);
564 28 : if (d == 0) retgc_const(av, cgetg(1, t_VEC));
565 28 : if (d < 0) pari_err_ROOTS0("zncoppersmith");
566 28 : if (B && typ(B) != t_INT) B = gceil(B);
567 28 : if (abscmpiu(X, X_SMALL) <= 0)
568 7 : return gc_upto(av, do_exhaustive(P, N, itos(X), B));
569 :
570 21 : if (B && equalii(B,N)) B = NULL;
571 21 : if (B) bnd = 1; /* bnd-hack is only for the case B = N */
572 21 : cP = gel(P,d+2);
573 21 : if (!gequal1(cP))
574 : {
575 : GEN r, z;
576 14 : gel(P,d+2) = cP = bezout(cP, N, &z, &r);
577 35 : for (j = 0; j < d; j++) gel(P,j+2) = Fp_mul(gel(P,j+2), z, N);
578 14 : if (!is_pm1(cP))
579 : {
580 7 : P = Q_primitive_part(P, &cP);
581 7 : if (cP) { N = diviiexact(N,cP); B = gceil(gdiv(B, cP)); }
582 : }
583 : }
584 21 : if (DEBUGLEVEL >= 2) err_printf("Modified P: %Ps\n", P);
585 :
586 21 : choose_params(P, N, X, B, &delta, &t);
587 14 : if (DEBUGLEVEL >= 2)
588 0 : err_printf("Init: trying delta = %d, t = %d\n", delta, t);
589 : for(;;)
590 : {
591 14 : dim = d * delta + t;
592 : /* TODO: In case of failure do not recompute the full vector */
593 14 : Xpowers = (GEN*)new_chunk(dim + 1);
594 14 : Xpowers[0] = gen_1;
595 217 : for (j = 1; j <= dim; j++) Xpowers[j] = mulii(Xpowers[j-1], X);
596 :
597 : /* TODO: in case of failure, use the part of the matrix already computed */
598 14 : M = zeromatcopy(dim,dim);
599 :
600 : /* Rows of M correspond to the polynomials
601 : * N^delta, N^delta Xi, ... N^delta (Xi)^d-1,
602 : * N^(delta-1)P(Xi), N^(delta-1)XiP(Xi), ... N^(delta-1)P(Xi)(Xi)^d-1,
603 : * ...
604 : * P(Xi)^delta, XiP(Xi)^delta, ..., P(Xi)^delta(Xi)^t-1 */
605 42 : for (j = 1; j <= d; j++) gcoeff(M, j, j) = gel(Xpowers,j-1);
606 :
607 : /* P-part */
608 14 : if (delta) row = d + 1; else row = 0;
609 :
610 14 : Q = P;
611 70 : for (i = 1; i < delta; i++)
612 : {
613 182 : for (j = 0; j < d; j++,row++)
614 1239 : for (l = j + 1; l <= row; l++)
615 1113 : gcoeff(M, l, row) = mulii(Xpowers[l-1], gel(Q,l-j+1));
616 56 : Q = ZX_mul(Q, P);
617 : }
618 63 : for (j = 0; j < t; row++, j++)
619 490 : for (l = j + 1; l <= row; l++)
620 441 : gcoeff(M, l, row) = mulii(Xpowers[l-1], gel(Q,l-j+1));
621 :
622 : /* N-part */
623 14 : row = dim - t; N0 = N;
624 84 : while (row >= 1)
625 : {
626 224 : for (j = 0; j < d; j++,row--)
627 1421 : for (l = 1; l <= row; l++)
628 1267 : gcoeff(M, l, row) = mulii(gmael(M, row, l), N0);
629 70 : if (row >= 1) N0 = mulii(N0, N);
630 : }
631 : /* Z is the upper bound for the L^1 norm of the polynomial,
632 : ie. N^delta if B = N, B^delta otherwise */
633 14 : if (B) Z = powiu(B, delta); else Z = N0;
634 :
635 14 : if (DEBUGLEVEL >= 2)
636 : {
637 0 : if (DEBUGLEVEL >= 6) err_printf("Matrix to be reduced:\n%Ps\n", M);
638 0 : err_printf("Entering LLL\nbitsize bound: %ld\n", expi(Z));
639 0 : err_printf("expected shvector bitsize: %ld\n", expi(ZM_det_triangular(M))/dim);
640 : }
641 :
642 14 : sh = ZM_lll(M, 0.75, LLL_INPLACE);
643 : /* Take the first vector if it is non constant */
644 14 : short_pol = gel(sh,1);
645 14 : if (ZV_isscalar(short_pol)) short_pol = gel(sh, 2);
646 :
647 14 : nsp = gen_0;
648 217 : for (j = 1; j <= dim; j++) nsp = addii(nsp, absi_shallow(gel(short_pol,j)));
649 :
650 14 : if (DEBUGLEVEL >= 2)
651 : {
652 0 : err_printf("Candidate: %Ps\n", short_pol);
653 0 : err_printf("bitsize Norm: %ld\n", expi(nsp));
654 0 : err_printf("bitsize bound: %ld\n", expi(mului(bnd, Z)));
655 : }
656 14 : if (cmpii(nsp, mului(bnd, Z)) < 0) break; /* SUCCESS */
657 :
658 : /* Failed with the precomputed or supplied value */
659 0 : if (++t == d) { delta++; t = 1; }
660 0 : if (DEBUGLEVEL >= 2)
661 0 : err_printf("Increasing dim, delta = %d t = %d\n", delta, t);
662 : }
663 14 : bnd = itos(divii(nsp, Z)) + 1;
664 :
665 14 : while (!signe(gel(short_pol,dim))) dim--;
666 :
667 14 : R = cgetg(dim + 2, t_POL); R[1] = P[1];
668 217 : for (j = 1; j <= dim; j++)
669 203 : gel(R,j+1) = diviiexact(gel(short_pol,j), Xpowers[j-1]);
670 14 : gel(R,2) = subii(gel(R,2), mului(bnd - 1, N0));
671 :
672 14 : sol = cgetg(1, t_VEC);
673 84 : for (i = -bnd + 1; i < bnd; i++)
674 : {
675 70 : GEN r = nfrootsQ(R);
676 70 : if (DEBUGLEVEL >= 2) err_printf("Roots: %Ps\n", r);
677 91 : for (j = 1; j < lg(r); j++)
678 : {
679 21 : GEN z = gel(r,j);
680 21 : if (typ(z) == t_INT && sol_OK(FpX_eval(P,z,N), N,B))
681 14 : sol = shallowconcat(sol, z);
682 : }
683 70 : if (i < bnd) gel(R,2) = addii(gel(R,2), Z);
684 : }
685 14 : return gc_upto(av, ZV_sort_uniq(sol));
686 : }
687 :
688 : /********************************************************************/
689 : /** **/
690 : /** LINEAR & ALGEBRAIC DEPENDENCE **/
691 : /** **/
692 : /********************************************************************/
693 :
694 : static int
695 8123 : real_indep(GEN re, GEN im, long bit)
696 : {
697 8123 : GEN d = gsub(gmul(gel(re,1),gel(im,2)), gmul(gel(re,2),gel(im,1)));
698 8123 : return (!gequal0(d) && gexpo(d) > - bit);
699 : }
700 :
701 : GEN
702 15302 : lindepfull_bit(GEN x, long bit)
703 : {
704 15302 : long lx = lg(x), ly, i, j;
705 : GEN re, im, M;
706 :
707 15302 : if (! is_vec_t(typ(x))) pari_err_TYPE("lindep2",x);
708 15302 : if (lx <= 2)
709 : {
710 21 : if (lx == 2 && gequal0(x)) return matid(1);
711 14 : return NULL;
712 : }
713 15281 : re = real_i(x);
714 15281 : im = imag_i(x);
715 : /* independent over R ? */
716 15281 : if (lx == 3 && real_indep(re,im,bit)) return NULL;
717 15267 : if (gequal0(im)) im = NULL;
718 15267 : ly = im? lx+2: lx+1;
719 15267 : M = cgetg(lx,t_MAT);
720 60701 : for (i=1; i<lx; i++)
721 : {
722 45434 : GEN c = cgetg(ly,t_COL); gel(M,i) = c;
723 209394 : for (j=1; j<lx; j++) gel(c,j) = gen_0;
724 45434 : gel(c,i) = gen_1;
725 45434 : gel(c,lx) = gtrunc2n(gel(re,i), bit);
726 45434 : if (im) gel(c,lx+1) = gtrunc2n(gel(im,i), bit);
727 : }
728 15267 : return ZM_lll(M, 0.99, LLL_INPLACE);
729 : }
730 : GEN
731 3311 : lindep_bit(GEN x, long bit)
732 : {
733 3311 : pari_sp av = avma;
734 3311 : GEN v, M = lindepfull_bit(x,bit);
735 3311 : if (!M) retgc_const(av, cgetg(1, t_COL));
736 3283 : v = gel(M,1); setlg(v, lg(M));
737 3283 : return gc_GEN(av, v);
738 : }
739 : /* deprecated */
740 : GEN
741 112 : lindep2(GEN x, long dig)
742 : {
743 : long bit;
744 112 : if (dig < 0) pari_err_DOMAIN("lindep2", "accuracy", "<", gen_0, stoi(dig));
745 112 : if (dig) bit = (long) (dig/LOG10_2);
746 : else
747 : {
748 98 : bit = gprecision(x);
749 98 : if (!bit)
750 : {
751 35 : x = Q_primpart(x); /* left on stack */
752 35 : bit = 32 + gexpo(x);
753 : }
754 : else
755 63 : bit = (long)prec2nbits_mul(bit, 0.8);
756 : }
757 112 : return lindep_bit(x, bit);
758 : }
759 :
760 : /* x is a vector of elts of a p-adic field */
761 : GEN
762 28 : lindep_padic(GEN x)
763 : {
764 28 : long i, j, prec = LONG_MAX, nx = lg(x)-1, v;
765 28 : pari_sp av = avma;
766 28 : GEN p = NULL, pn, m, a;
767 :
768 28 : if (nx < 2) return cgetg(1,t_COL);
769 147 : for (i=1; i<=nx; i++)
770 : {
771 119 : GEN c = gel(x,i), q;
772 119 : if (typ(c) != t_PADIC) continue;
773 :
774 91 : j = precp(c); if (j < prec) prec = j;
775 91 : q = padic_p(c);
776 91 : if (!p) p = q; else if (!equalii(p, q)) pari_err_MODULUS("lindep_padic", p, q);
777 : }
778 28 : if (!p) pari_err_TYPE("lindep_padic [not a p-adic vector]",x);
779 28 : v = gvaluation(x,p); pn = powiu(p,prec);
780 28 : if (v) x = gmul(x, powis(p, -v));
781 28 : x = RgV_to_FpV(x, pn);
782 :
783 28 : a = negi(gel(x,1));
784 28 : m = cgetg(nx,t_MAT);
785 119 : for (i=1; i<nx; i++)
786 : {
787 91 : GEN c = zerocol(nx);
788 91 : gel(c,1+i) = a;
789 91 : gel(c,1) = gel(x,i+1);
790 91 : gel(m,i) = c;
791 : }
792 28 : m = ZM_lll(ZM_hnfmodid(m, pn), 0.99, LLL_INPLACE);
793 28 : return gc_GEN(av, gel(m,1));
794 : }
795 : /* x is a vector of t_POL/t_SER */
796 : GEN
797 77 : lindep_Xadic(GEN x)
798 : {
799 77 : long i, prec = LONG_MAX, deg = 0, lx = lg(x), vx, v;
800 77 : pari_sp av = avma;
801 : GEN m;
802 :
803 77 : if (lx == 1) return cgetg(1,t_COL);
804 77 : vx = gvar(x);
805 77 : if (gequal0(x)) return col_ei(lx-1,1);
806 70 : v = gvaluation(x, pol_x(vx));
807 70 : if (!v) x = shallowcopy(x);
808 0 : else if (v > 0) x = gdiv(x, pol_xn(v, vx));
809 0 : else x = gmul(x, pol_xn(-v, vx));
810 : /* all t_SER have valuation >= 0 */
811 735 : for (i=1; i<lx; i++)
812 : {
813 665 : GEN c = gel(x,i);
814 665 : if (gvar(c) != vx) { gel(x,i) = scalarpol_shallow(c, vx); continue; }
815 658 : switch(typ(c))
816 : {
817 231 : case t_POL: deg = maxss(deg, degpol(c)); break;
818 0 : case t_RFRAC: pari_err_TYPE("lindep_Xadic", c);
819 427 : case t_SER:
820 427 : prec = minss(prec, valser(c)+lg(c)-2);
821 427 : gel(x,i) = ser2rfrac_i(c);
822 : }
823 : }
824 70 : if (prec == LONG_MAX) prec = deg+1;
825 70 : m = RgXV_to_RgM(x, prec);
826 70 : return gc_upto(av, deplin(m));
827 : }
828 : static GEN
829 35 : vec_lindep(GEN x)
830 : {
831 35 : pari_sp av = avma;
832 35 : long i, l = lg(x); /* > 1 */
833 35 : long t = typ(gel(x,1)), h = lg(gel(x,1));
834 35 : GEN m = cgetg(l, t_MAT);
835 126 : for (i = 1; i < l; i++)
836 : {
837 98 : GEN y = gel(x,i);
838 98 : if (lg(y) != h || typ(y) != t) pari_err_TYPE("lindep",x);
839 91 : if (t != t_COL) y = shallowtrans(y); /* Sigh */
840 91 : gel(m,i) = y;
841 : }
842 28 : return gc_upto(av, deplin(m));
843 : }
844 :
845 : GEN
846 0 : lindep(GEN x) { return lindep2(x, 0); }
847 :
848 : GEN
849 434 : lindep0(GEN x,long bit)
850 : {
851 434 : long i, tx = typ(x);
852 434 : if (tx == t_MAT) return deplin(x);
853 147 : if (! is_vec_t(tx)) pari_err_TYPE("lindep",x);
854 441 : for (i = 1; i < lg(x); i++)
855 357 : switch(typ(gel(x,i)))
856 : {
857 7 : case t_PADIC: return lindep_padic(x);
858 21 : case t_POL:
859 : case t_RFRAC:
860 21 : case t_SER: return lindep_Xadic(x);
861 35 : case t_VEC:
862 35 : case t_COL: return vec_lindep(x);
863 : }
864 84 : return lindep2(x, bit);
865 : }
866 :
867 : GEN
868 77 : algdep0(GEN x, long n, long bit)
869 : {
870 77 : long tx = typ(x), i;
871 : pari_sp av;
872 : GEN y;
873 :
874 77 : if (! is_scalar_t(tx)) pari_err_TYPE("algdep0",x);
875 77 : if (tx == t_POLMOD)
876 : {
877 14 : av = avma; y = minpoly(x, 0);
878 14 : return (degpol(y) > n)? gc_const(av, gen_1): y;
879 : }
880 63 : if (gequal0(x)) return pol_x(0);
881 63 : 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 49 : av = avma; y = cgetg(n+2,t_COL);
888 49 : gel(y,1) = gen_1;
889 49 : gel(y,2) = x; /* n >= 1 */
890 210 : for (i=3; i<=n+1; i++) gel(y,i) = gmul(gel(y,i-1),x);
891 49 : if (typ(x) == t_PADIC)
892 21 : y = lindep_padic(y);
893 : else
894 28 : y = lindep2(y, bit);
895 49 : if (lg(y) == 1) pari_err(e_DOMAIN,"algdep", "degree(x)",">", stoi(n), x);
896 49 : y = RgV_to_RgX(y, 0);
897 49 : if (signe(leading_coeff(y)) > 0) return gc_GEN(av, y);
898 14 : return gc_upto(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 : static GEN
908 56 : sertomat(GEN S, long p, long r, long vy)
909 : {
910 : long n, m;
911 56 : GEN v = cgetg(r*p+1, t_VEC); /* v[r*n+m+1] = s^n * y^m */
912 : /* n = 0 */
913 245 : for (m = 0; m < r; m++) gel(v, m+1) = pol_xn(m, vy);
914 175 : for(n=1; n < p; n++)
915 546 : for (m = 0; m < r; m++)
916 : {
917 427 : GEN c = gel(S,n);
918 427 : if (m)
919 : {
920 308 : c = shallowcopy(c);
921 308 : setvalser(c, valser(c) + m);
922 : }
923 427 : gel(v, r*n + m + 1) = c;
924 : }
925 56 : return v;
926 : }
927 :
928 : GEN
929 42 : seralgdep(GEN s, long p, long r)
930 : {
931 42 : pari_sp av = avma;
932 : long vy, i, n, prec;
933 : GEN S, v, D;
934 :
935 42 : if (typ(s) != t_SER) pari_err_TYPE("seralgdep",s);
936 42 : if (p <= 0) pari_err_DOMAIN("seralgdep", "p", "<=", gen_0, stoi(p));
937 42 : if (r < 0) pari_err_DOMAIN("seralgdep", "r", "<", gen_0, stoi(r));
938 42 : if (is_bigint(addiu(muluu(p, r), 1))) pari_err_OVERFLOW("seralgdep");
939 42 : vy = varn(s);
940 42 : if (!vy) pari_err_PRIORITY("seralgdep", s, ">", 0);
941 42 : r++; p++;
942 42 : prec = valser(s) + lg(s)-2;
943 42 : if (r > prec) r = prec;
944 42 : S = cgetg(p+1, t_VEC); gel(S, 1) = s;
945 119 : for (i = 2; i <= p; i++) gel(S,i) = gmul(gel(S,i-1), s);
946 42 : v = sertomat(S, p, r, vy);
947 42 : D = lindep_Xadic(v);
948 42 : if (lg(D) == 1) { set_avma(av); return gen_0; }
949 35 : v = cgetg(p+1, t_VEC);
950 133 : for (n = 0; n < p; n++)
951 98 : gel(v, n+1) = RgV_to_RgX(vecslice(D, r*n+1, r*n+r), vy);
952 35 : return gc_GEN(av, RgV_to_RgX(v, 0));
953 : }
954 :
955 : GEN
956 14 : serdiffdep(GEN s, long p, long r)
957 : {
958 14 : pari_sp av = avma;
959 : long vy, i, n, prec;
960 : GEN P, S, v, D;
961 :
962 14 : if (typ(s) != t_SER) pari_err_TYPE("serdiffdep",s);
963 14 : if (p <= 0) pari_err_DOMAIN("serdiffdep", "p", "<=", gen_0, stoi(p));
964 14 : if (r < 0) pari_err_DOMAIN("serdiffdep", "r", "<", gen_0, stoi(r));
965 14 : if (is_bigint(addiu(muluu(p, r), 1))) pari_err_OVERFLOW("serdiffdep");
966 14 : vy = varn(s);
967 14 : if (!vy) pari_err_PRIORITY("serdiffdep", s, ">", 0);
968 14 : r++; p++;
969 14 : prec = valser(s) + lg(s)-2;
970 14 : if (r > prec) r = prec;
971 14 : S = cgetg(p+1, t_VEC); gel(S, 1) = s;
972 56 : for (i = 2; i <= p; i++) gel(S,i) = derivser(gel(S,i-1));
973 14 : v = sertomat(S, p, r, vy);
974 14 : D = lindep_Xadic(v);
975 14 : if (lg(D) == 1) { set_avma(av); return gen_0; }
976 14 : P = RgV_to_RgX(vecslice(D, 1, r), vy);
977 14 : v = cgetg(p, t_VEC);
978 56 : for (n = 1; n < p; n++)
979 42 : gel(v, n) = RgV_to_RgX(vecslice(D, r*n+1, r*n+r), vy);
980 14 : return gc_GEN(av, mkvec2(RgV_to_RgX(v, 0), gneg(P)));
981 : }
982 :
983 : /* FIXME: could precompute ZM_lll attached to V[2..] */
984 : static GEN
985 11991 : lindepcx(GEN V, long bit)
986 : {
987 11991 : GEN Vr = real_i(V), Vi = imag_i(V);
988 11991 : long d = gexpo(Vr) - gexpo(Vi);
989 11991 : if (d < -bit) V = Vi;
990 11991 : else if (d > bit) V = Vr;
991 11991 : return lindepfull_bit(V, bit);
992 : }
993 : /* c floating point t_REAL or t_COMPLEX, T ZX, recognize in Q[x]/(T).
994 : * V helper vector (containing complex roots of T), MODIFIED */
995 : static GEN
996 11991 : cx_bestapprnf(GEN c, GEN T, GEN V, long bit)
997 : {
998 11991 : GEN M, a, v = NULL;
999 : long i, l;
1000 11991 : gel(V,1) = gneg(c); M = lindepcx(V, bit);
1001 11991 : if (!M) pari_err(e_MISC, "cannot rationalize coeff in bestapprnf");
1002 11991 : l = lg(M); a = NULL;
1003 11991 : for (i = 1; i < l; i ++) { v = gel(M,i); a = gel(v,1); if (signe(a)) break; }
1004 11991 : v = RgC_Rg_div(vecslice(v, 2, lg(M)-1), a);
1005 11991 : if (!T) return gel(v,1);
1006 4830 : v = RgV_to_RgX(v, varn(T)); l = lg(v);
1007 4830 : if (l == 2) return gen_0;
1008 4165 : if (l == 3) return gel(v,2);
1009 3668 : return mkpolmod(v, T);
1010 : }
1011 : static GEN
1012 14784 : bestapprnf_i(GEN x, GEN T, GEN V, long bit)
1013 : {
1014 14784 : long i, l, tx = typ(x);
1015 : GEN z;
1016 14784 : switch (tx)
1017 : {
1018 1505 : case t_INT: case t_FRAC: return x;
1019 11991 : case t_REAL: case t_COMPLEX: return cx_bestapprnf(x, T, V, bit);
1020 0 : case t_POLMOD: if (RgX_equal(gel(x,1),T)) return x;
1021 0 : break;
1022 1288 : case t_POL: case t_SER: case t_VEC: case t_COL: case t_MAT:
1023 1288 : l = lg(x); z = cgetg(l, tx);
1024 1974 : for (i = 1; i < lontyp[tx]; i++) z[i] = x[i];
1025 13993 : for (; i < l; i++) gel(z,i) = bestapprnf_i(gel(x,i), T, V, bit);
1026 1288 : return z;
1027 : }
1028 0 : pari_err_TYPE("mfcxtoQ", x);
1029 : return NULL;/*LCOV_EXCL_LINE*/
1030 : }
1031 :
1032 : GEN
1033 2163 : bestapprnf(GEN x, GEN T, GEN roT, long prec)
1034 : {
1035 2163 : pari_sp av = avma;
1036 2163 : long tx = typ(x), dT = 1, bit;
1037 : GEN V;
1038 :
1039 2163 : if (T)
1040 : {
1041 1610 : if (typ(T) != t_POL) T = nf_get_pol(checknf(T));
1042 1610 : else if (!RgX_is_ZX(T)) pari_err_TYPE("bestapprnf", T);
1043 1610 : dT = degpol(T);
1044 : }
1045 2163 : if (is_rational_t(tx)) return gcopy(x);
1046 2079 : if (tx == t_POLMOD)
1047 : {
1048 0 : if (!T || !RgX_equal(T, gel(x,1))) pari_err_TYPE("bestapprnf",x);
1049 0 : return gcopy(x);
1050 : }
1051 :
1052 2079 : if (roT)
1053 : {
1054 644 : long l = gprecision(roT);
1055 644 : switch(typ(roT))
1056 : {
1057 644 : case t_INT: case t_FRAC: case t_REAL: case t_COMPLEX: break;
1058 0 : default: pari_err_TYPE("bestapprnf", roT);
1059 : }
1060 644 : if (prec < l) prec = l;
1061 : }
1062 1435 : else if (!T)
1063 525 : roT = gen_1;
1064 : else
1065 : {
1066 910 : long n = poliscyclo(T); /* cyclotomic is an important special case */
1067 910 : roT = n? rootsof1u_cx(n,prec): gel(QX_complex_roots(T,prec), 1);
1068 : }
1069 2079 : V = vec_prepend(gpowers(roT, dT-1), NULL);
1070 2079 : bit = prec2nbits_mul(prec, 0.8);
1071 2079 : return gc_GEN(av, bestapprnf_i(x, T, V, bit));
1072 : }
1073 :
1074 : /********************************************************************/
1075 : /** **/
1076 : /** MINIM **/
1077 : /** **/
1078 : /********************************************************************/
1079 : void
1080 124212 : minim_alloc(long n, double ***q, GEN *x, double **y, double **z, double **v)
1081 : {
1082 124212 : long i, s = n * sizeof(double);
1083 :
1084 124212 : *x = cgetg(n, t_VECSMALL);
1085 124211 : *q = (double**) new_chunk(n);
1086 124212 : *y = (double*) stack_malloc_align(s, sizeof(double));
1087 124212 : *z = (double*) stack_malloc_align(s, sizeof(double));
1088 124212 : *v = (double*) stack_malloc_align(s, sizeof(double));
1089 535774 : for (i=1; i<n; i++) (*q)[i] = (double*) stack_malloc_align(s, sizeof(double));
1090 124214 : }
1091 :
1092 : static void
1093 70 : cvp_alloc(long n, double **t, double **tpre)
1094 : {
1095 70 : long s = n * sizeof(double);
1096 70 : *t = (double*) stack_malloc_align(s, sizeof(double));
1097 70 : *tpre = (double*) stack_malloc_align(s, sizeof(double));
1098 70 : }
1099 :
1100 : static GEN
1101 5502 : ZC_canon(GEN V)
1102 : {
1103 5502 : long l = lg(V), j, s;
1104 11242 : for (j = 1; j < l; j++)
1105 11242 : if ((s = signe(gel(V,j)))) return s < 0? ZC_neg(V): V;
1106 0 : return V;
1107 : }
1108 : static GEN
1109 5502 : ZM_zc_mul_canon(GEN u, GEN x) { return ZC_canon(ZM_zc_mul(u,x)); }
1110 : static GEN
1111 240366 : ZM_zc_mul_canon_zm(GEN u, GEN x)
1112 : {
1113 240366 : pari_sp av = avma;
1114 240366 : GEN y = ZV_to_zv(ZM_zc_mul(u,x));
1115 240366 : zv_canon_inplace(y); return gc_upto(av, y);
1116 : }
1117 :
1118 : struct qfvec
1119 : {
1120 : GEN a, r, u;
1121 : };
1122 :
1123 : static void
1124 0 : err_minim(GEN a)
1125 : {
1126 0 : pari_err_DOMAIN("minim0","form","is not",
1127 : strtoGENstr("positive definite"),a);
1128 0 : }
1129 :
1130 : static GEN
1131 902 : minim_lll(GEN a, GEN *u)
1132 : {
1133 902 : *u = lllgramint(a);
1134 902 : if (lg(*u) != lg(a)) err_minim(a);
1135 902 : return qf_ZM_apply(a,*u);
1136 : }
1137 :
1138 : static void
1139 902 : forqfvec_init_dolll(struct qfvec *qv, GEN *pa, long dolll)
1140 : {
1141 902 : GEN r, u, a = *pa;
1142 902 : if (!dolll) u = NULL;
1143 : else
1144 : {
1145 860 : if (typ(a) != t_MAT || !RgM_is_ZM(a)) pari_err_TYPE("qfminim",a);
1146 860 : a = *pa = minim_lll(a, &u);
1147 : }
1148 902 : qv->a = RgM_gtofp(a, DEFAULTPREC);
1149 902 : r = qfgaussred_positive(qv->a);
1150 902 : if (!r)
1151 : {
1152 0 : r = qfgaussred_positive(a); /* exact computation */
1153 0 : if (!r) err_minim(a);
1154 0 : r = RgM_gtofp(r, DEFAULTPREC);
1155 : }
1156 902 : qv->r = r;
1157 902 : qv->u = u;
1158 902 : }
1159 :
1160 : static void
1161 42 : forqfvec_init(struct qfvec *qv, GEN a)
1162 42 : { forqfvec_init_dolll(qv, &a, 1); }
1163 :
1164 : static void
1165 42 : forqfvec_i(void *E, long (*fun)(void *, GEN, GEN, double), struct qfvec *qv, GEN BORNE)
1166 : {
1167 42 : GEN x, a = qv->a, r = qv->r, u = qv->u;
1168 42 : long n = lg(a)-1, i, j, k;
1169 : double p,BOUND,*v,*y,*z,**q;
1170 42 : const double eps = 1e-10;
1171 42 : if (!BORNE) BORNE = gen_0;
1172 : else
1173 : {
1174 28 : BORNE = gfloor(BORNE);
1175 28 : if (typ(BORNE) != t_INT) pari_err_TYPE("minim0",BORNE);
1176 35 : if (signe(BORNE) <= 0) return;
1177 : }
1178 35 : if (n == 0) return;
1179 28 : minim_alloc(n+1, &q, &x, &y, &z, &v);
1180 98 : for (j=1; j<=n; j++)
1181 : {
1182 70 : v[j] = rtodbl(gcoeff(r,j,j));
1183 133 : for (i=1; i<j; i++) q[i][j] = rtodbl(gcoeff(r,i,j));
1184 : }
1185 :
1186 28 : if (gequal0(BORNE))
1187 : {
1188 : double c;
1189 14 : p = rtodbl(gcoeff(a,1,1));
1190 42 : for (i=2; i<=n; i++) { c = rtodbl(gcoeff(a,i,i)); if (c < p) p = c; }
1191 14 : BORNE = roundr(dbltor(p));
1192 : }
1193 : else
1194 14 : p = gtodouble(BORNE);
1195 28 : BOUND = p * (1 + eps);
1196 28 : if (BOUND > (double)ULONG_MAX || (ulong)BOUND != (ulong)p)
1197 7 : pari_err_PREC("forqfvec");
1198 :
1199 21 : k = n; y[n] = z[n] = 0;
1200 21 : x[n] = (long)sqrt(BOUND/v[n]);
1201 56 : for(;;x[1]--)
1202 : {
1203 : do
1204 : {
1205 140 : if (k>1)
1206 : {
1207 84 : long l = k-1;
1208 84 : z[l] = 0;
1209 245 : for (j=k; j<=n; j++) z[l] += q[l][j]*x[j];
1210 84 : p = (double)x[k] + z[k];
1211 84 : y[l] = y[k] + p*p*v[k];
1212 84 : x[l] = (long)floor(sqrt((BOUND-y[l])/v[l])-z[l]);
1213 84 : k = l;
1214 : }
1215 : for(;;)
1216 : {
1217 189 : p = (double)x[k] + z[k];
1218 189 : if (y[k] + p*p*v[k] <= BOUND) break;
1219 49 : k++; x[k]--;
1220 : }
1221 140 : } while (k > 1);
1222 77 : if (! x[1] && y[1]<=eps) break;
1223 :
1224 56 : p = (double)x[1] + z[1]; p = y[1] + p*p*v[1]; /* norm(x) */
1225 56 : if (fun(E, u, x, p)) break;
1226 : }
1227 : }
1228 :
1229 : void
1230 0 : forqfvec(void *E, long (*fun)(void *, GEN, GEN, double), GEN a, GEN BORNE)
1231 : {
1232 0 : pari_sp av = avma;
1233 : struct qfvec qv;
1234 0 : forqfvec_init(&qv, a);
1235 0 : forqfvec_i(E, fun, &qv, BORNE);
1236 0 : set_avma(av);
1237 0 : }
1238 :
1239 : struct qfvecwrap
1240 : {
1241 : void *E;
1242 : long (*fun)(void *, GEN);
1243 : };
1244 :
1245 : static long
1246 56 : forqfvec_wrap(void *E, GEN u, GEN x, double d)
1247 : {
1248 56 : pari_sp av = avma;
1249 56 : struct qfvecwrap *W = (struct qfvecwrap *) E;
1250 : (void) d;
1251 56 : return gc_long(av, W->fun(W->E, ZM_zc_mul_canon(u, x)));
1252 : }
1253 :
1254 : void
1255 42 : forqfvec1(void *E, long (*fun)(void *, GEN), GEN a, GEN BORNE)
1256 : {
1257 42 : pari_sp av = avma;
1258 : struct qfvecwrap wr;
1259 : struct qfvec qv;
1260 42 : wr.E = E; wr.fun = fun;
1261 42 : forqfvec_init(&qv, a);
1262 42 : forqfvec_i((void*) &wr, forqfvec_wrap, &qv, BORNE);
1263 35 : set_avma(av);
1264 35 : }
1265 :
1266 : void
1267 42 : forqfvec0(GEN a, GEN BORNE, GEN code)
1268 42 : { EXPRVOID_WRAP(code, forqfvec1(EXPR_ARGVOID, a, BORNE)) }
1269 :
1270 : enum { min_ALL = 0, min_FIRST, min_VECSMALL, min_VECSMALL2 };
1271 :
1272 : static int
1273 923 : stockmax_init(const char *fun, GEN STOCKMAX, long *maxrank)
1274 : {
1275 923 : long r = 200;
1276 923 : if (!STOCKMAX) { *maxrank = 200; return 1; }
1277 511 : STOCKMAX = gfloor(STOCKMAX);
1278 511 : if (typ(STOCKMAX) != t_INT) pari_err_TYPE(fun, STOCKMAX);
1279 511 : r = itos(STOCKMAX);
1280 511 : if (r < 0)
1281 : {
1282 0 : char *e = stack_strcat(fun, "[negative number of vectors]");
1283 0 : pari_err_TYPE(e, STOCKMAX);
1284 : }
1285 511 : *maxrank = r; return 0;
1286 : }
1287 :
1288 : /* Minimal vectors for the integral definite quadratic form: a.
1289 : * Result u:
1290 : * u[1]= Number of vectors of square norm <= BORNE
1291 : * u[2]= maximum norm found
1292 : * u[3]= list of vectors found (at most STOCKMAX, unless NULL)
1293 : *
1294 : * If BORNE = NULL: Minimal nonzero vectors.
1295 : * flag = min_ALL, as above
1296 : * flag = min_FIRST, exits when first suitable vector is found.
1297 : * flag = min_VECSMALL, return a t_VECSMALL of (half) the number of vectors
1298 : * for each norm
1299 : * flag = min_VECSMALL2, same but count only vectors with even norm, and shift
1300 : * the answer */
1301 : static GEN
1302 847 : minim0_dolll(GEN a, GEN BORNE, GEN STOCKMAX, long flag, long dolll)
1303 : {
1304 : GEN x, u, r, L, gnorme;
1305 847 : long n = lg(a)-1, i, j, k, s, maxrank, sBORNE;
1306 847 : pari_sp av = avma, av1;
1307 : double p,maxnorm,BOUND,*v,*y,*z,**q;
1308 847 : const double eps = 1e-10;
1309 : int stockall;
1310 : struct qfvec qv;
1311 :
1312 847 : if (!BORNE)
1313 56 : sBORNE = 0;
1314 : else
1315 : {
1316 791 : BORNE = gfloor(BORNE);
1317 791 : if (typ(BORNE) != t_INT) pari_err_TYPE("minim0",BORNE);
1318 791 : if (is_bigint(BORNE)) pari_err_PREC( "qfminim");
1319 790 : sBORNE = itos(BORNE); set_avma(av);
1320 790 : if (sBORNE < 0) sBORNE = 0;
1321 : }
1322 846 : stockall = stockmax_init("minim0", STOCKMAX, &maxrank);
1323 :
1324 846 : switch(flag)
1325 : {
1326 462 : case min_VECSMALL:
1327 : case min_VECSMALL2:
1328 462 : if (sBORNE <= 0) return cgetg(1, t_VECSMALL);
1329 434 : L = zero_zv(sBORNE);
1330 434 : if (flag == min_VECSMALL2) sBORNE <<= 1;
1331 434 : if (n == 0) return L;
1332 434 : break;
1333 35 : case min_FIRST:
1334 35 : if (n == 0 || (!sBORNE && BORNE)) return cgetg(1,t_VEC);
1335 21 : L = NULL; /* gcc -Wall */
1336 21 : break;
1337 349 : case min_ALL:
1338 349 : if (n == 0 || (!sBORNE && BORNE))
1339 14 : retmkvec3(gen_0, gen_0, cgetg(1, t_MAT));
1340 335 : L = new_chunk(1+maxrank);
1341 335 : break;
1342 0 : default:
1343 0 : return NULL;
1344 : }
1345 790 : minim_alloc(n+1, &q, &x, &y, &z, &v);
1346 :
1347 790 : forqfvec_init_dolll(&qv, &a, dolll);
1348 790 : av1 = avma;
1349 790 : r = qv.r;
1350 790 : u = qv.u;
1351 5912 : for (j=1; j<=n; j++)
1352 : {
1353 5122 : v[j] = rtodbl(gcoeff(r,j,j));
1354 29579 : for (i=1; i<j; i++) q[i][j] = rtodbl(gcoeff(r,i,j)); /* |.| <= 1/2 */
1355 : }
1356 :
1357 790 : if (sBORNE) maxnorm = 0.;
1358 : else
1359 : {
1360 56 : GEN B = gcoeff(a,1,1);
1361 56 : long t = 1;
1362 616 : for (i=2; i<=n; i++)
1363 : {
1364 560 : GEN c = gcoeff(a,i,i);
1365 560 : if (cmpii(c, B) < 0) { B = c; t = i; }
1366 : }
1367 56 : if (flag == min_FIRST) return gc_GEN(av, mkvec2(B, gel(u,t)));
1368 49 : maxnorm = -1.; /* don't update maxnorm */
1369 49 : if (is_bigint(B)) return NULL;
1370 48 : sBORNE = itos(B);
1371 : }
1372 782 : BOUND = sBORNE * (1 + eps);
1373 782 : if ((long)BOUND != sBORNE) return NULL;
1374 :
1375 770 : s = 0;
1376 770 : k = n; y[n] = z[n] = 0;
1377 770 : x[n] = (long)sqrt(BOUND/v[n]);
1378 1223264 : for(;;x[1]--)
1379 : {
1380 : do
1381 : {
1382 2245614 : if (k>1)
1383 : {
1384 1022259 : long l = k-1;
1385 1022259 : z[l] = 0;
1386 11756080 : for (j=k; j<=n; j++) z[l] += q[l][j]*x[j];
1387 1022259 : p = (double)x[k] + z[k];
1388 1022259 : y[l] = y[k] + p*p*v[k];
1389 1022259 : x[l] = (long)floor(sqrt((BOUND-y[l])/v[l])-z[l]);
1390 1022259 : k = l;
1391 : }
1392 : for(;;)
1393 : {
1394 3263729 : p = (double)x[k] + z[k];
1395 3263729 : if (y[k] + p*p*v[k] <= BOUND) break;
1396 1018115 : k++; x[k]--;
1397 : }
1398 : }
1399 2245614 : while (k > 1);
1400 1224034 : if (! x[1] && y[1]<=eps) break;
1401 :
1402 1223271 : p = (double)x[1] + z[1];
1403 1223271 : p = y[1] + p*p*v[1]; /* norm(x) */
1404 1223271 : if (maxnorm >= 0)
1405 : {
1406 1220723 : if (p > maxnorm) maxnorm = p;
1407 : }
1408 : else
1409 : { /* maxnorm < 0 : only look for minimal vectors */
1410 2548 : pari_sp av2 = avma;
1411 2548 : gnorme = roundr(dbltor(p));
1412 2548 : if (cmpis(gnorme, sBORNE) >= 0) set_avma(av2);
1413 : else
1414 : {
1415 14 : sBORNE = itos(gnorme); set_avma(av1);
1416 14 : BOUND = sBORNE * (1+eps);
1417 14 : L = new_chunk(maxrank+1);
1418 14 : s = 0;
1419 : }
1420 : }
1421 1223271 : s++;
1422 :
1423 1223271 : switch(flag)
1424 : {
1425 7 : case min_FIRST:
1426 7 : if (dolll) x = ZM_zc_mul_canon(u,x);
1427 7 : return gc_GEN(av, mkvec2(roundr(dbltor(p)), x));
1428 :
1429 248241 : case min_ALL:
1430 248241 : if (s > maxrank && stockall) /* overflow */
1431 : {
1432 490 : long maxranknew = maxrank << 1;
1433 490 : GEN Lnew = new_chunk(1 + maxranknew);
1434 344890 : for (i=1; i<=maxrank; i++) Lnew[i] = L[i];
1435 490 : L = Lnew; maxrank = maxranknew;
1436 : }
1437 248241 : if (s<=maxrank) gel(L,s) = leafcopy(x);
1438 248241 : break;
1439 :
1440 39200 : case min_VECSMALL:
1441 39200 : { ulong norm = (ulong)(p + 0.5); L[norm]++; }
1442 39200 : break;
1443 :
1444 935823 : case min_VECSMALL2:
1445 935823 : { ulong norm = (ulong)(p + 0.5); if (!odd(norm)) L[norm>>1]++; }
1446 935823 : break;
1447 :
1448 : }
1449 : }
1450 763 : switch(flag)
1451 : {
1452 7 : case min_FIRST:
1453 7 : retgc_const(av, cgetg(1, t_VEC));
1454 434 : case min_VECSMALL:
1455 : case min_VECSMALL2:
1456 434 : set_avma((pari_sp)L); return L;
1457 : }
1458 322 : r = (maxnorm >= 0) ? roundr(dbltor(maxnorm)): stoi(sBORNE);
1459 322 : k = minss(s,maxrank);
1460 322 : L[0] = evaltyp(t_MAT) | evallg(k + 1);
1461 322 : if (dolll)
1462 246092 : for (j=1; j<=k; j++)
1463 245805 : gel(L,j) = dolll==1 ? ZM_zc_mul_canon(u, gel(L,j))
1464 245805 : : ZM_zc_mul_canon_zm(u, gel(L,j));
1465 322 : return gc_GEN(av, mkvec3(stoi(s<<1), r, L));
1466 : }
1467 :
1468 : /* Closest vectors for the integral definite quadratic form: a.
1469 : * Code bases on minim0_dolll
1470 : * Result u:
1471 : * u[1]= Number of closest vectors of square distance <= BORNE
1472 : * u[2]= maximum squared distance found
1473 : * u[3]= list of vectors found (at most STOCKMAX, unless NULL)
1474 : *
1475 : * If BORNE = NULL or <= 0.: returns closest vectors.
1476 : * flag = min_ALL, as above
1477 : * flag = min_FIRST, exits when first suitable vector is found.
1478 : */
1479 : static GEN
1480 91 : cvp0_dolll(GEN a, GEN target, GEN BORNE, GEN STOCKMAX, long flag, long dolll)
1481 : {
1482 : GEN x, u, r, L;
1483 91 : long n = lg(a)-1, i, j, k, s, maxrank;
1484 91 : pari_sp av = avma, av1;
1485 : double p,maxnorm,BOUND,*v,*y,*z,*tt,**q, *tpre, sBORNE;
1486 91 : const double eps = 1e-10;
1487 : int stockall;
1488 : struct qfvec qv;
1489 91 : int done = 0;
1490 :
1491 91 : if (!is_vec_t(typ(target))) pari_err_TYPE("cvp0",target);
1492 91 : if (n != lg(target)-1) pari_err_TYPE("cvp0 [different dimensions]",target);
1493 77 : if (!BORNE)
1494 0 : sBORNE = 0.;
1495 : else
1496 : {
1497 77 : if (!is_real_t(typ(BORNE))) pari_err_TYPE("cvp0",BORNE);
1498 77 : sBORNE = gtodouble(BORNE);
1499 77 : if (sBORNE < 0.) sBORNE = 0.;
1500 : }
1501 77 : stockall = stockmax_init("cvp0", STOCKMAX, &maxrank);
1502 :
1503 77 : L = (flag==min_ALL) ? new_chunk(1+maxrank) : NULL;
1504 77 : if (n == 0)
1505 : {
1506 7 : if (flag==min_ALL) retmkvec3(gen_0, gen_0, cgetg(1, t_MAT));
1507 0 : return cgetg(1,t_VEC);
1508 : }
1509 :
1510 70 : minim_alloc(n+1, &q, &x, &y, &z, &v);
1511 70 : cvp_alloc(n+1, &tt, &tpre);
1512 :
1513 70 : forqfvec_init_dolll(&qv, &a, dolll);
1514 70 : av1 = avma;
1515 70 : r = qv.r;
1516 70 : u = qv.u;
1517 392 : for (j=1; j<=n; j++)
1518 : {
1519 322 : v[j] = rtodbl(gcoeff(r,j,j));
1520 1729 : for (i=1; i<j; i++) q[i][j] = rtodbl(gcoeff(r,i,j)); /* |.| <= 1/2 */
1521 : }
1522 :
1523 70 : if (dolll)
1524 : {
1525 70 : GEN tv = RgM_RgC_mul(ZM_inv(u, NULL), target);
1526 392 : for (j=1; j<=n; j++) tt[j] = gtodouble(gel(tv, j));
1527 : } else
1528 0 : for (j=1; j<=n; j++) tt[j] = gtodouble(gel(target, j));
1529 : /* precompute contribution of tt to z[l] */
1530 392 : for(k=1; k <= n; k++)
1531 : {
1532 322 : tpre[k] = -tt[k];
1533 1729 : for(j=k+1; j<=n; j++) tpre[k] -= q[k][j] * tt[j];
1534 : }
1535 :
1536 70 : if (sBORNE) maxnorm = 0.;
1537 : else
1538 : {
1539 28 : GEN B = gcoeff(a,1,1);
1540 112 : for (i = 2; i <= n; i++) B = addii(B, gcoeff(a,i,i));
1541 28 : maxnorm = -1.; /* don't update maxnorm */
1542 28 : if (is_bigint(B)) return NULL;
1543 28 : sBORNE = 0.;
1544 140 : for(i=1; i<=n; i++) sBORNE += v[i];
1545 : }
1546 70 : BOUND = sBORNE * (1 + eps);
1547 :
1548 70 : s = 0;
1549 70 : k = n; y[n] = 0;
1550 70 : z[n] = tpre[n];
1551 70 : x[n] = (long)floor(sqrt(BOUND/v[n])-z[n]);
1552 889 : for(;;x[1]--)
1553 : {
1554 : do
1555 : {
1556 8582 : if (k>1)
1557 : {
1558 7665 : long l = k-1;
1559 7665 : z[l] = tpre[l];
1560 61488 : for (j=k; j<=n; j++) z[l] += q[l][j]*x[j];
1561 7665 : p = (double)x[k] + z[k];
1562 7665 : y[l] = y[k] + p*p*v[k];
1563 7665 : x[l] = (long)floor(sqrt((BOUND-y[l])/v[l])-z[l]);
1564 7665 : k = l;
1565 : }
1566 : for(;;)
1567 : {
1568 16247 : p = (double)x[k] + z[k];
1569 16247 : if (y[k] + p*p*v[k] <= BOUND) break;
1570 7735 : if (k >= n) { done = 1; break; }
1571 7665 : k++; x[k]--;
1572 : }
1573 : }
1574 8582 : while (k > 1 && !done);
1575 959 : if (done) break;
1576 :
1577 889 : p = (double)x[1] + z[1];
1578 889 : p = y[1] + p*p*v[1]; /* norm(x-target) */
1579 889 : if (maxnorm >= 0)
1580 : {
1581 175 : if (p > maxnorm) maxnorm = p;
1582 : }
1583 : else
1584 : { /* maxnorm < 0 : only look for closest vectors */
1585 714 : if (p * (1+10*eps) < sBORNE) {
1586 231 : sBORNE = p; set_avma(av1);
1587 231 : BOUND = sBORNE * (1+eps);
1588 231 : L = new_chunk(maxrank+1);
1589 231 : s = 0;
1590 : }
1591 : }
1592 889 : s++;
1593 :
1594 889 : switch(flag)
1595 : {
1596 0 : case min_FIRST:
1597 0 : if (dolll) x = ZM_zc_mul(u,x);
1598 0 : return gc_GEN(av, mkvec2(dbltor(p), x));
1599 :
1600 889 : case min_ALL:
1601 889 : if (s > maxrank && stockall) /* overflow */
1602 : {
1603 0 : long maxranknew = maxrank << 1;
1604 0 : GEN Lnew = new_chunk(1 + maxranknew);
1605 0 : for (i=1; i<=maxrank; i++) Lnew[i] = L[i];
1606 0 : L = Lnew; maxrank = maxranknew;
1607 : }
1608 889 : if (s<=maxrank) gel(L,s) = leafcopy(x);
1609 889 : break;
1610 : }
1611 : }
1612 70 : switch(flag)
1613 : {
1614 0 : case min_FIRST:
1615 0 : retgc_const(av, cgetg(1, t_VEC));
1616 : }
1617 70 : r = (maxnorm >= 0) ? dbltor(maxnorm): dbltor(sBORNE);
1618 70 : k = minss(s,maxrank);
1619 70 : L[0] = evaltyp(t_MAT) | evallg(k + 1);
1620 322 : for (j=1; j<=k; j++)
1621 252 : gel(L,j) = dolll==1 ? ZM_zc_mul(u, gel(L,j))
1622 252 : : zc_to_ZC(gel(L,j));
1623 70 : return gc_GEN(av, mkvec3(stoi(s), r, L));
1624 : }
1625 :
1626 : static GEN
1627 553 : minim0(GEN a, GEN BORNE, GEN STOCKMAX, long flag)
1628 : {
1629 553 : GEN v = minim0_dolll(a, BORNE, STOCKMAX, flag, 1);
1630 552 : if (!v) pari_err_PREC("qfminim");
1631 546 : return v;
1632 : }
1633 :
1634 : static GEN
1635 91 : cvp0(GEN a, GEN target, GEN BORNE, GEN STOCKMAX, long flag)
1636 : {
1637 91 : GEN v = cvp0_dolll(a, target, BORNE, STOCKMAX, flag, 1);
1638 77 : if (!v) pari_err_PREC("qfcvp");
1639 77 : return v;
1640 : }
1641 :
1642 : static GEN
1643 252 : minim0_zm(GEN a, GEN BORNE, GEN STOCKMAX, long flag)
1644 : {
1645 252 : GEN v = minim0_dolll(a, BORNE, STOCKMAX, flag, 2);
1646 252 : if (!v) pari_err_PREC("qfminim");
1647 252 : return v;
1648 : }
1649 :
1650 : GEN
1651 462 : qfrep0(GEN a, GEN borne, long flag)
1652 462 : { return minim0(a, borne, gen_0, (flag & 1)? min_VECSMALL2: min_VECSMALL); }
1653 :
1654 : GEN
1655 133 : qfminim0(GEN a, GEN borne, GEN stockmax, long flag, long prec)
1656 : {
1657 133 : switch(flag)
1658 : {
1659 49 : case 0: return minim0(a,borne,stockmax,min_ALL);
1660 35 : case 1: return minim0(a,borne,gen_0 ,min_FIRST);
1661 49 : case 2:
1662 : {
1663 49 : long maxnum = -1;
1664 49 : if (typ(a) != t_MAT) pari_err_TYPE("qfminim",a);
1665 49 : if (stockmax) {
1666 14 : if (typ(stockmax) != t_INT) pari_err_TYPE("qfminim",stockmax);
1667 14 : maxnum = itos(stockmax);
1668 : }
1669 49 : a = fincke_pohst(a,borne,maxnum,prec,NULL);
1670 42 : if (!a) pari_err_PREC("qfminim");
1671 42 : return a;
1672 : }
1673 0 : default: pari_err_FLAG("qfminim");
1674 : }
1675 : return NULL; /* LCOV_EXCL_LINE */
1676 : }
1677 :
1678 :
1679 : GEN
1680 91 : qfcvp0(GEN a, GEN target, GEN borne, GEN stockmax, long flag)
1681 : {
1682 91 : switch(flag)
1683 : {
1684 91 : case 0: return cvp0(a,target,borne,stockmax,min_ALL);
1685 0 : case 1: return cvp0(a,target,borne,gen_0 ,min_FIRST);
1686 : /* case 2:
1687 : TODO: more robust finke_pohst enumeration */
1688 0 : default: pari_err_FLAG("qfcvp");
1689 : }
1690 : return NULL; /* LCOV_EXCL_LINE */
1691 : }
1692 :
1693 : GEN
1694 7 : minim(GEN a, GEN borne, GEN stockmax)
1695 7 : { return minim0(a,borne,stockmax,min_ALL); }
1696 :
1697 : GEN
1698 252 : minim_zm(GEN a, GEN borne, GEN stockmax)
1699 252 : { return minim0_zm(a,borne,stockmax,min_ALL); }
1700 :
1701 : GEN
1702 42 : minim_raw(GEN a, GEN BORNE, GEN STOCKMAX)
1703 42 : { return minim0_dolll(a, BORNE, STOCKMAX, min_ALL, 0); }
1704 :
1705 : GEN
1706 0 : minim2(GEN a, GEN borne, GEN stockmax)
1707 0 : { return minim0(a,borne,stockmax,min_FIRST); }
1708 :
1709 : /* If V depends linearly from the columns of the matrix, return 0.
1710 : * Otherwise, update INVP and L and return 1. No GC. */
1711 : static int
1712 1652 : addcolumntomatrix(GEN V, GEN invp, GEN L)
1713 : {
1714 1652 : long i,j,k, n = lg(invp);
1715 1652 : GEN a = cgetg(n, t_COL), ak = NULL, mak;
1716 :
1717 84231 : for (k = 1; k < n; k++)
1718 83706 : if (!L[k])
1719 : {
1720 27902 : ak = RgMrow_zc_mul(invp, V, k);
1721 27902 : if (!gequal0(ak)) break;
1722 : }
1723 1652 : if (k == n) return 0;
1724 1127 : L[k] = 1;
1725 1127 : mak = gneg_i(ak);
1726 43253 : for (i=k+1; i<n; i++)
1727 42126 : gel(a,i) = gdiv(RgMrow_zc_mul(invp, V, i), mak);
1728 43883 : for (j=1; j<=k; j++)
1729 : {
1730 42756 : GEN c = gel(invp,j), ck = gel(c,k);
1731 42756 : if (gequal0(ck)) continue;
1732 8757 : gel(c,k) = gdiv(ck, ak);
1733 8757 : if (j==k)
1734 43253 : for (i=k+1; i<n; i++)
1735 42126 : gel(c,i) = gmul(gel(a,i), ck);
1736 : else
1737 184814 : for (i=k+1; i<n; i++)
1738 177184 : gel(c,i) = gadd(gel(c,i), gmul(gel(a,i), ck));
1739 : }
1740 1127 : return 1;
1741 : }
1742 :
1743 : GEN
1744 42 : qfperfection(GEN a)
1745 : {
1746 42 : pari_sp av = avma;
1747 : GEN u, L;
1748 42 : long r, s, k, l, n = lg(a)-1;
1749 :
1750 42 : if (!n) return gen_0;
1751 42 : if (typ(a) != t_MAT || !RgM_is_ZM(a)) pari_err_TYPE("qfperfection",a);
1752 42 : a = minim_lll(a, &u);
1753 42 : L = minim_raw(a,NULL,NULL);
1754 42 : r = (n*(n+1)) >> 1;
1755 42 : if (L)
1756 : {
1757 : GEN D, V, invp;
1758 35 : L = gel(L, 3); l = lg(L);
1759 35 : if (l == 2) { set_avma(av); return gen_1; }
1760 : /* |L[i]|^2 fits into a long for all i */
1761 21 : D = zero_zv(r);
1762 21 : V = cgetg(r+1, t_VECSMALL);
1763 21 : invp = matid(r);
1764 21 : s = 0;
1765 1659 : for (k = 1; k < l; k++)
1766 : {
1767 1652 : pari_sp av2 = avma;
1768 1652 : GEN x = gel(L,k);
1769 : long i, j, I;
1770 21098 : for (i = I = 1; i<=n; i++)
1771 145278 : for (j=i; j<=n; j++,I++) V[I] = x[i]*x[j];
1772 1652 : if (!addcolumntomatrix(V,invp,D)) set_avma(av2);
1773 1127 : else if (++s == r) break;
1774 : }
1775 : }
1776 : else
1777 : {
1778 : GEN M;
1779 7 : L = fincke_pohst(a,NULL,-1, DEFAULTPREC, NULL);
1780 7 : if (!L) pari_err_PREC("qfminim");
1781 7 : L = gel(L, 3); l = lg(L);
1782 7 : if (l == 2) { set_avma(av); return gen_1; }
1783 7 : M = cgetg(l, t_MAT);
1784 959 : for (k = 1; k < l; k++)
1785 : {
1786 952 : GEN x = gel(L,k), c = cgetg(r+1, t_COL);
1787 : long i, I, j;
1788 16184 : for (i = I = 1; i<=n; i++)
1789 144704 : for (j=i; j<=n; j++,I++) gel(c,I) = mulii(gel(x,i), gel(x,j));
1790 952 : gel(M,k) = c;
1791 : }
1792 7 : s = ZM_rank(M);
1793 : }
1794 28 : return gc_utoipos(av, s);
1795 : }
1796 :
1797 : static GEN
1798 141 : clonefill(GEN S, long s, long t)
1799 : { /* initialize to dummy values */
1800 141 : GEN T = S, dummy = cgetg(1, t_STR);
1801 : long i;
1802 310917 : for (i = s+1; i <= t; i++) gel(S,i) = dummy;
1803 141 : S = gclone(S); if (isclone(T)) gunclone(T);
1804 141 : return S;
1805 : }
1806 :
1807 : /* increment ZV x, by incrementing cell of index k. Initial value x0[k] was
1808 : * chosen to minimize qf(x) for given x0[1..k-1] and x0[k+1,..] = 0
1809 : * The last nonzero entry must be positive and goes through x0[k]+1,2,3,...
1810 : * Others entries go through: x0[k]+1,-1,2,-2,...*/
1811 : INLINE void
1812 2952918 : step(GEN x, GEN y, GEN inc, long k)
1813 : {
1814 2952918 : if (!signe(gel(y,k))) /* x[k+1..] = 0 */
1815 160814 : gel(x,k) = addiu(gel(x,k), 1); /* leading coeff > 0 */
1816 : else
1817 : {
1818 2792104 : long i = inc[k];
1819 2792104 : gel(x,k) = addis(gel(x,k), i),
1820 2792120 : inc[k] = (i > 0)? -1-i: 1-i;
1821 : }
1822 2952933 : }
1823 :
1824 : /* 1 if we are "sure" that x < y, up to few rounding errors, i.e.
1825 : * x < y - epsilon. More precisely :
1826 : * - sign(x - y) < 0
1827 : * - lgprec(x-y) > 3 || expo(x - y) - expo(x) > -24 */
1828 : static int
1829 1216469 : mplessthan(GEN x, GEN y)
1830 : {
1831 1216469 : pari_sp av = avma;
1832 1216469 : GEN z = mpsub(x, y);
1833 1216468 : set_avma(av);
1834 1216467 : if (typ(z) == t_INT) return (signe(z) < 0);
1835 1216467 : if (signe(z) >= 0) return 0;
1836 22160 : if (realprec(z) > LOWDEFAULTPREC) return 1;
1837 22160 : return ( expo(z) - mpexpo(x) > -24 );
1838 : }
1839 :
1840 : /* 1 if we are "sure" that x > y, up to few rounding errors, i.e.
1841 : * x > y + epsilon */
1842 : static int
1843 4621599 : mpgreaterthan(GEN x, GEN y)
1844 : {
1845 4621599 : pari_sp av = avma;
1846 4621599 : GEN z = mpsub(x, y);
1847 4621621 : set_avma(av);
1848 4621645 : if (typ(z) == t_INT) return (signe(z) > 0);
1849 4621645 : if (signe(z) <= 0) return 0;
1850 2689971 : if (realprec(z) > LOWDEFAULTPREC) return 1;
1851 476044 : return ( expo(z) - mpexpo(x) > -24 );
1852 : }
1853 :
1854 : /* x a t_INT, y t_INT or t_REAL */
1855 : INLINE GEN
1856 1228563 : mulimp(GEN x, GEN y)
1857 : {
1858 1228563 : if (typ(y) == t_INT) return mulii(x,y);
1859 1228563 : return signe(x) ? mulir(x,y): gen_0;
1860 : }
1861 : /* x + y*z, x,z two mp's, y a t_INT */
1862 : INLINE GEN
1863 13538834 : addmulimp(GEN x, GEN y, GEN z)
1864 : {
1865 13538834 : if (!signe(y)) return x;
1866 5831089 : if (typ(z) == t_INT) return mpadd(x, mulii(y, z));
1867 5831089 : return mpadd(x, mulir(y, z));
1868 : }
1869 :
1870 : /* yk + vk * (xk + zk)^2 */
1871 : static GEN
1872 5780587 : norm_aux(GEN xk, GEN yk, GEN zk, GEN vk)
1873 : {
1874 5780587 : GEN t = mpadd(xk, zk);
1875 5780574 : if (typ(t) == t_INT) { /* probably gen_0, avoid loss of accuracy */
1876 306221 : yk = addmulimp(yk, sqri(t), vk);
1877 : } else {
1878 5474353 : yk = mpadd(yk, mpmul(sqrr(t), vk));
1879 : }
1880 5780527 : return yk;
1881 : }
1882 : /* yk + vk * (xk + zk)^2 < B + epsilon */
1883 : static int
1884 4167491 : check_bound(GEN B, GEN xk, GEN yk, GEN zk, GEN vk)
1885 : {
1886 4167491 : pari_sp av = avma;
1887 4167491 : int f = mpgreaterthan(norm_aux(xk,yk,zk,vk), B);
1888 4167486 : return gc_bool(av, !f);
1889 : }
1890 :
1891 : /* q(k-th canonical basis vector), where q is given in Cholesky form
1892 : * q(x) = sum_{i = 1}^n q[i,i] (x[i] + sum_{j > i} q[i,j] x[j])^2.
1893 : * Namely q(e_k) = q[k,k] + sum_{i < k} q[i,i] q[i,k]^2
1894 : * Assume 1 <= k <= n. */
1895 : static GEN
1896 182 : cholesky_norm_ek(GEN q, long k)
1897 : {
1898 182 : GEN t = gcoeff(q,k,k);
1899 : long i;
1900 1484 : for (i = 1; i < k; i++) t = norm_aux(gen_0, t, gcoeff(q,i,k), gcoeff(q,i,i));
1901 182 : return t;
1902 : }
1903 :
1904 : /* q is the Cholesky decomposition of a quadratic form
1905 : * Enumerate vectors whose norm is less than BORNE (Algo 2.5.7),
1906 : * minimal vectors if BORNE = NULL (implies check = NULL).
1907 : * If (check != NULL) consider only vectors passing the check, and assumes
1908 : * we only want the smallest possible vectors */
1909 : static GEN
1910 14713 : smallvectors(GEN q, GEN BORNE, long maxnum, FP_chk_fun *CHECK)
1911 : {
1912 14713 : long N = lg(q), n = N-1, i, j, k, s, stockmax, checkcnt = 1;
1913 : pari_sp av, av1;
1914 : GEN inc, S, x, y, z, v, p1, alpha, norms;
1915 : GEN norme1, normax1, borne1, borne2;
1916 14713 : GEN (*check)(void *,GEN) = CHECK? CHECK->f: NULL;
1917 14713 : void *data = CHECK? CHECK->data: NULL;
1918 14713 : const long skipfirst = CHECK? CHECK->skipfirst: 0;
1919 14713 : const int stockall = (maxnum == -1);
1920 :
1921 14713 : alpha = dbltor(0.95);
1922 14713 : normax1 = gen_0;
1923 :
1924 14713 : v = cgetg(N,t_VEC);
1925 14713 : inc = const_vecsmall(n, 1);
1926 :
1927 14713 : av = avma;
1928 14713 : stockmax = stockall? 2000: maxnum;
1929 14713 : norms = cgetg(check?(stockmax+1): 1,t_VEC); /* unused if (!check) */
1930 14713 : S = cgetg(stockmax+1,t_VEC);
1931 14713 : x = cgetg(N,t_COL);
1932 14713 : y = cgetg(N,t_COL);
1933 14713 : z = cgetg(N,t_COL);
1934 97807 : for (i=1; i<N; i++) {
1935 83094 : gel(v,i) = gcoeff(q,i,i);
1936 83094 : gel(x,i) = gel(y,i) = gel(z,i) = gen_0;
1937 : }
1938 14713 : if (BORNE)
1939 : {
1940 14692 : borne1 = BORNE;
1941 14692 : if (gsigne(borne1) <= 0) retmkvec3(gen_0, gen_0, cgetg(1,t_MAT));
1942 14678 : if (typ(borne1) != t_REAL)
1943 : {
1944 : long prec;
1945 419 : prec = nbits2prec(gexpo(borne1) + 10);
1946 419 : borne1 = gtofp(borne1, maxss(prec, DEFAULTPREC));
1947 : }
1948 : }
1949 : else
1950 : {
1951 21 : borne1 = gcoeff(q,1,1);
1952 203 : for (i=2; i<N; i++)
1953 : {
1954 182 : GEN b = cholesky_norm_ek(q, i);
1955 182 : if (gcmp(b, borne1) < 0) borne1 = b;
1956 : }
1957 : /* borne1 = norm of smallest basis vector */
1958 : }
1959 14699 : borne2 = mulrr(borne1,alpha);
1960 14699 : if (DEBUGLEVEL>2)
1961 0 : err_printf("smallvectors looking for norm < %P.4G\n",borne1);
1962 14699 : s = 0; k = n;
1963 383967 : for(;; step(x,y,inc,k)) /* main */
1964 : { /* x (supposedly) small vector, ZV.
1965 : * For all t >= k, we have
1966 : * z[t] = sum_{j > t} q[t,j] * x[j]
1967 : * y[t] = sum_{i > t} q[i,i] * (x[i] + z[i])^2
1968 : * = 0 <=> x[i]=0 for all i>t */
1969 : do
1970 : {
1971 1612527 : int skip = 0;
1972 1612527 : if (k > 1)
1973 : {
1974 1228563 : long l = k-1;
1975 1228563 : av1 = avma;
1976 1228563 : p1 = mulimp(gel(x,k), gcoeff(q,l,k));
1977 14461190 : for (j=k+1; j<N; j++) p1 = addmulimp(p1, gel(x,j), gcoeff(q,l,j));
1978 1228559 : gel(z,l) = gc_leaf(av1,p1);
1979 :
1980 1228564 : av1 = avma;
1981 1228564 : p1 = norm_aux(gel(x,k), gel(y,k), gel(z,k), gel(v,k));
1982 1228562 : gel(y,l) = gc_leaf(av1, p1);
1983 : /* skip the [x_1,...,x_skipfirst,0,...,0] */
1984 1228563 : if ((l <= skipfirst && !signe(gel(y,skipfirst)))
1985 1228563 : || mplessthan(borne1, gel(y,l))) skip = 1;
1986 : else /* initial value, minimizing (x[l] + z[l])^2, hence qf(x) for
1987 : the given x[1..l-1] */
1988 1214577 : gel(x,l) = mpround( mpneg(gel(z,l)) );
1989 1228563 : k = l;
1990 : }
1991 1228558 : for(;; step(x,y,inc,k))
1992 : { /* at most 2n loops */
1993 2841085 : if (!skip)
1994 : {
1995 2827100 : if (check_bound(borne1, gel(x,k),gel(y,k),gel(z,k),gel(v,k))) break;
1996 1340396 : step(x,y,inc,k);
1997 1340410 : if (check_bound(borne1, gel(x,k),gel(y,k),gel(z,k),gel(v,k))) break;
1998 : }
1999 1243257 : skip = 0; inc[k] = 1;
2000 1243257 : if (++k > n) goto END;
2001 : }
2002 :
2003 1597838 : if (gc_needed(av,2))
2004 : {
2005 15 : if(DEBUGMEM>1) pari_warn(warnmem,"smallvectors");
2006 15 : if (stockmax) S = clonefill(S, s, stockmax);
2007 15 : if (check) {
2008 15 : GEN dummy = cgetg(1, t_STR);
2009 9629 : for (i=s+1; i<=stockmax; i++) gel(norms,i) = dummy;
2010 : }
2011 15 : (void)gc_all(av,7,&x,&y,&z,&normax1,&borne1,&borne2,&norms);
2012 : }
2013 : }
2014 1597838 : while (k > 1);
2015 383969 : if (!signe(gel(x,1)) && !signe(gel(y,1))) continue; /* exclude 0 */
2016 :
2017 383254 : av1 = avma;
2018 383254 : norme1 = norm_aux(gel(x,1),gel(y,1),gel(z,1),gel(v,1));
2019 383252 : if (mpgreaterthan(norme1,borne1)) { set_avma(av1); continue; /* main */ }
2020 :
2021 383251 : norme1 = gc_leaf(av1,norme1);
2022 383251 : if (check)
2023 : {
2024 314665 : if (checkcnt < 5 && mpcmp(norme1, borne2) < 0)
2025 : {
2026 4416 : if (!check(data,x)) { checkcnt++ ; continue; /* main */}
2027 476 : if (DEBUGLEVEL>4) err_printf("New bound: %Ps", norme1);
2028 476 : borne1 = norme1;
2029 476 : borne2 = mulrr(borne1, alpha);
2030 476 : s = 0; checkcnt = 0;
2031 : }
2032 : }
2033 : else
2034 : {
2035 68586 : if (!BORNE) /* find minimal vectors */
2036 : {
2037 1890 : if (mplessthan(norme1, borne1))
2038 : { /* strictly smaller vector than previously known */
2039 0 : borne1 = norme1; /* + epsilon */
2040 0 : s = 0;
2041 : }
2042 : }
2043 : else
2044 66696 : if (mpcmp(norme1,normax1) > 0) normax1 = norme1;
2045 : }
2046 379311 : if (++s > stockmax) continue; /* too many vectors: no longer remember */
2047 378380 : if (check) gel(norms,s) = norme1;
2048 378380 : gel(S,s) = leafcopy(x);
2049 378381 : if (s != stockmax) continue; /* still room, get next vector */
2050 :
2051 126 : if (check)
2052 : { /* overflow, eliminate vectors failing "check" */
2053 105 : pari_sp av2 = avma;
2054 : long imin, imax;
2055 105 : GEN per = indexsort(norms), S2 = cgetg(stockmax+1, t_VEC);
2056 105 : if (DEBUGLEVEL>2) err_printf("sorting... [%ld elts]\n",s);
2057 : /* let N be the minimal norm so far for x satisfying 'check'. Keep
2058 : * all elements of norm N */
2059 26593 : for (i = 1; i <= s; i++)
2060 : {
2061 26586 : long k = per[i];
2062 26586 : if (check(data,gel(S,k))) { borne1 = gel(norms,k); break; }
2063 : }
2064 105 : imin = i;
2065 20943 : for (; i <= s; i++)
2066 20922 : if (mpgreaterthan(gel(norms,per[i]), borne1)) break;
2067 105 : imax = i;
2068 20943 : for (i=imin, s=0; i < imax; i++) gel(S2,++s) = gel(S,per[i]);
2069 20943 : for (i = 1; i <= s; i++) gel(S,i) = gel(S2,i);
2070 105 : set_avma(av2);
2071 105 : if (s) { borne2 = mulrr(borne1, alpha); checkcnt = 0; }
2072 105 : if (!stockall) continue;
2073 105 : if (s > stockmax/2) stockmax <<= 1;
2074 105 : norms = cgetg(stockmax+1, t_VEC);
2075 20943 : for (i = 1; i <= s; i++) gel(norms,i) = borne1;
2076 : }
2077 : else
2078 : {
2079 21 : if (!stockall && BORNE) goto END;
2080 21 : if (!stockall) continue;
2081 21 : stockmax <<= 1;
2082 : }
2083 :
2084 : {
2085 126 : GEN Snew = clonefill(vec_lengthen(S,stockmax), s, stockmax);
2086 126 : if (isclone(S)) gunclone(S);
2087 126 : S = Snew;
2088 : }
2089 : }
2090 14699 : END:
2091 14699 : if (s < stockmax) stockmax = s;
2092 14699 : if (check)
2093 : {
2094 : GEN per, alph, pols, p;
2095 14671 : if (DEBUGLEVEL>2) err_printf("final sort & check...\n");
2096 14671 : setlg(norms,stockmax+1); per = indexsort(norms);
2097 14671 : alph = cgetg(stockmax+1,t_VEC);
2098 14671 : pols = cgetg(stockmax+1,t_VEC);
2099 84520 : for (j=0,i=1; i<=stockmax; i++)
2100 : {
2101 70117 : long t = per[i];
2102 70117 : GEN N = gel(norms,t);
2103 70117 : if (j && mpgreaterthan(N, borne1)) break;
2104 69849 : if ((p = check(data,gel(S,t))))
2105 : {
2106 55908 : if (!j) borne1 = N;
2107 55908 : j++;
2108 55908 : gel(pols,j) = p;
2109 55908 : gel(alph,j) = gel(S,t);
2110 : }
2111 : }
2112 14671 : setlg(pols,j+1);
2113 14671 : setlg(alph,j+1);
2114 14671 : if (stockmax && isclone(S)) { alph = gcopy(alph); gunclone(S); }
2115 14671 : return mkvec2(pols, alph);
2116 : }
2117 28 : if (stockmax)
2118 : {
2119 21 : setlg(S,stockmax+1);
2120 21 : settyp(S,t_MAT);
2121 21 : if (isclone(S)) { p1 = S; S = gcopy(S); gunclone(p1); }
2122 : }
2123 : else
2124 7 : S = cgetg(1,t_MAT);
2125 28 : return mkvec3(utoi(s<<1), borne1, S);
2126 : }
2127 :
2128 : /* solve q(x) = x~.a.x <= bound, a > 0.
2129 : * If check is non-NULL keep x only if check(x).
2130 : * If a is a vector, assume a[1] is the LLL-reduced Cholesky form of q */
2131 : GEN
2132 14734 : fincke_pohst(GEN a, GEN B0, long stockmax, long PREC, FP_chk_fun *CHECK)
2133 : {
2134 14734 : pari_sp av = avma;
2135 : VOLATILE long i,j,l;
2136 14734 : VOLATILE GEN r,rinv,rinvtrans,u,v,res,z,vnorm,rperm,perm,uperm, bound = B0;
2137 :
2138 14734 : if (typ(a) == t_VEC)
2139 : {
2140 14266 : r = gel(a,1);
2141 14266 : u = NULL;
2142 : }
2143 : else
2144 : {
2145 468 : long prec = PREC;
2146 468 : l = lg(a);
2147 468 : if (l == 1)
2148 : {
2149 7 : if (CHECK) pari_err_TYPE("fincke_pohst [dimension 0]", a);
2150 7 : retmkvec3(gen_0, gen_0, cgetg(1,t_MAT));
2151 : }
2152 461 : u = lllfp(a, 0.75, LLL_GRAM | LLL_IM);
2153 454 : if (!u || lg(u) != lg(a)) return gc_NULL(av);
2154 454 : r = qf_RgM_apply(a,u);
2155 454 : i = gprecision(r);
2156 454 : if (i)
2157 412 : prec = i;
2158 : else {
2159 42 : prec = DEFAULTPREC + nbits2extraprec(gexpo(r));
2160 42 : if (prec < PREC) prec = PREC;
2161 : }
2162 454 : if (DEBUGLEVEL>2) err_printf("first LLL: prec = %ld\n", prec);
2163 454 : r = qfgaussred_positive(r);
2164 454 : if (!r) return gc_NULL(av);
2165 1984 : for (i=1; i<l; i++)
2166 : {
2167 1530 : GEN s = gsqrt(gcoeff(r,i,i), prec);
2168 1530 : gcoeff(r,i,i) = s;
2169 4236 : for (j=i+1; j<l; j++) gcoeff(r,i,j) = gmul(s, gcoeff(r,i,j));
2170 : }
2171 : }
2172 : /* now r~ * r = a in LLL basis */
2173 14720 : rinv = RgM_inv_upper(r);
2174 14720 : if (!rinv) return gc_NULL(av);
2175 14720 : rinvtrans = shallowtrans(rinv);
2176 14720 : if (DEBUGLEVEL>2)
2177 0 : err_printf("Fincke-Pohst, final LLL: prec = %ld\n", gprecision(rinvtrans));
2178 14720 : v = lll(rinvtrans);
2179 14720 : if (lg(v) != lg(rinvtrans)) return gc_NULL(av);
2180 :
2181 14720 : rinvtrans = RgM_mul(rinvtrans, v);
2182 14720 : v = ZM_inv(shallowtrans(v),NULL);
2183 14720 : r = RgM_mul(r,v);
2184 14720 : u = u? ZM_mul(u,v): v;
2185 :
2186 14720 : l = lg(r);
2187 14720 : vnorm = cgetg(l,t_VEC);
2188 97841 : for (j=1; j<l; j++) gel(vnorm,j) = gnorml2(gel(rinvtrans,j));
2189 14720 : rperm = cgetg(l,t_MAT);
2190 14720 : uperm = cgetg(l,t_MAT); perm = indexsort(vnorm);
2191 97842 : for (i=1; i<l; i++) { uperm[l-i] = u[perm[i]]; rperm[l-i] = r[perm[i]]; }
2192 14720 : u = uperm;
2193 14720 : r = rperm; res = NULL;
2194 14720 : pari_CATCH(e_PREC) { }
2195 : pari_TRY {
2196 : GEN q;
2197 14720 : if (CHECK && CHECK->f_init) bound = CHECK->f_init(CHECK, r, u);
2198 14713 : q = gaussred_from_QR(r, gprecision(vnorm));
2199 14713 : if (q) res = smallvectors(q, bound, stockmax, CHECK);
2200 14713 : } pari_ENDCATCH;
2201 14720 : if (!res) return gc_NULL(av);
2202 14713 : if (CHECK)
2203 : {
2204 14671 : if (CHECK->f_post) res = CHECK->f_post(CHECK, res, u);
2205 14671 : return res;
2206 : }
2207 :
2208 42 : z = cgetg(4,t_VEC);
2209 42 : gel(z,1) = gcopy(gel(res,1));
2210 42 : gel(z,2) = gcopy(gel(res,2));
2211 42 : gel(z,3) = ZM_mul(u, gel(res,3)); return gc_upto(av,z);
2212 : }
|