Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
Download

Testing latest pari + WASM + node.js... and it works?! Wow.

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