Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
Download

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

28485 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
/** LINEAR ALGEBRA **/
18
/** (second part) **/
19
/** **/
20
/********************************************************************/
21
#include "pari.h"
22
#include "paripriv.h"
23
24
#define DEBUGLEVEL DEBUGLEVEL_mat
25
26
/*******************************************************************/
27
/* */
28
/* CHARACTERISTIC POLYNOMIAL */
29
/* */
30
/*******************************************************************/
31
32
static GEN
33
Flm_charpoly_i(GEN x, ulong p)
34
{
35
long lx = lg(x), r, i;
36
GEN H, y = cgetg(lx+1, t_VEC);
37
gel(y,1) = pol1_Flx(0); H = Flm_hess(x, p);
38
for (r = 1; r < lx; r++)
39
{
40
pari_sp av2 = avma;
41
ulong a = 1;
42
GEN z, b = zero_Flx(0);
43
for (i = r-1; i; i--)
44
{
45
a = Fl_mul(a, ucoeff(H,i+1,i), p);
46
if (!a) break;
47
b = Flx_add(b, Flx_Fl_mul(gel(y,i), Fl_mul(a,ucoeff(H,i,r),p), p), p);
48
}
49
z = Flx_sub(Flx_shift(gel(y,r), 1),
50
Flx_Fl_mul(gel(y,r), ucoeff(H,r,r), p), p);
51
/* (X - H[r,r])y[r] - b */
52
gel(y,r+1) = gerepileuptoleaf(av2, Flx_sub(z, b, p));
53
}
54
return gel(y,lx);
55
}
56
57
GEN
58
Flm_charpoly(GEN x, ulong p)
59
{
60
pari_sp av = avma;
61
return gerepileuptoleaf(av, Flm_charpoly_i(x,p));
62
}
63
64
GEN
65
FpM_charpoly(GEN x, GEN p)
66
{
67
pari_sp av = avma;
68
long lx, r, i;
69
GEN y, H;
70
71
if (lgefint(p) == 3)
72
{
73
ulong pp = p[2];
74
y = Flx_to_ZX(Flm_charpoly_i(ZM_to_Flm(x,pp), pp));
75
return gerepileupto(av, y);
76
}
77
lx = lg(x); y = cgetg(lx+1, t_VEC);
78
gel(y,1) = pol_1(0); H = FpM_hess(x, p);
79
for (r = 1; r < lx; r++)
80
{
81
pari_sp av2 = avma;
82
GEN z, a = gen_1, b = pol_0(0);
83
for (i = r-1; i; i--)
84
{
85
a = Fp_mul(a, gcoeff(H,i+1,i), p);
86
if (!signe(a)) break;
87
b = ZX_add(b, ZX_Z_mul(gel(y,i), Fp_mul(a,gcoeff(H,i,r),p)));
88
}
89
b = FpX_red(b, p);
90
z = FpX_sub(RgX_shift_shallow(gel(y,r), 1),
91
FpX_Fp_mul(gel(y,r), gcoeff(H,r,r), p), p);
92
z = FpX_sub(z,b,p);
93
if (r+1 == lx) { gel(y,lx) = z; break; }
94
gel(y,r+1) = gerepileupto(av2, z); /* (X - H[r,r])y[r] - b */
95
}
96
return gerepileupto(av, gel(y,lx));
97
}
98
99
GEN
100
charpoly0(GEN x, long v, long flag)
101
{
102
if (v<0) v = 0;
103
switch(flag)
104
{
105
case 0: return caradj(x,v,NULL);
106
case 1: return caract(x,v);
107
case 2: return carhess(x,v);
108
case 3: return carberkowitz(x,v);
109
case 4:
110
if (typ(x) != t_MAT) pari_err_TYPE("charpoly",x);
111
RgM_check_ZM(x, "charpoly");
112
x = ZM_charpoly(x); setvarn(x, v); return x;
113
case 5:
114
return charpoly(x, v);
115
}
116
pari_err_FLAG("charpoly");
117
return NULL; /* LCOV_EXCL_LINE */
118
}
119
120
/* characteristic pol. Easy cases. Return NULL in case it's not so easy. */
121
static GEN
122
easychar(GEN x, long v)
123
{
124
pari_sp av;
125
long lx;
126
GEN p1;
127
128
switch(typ(x))
129
{
130
case t_INT: case t_REAL: case t_INTMOD:
131
case t_FRAC: case t_PADIC:
132
p1=cgetg(4,t_POL);
133
p1[1]=evalsigne(1) | evalvarn(v);
134
gel(p1,2) = gneg(x); gel(p1,3) = gen_1;
135
return p1;
136
137
case t_COMPLEX: case t_QUAD:
138
p1 = cgetg(5,t_POL);
139
p1[1] = evalsigne(1) | evalvarn(v);
140
gel(p1,2) = gnorm(x); av = avma;
141
gel(p1,3) = gerepileupto(av, gneg(gtrace(x)));
142
gel(p1,4) = gen_1; return p1;
143
144
case t_FFELT: {
145
pari_sp ltop=avma;
146
p1 = FpX_to_mod(FF_charpoly(x), FF_p_i(x));
147
setvarn(p1,v); return gerepileupto(ltop,p1);
148
}
149
150
case t_POLMOD:
151
{
152
GEN A = gel(x,2), T = gel(x,1);
153
if (typ(A)==t_POL && RgX_is_QX(A) && RgX_is_ZX(T))
154
return QXQ_charpoly(A, T, v);
155
else
156
return RgXQ_charpoly(A, T, v);
157
}
158
case t_MAT:
159
lx=lg(x);
160
if (lx==1) return pol_1(v);
161
if (lgcols(x) != lx) break;
162
return NULL;
163
}
164
pari_err_TYPE("easychar",x);
165
return NULL; /* LCOV_EXCL_LINE */
166
}
167
/* compute charpoly by mapping to Fp first, return lift to Z */
168
static GEN
169
RgM_Fp_charpoly(GEN x, GEN p, long v)
170
{
171
GEN T;
172
if (lgefint(p) == 3)
173
{
174
ulong pp = itou(p);
175
T = Flm_charpoly_i(RgM_to_Flm(x, pp), pp);
176
T = Flx_to_ZX(T);
177
}
178
else
179
T = FpM_charpoly(RgM_to_FpM(x, p), p);
180
setvarn(T, v); return T;
181
}
182
GEN
183
charpoly(GEN x, long v)
184
{
185
GEN T, p = NULL;
186
long prec;
187
if ((T = easychar(x,v))) return T;
188
switch(RgM_type(x, &p,&T,&prec))
189
{
190
case t_INT:
191
T = ZM_charpoly(x); setvarn(T, v); break;
192
case t_INTMOD:
193
if (!BPSW_psp(p)) T = carberkowitz(x, v);
194
else
195
{
196
pari_sp av = avma;
197
T = RgM_Fp_charpoly(x,p,v);
198
T = gerepileupto(av, FpX_to_mod(T,p));
199
}
200
break;
201
case t_REAL:
202
case t_COMPLEX:
203
case t_PADIC:
204
T = carhess(x, v);
205
break;
206
default:
207
T = carberkowitz(x, v);
208
}
209
return T;
210
}
211
212
/* We possibly worked with an "invalid" polynomial p, satisfying
213
* varn(p) > gvar2(p). Fix this. */
214
static GEN
215
fix_pol(pari_sp av, GEN p)
216
{
217
long w = gvar2(p), v = varn(p);
218
if (w == v) pari_err_PRIORITY("charpoly", p, "=", w);
219
if (varncmp(w,v) < 0) p = gerepileupto(av, poleval(p, pol_x(v)));
220
return p;
221
}
222
GEN
223
caract(GEN x, long v)
224
{
225
pari_sp av = avma;
226
GEN T, C, x_k, Q;
227
long k, n;
228
229
if ((T = easychar(x,v))) return T;
230
231
n = lg(x)-1;
232
if (n == 1) return fix_pol(av, deg1pol(gen_1, gneg(gcoeff(x,1,1)), v));
233
234
x_k = pol_x(v); /* to be modified in place */
235
T = scalarpol(det(x), v); C = utoineg(n); Q = pol_x(v);
236
for (k=1; k<=n; k++)
237
{
238
GEN mk = utoineg(k), d;
239
gel(x_k,2) = mk;
240
d = det(RgM_Rg_add_shallow(x, mk));
241
T = RgX_add(RgX_mul(T, x_k), RgX_Rg_mul(Q, gmul(C, d)));
242
if (k == n) break;
243
244
Q = RgX_mul(Q, x_k);
245
C = diviuexact(mulsi(k-n,C), k+1); /* (-1)^k binomial(n,k) */
246
}
247
return fix_pol(av, RgX_Rg_div(T, mpfact(n)));
248
}
249
250
/* C = charpoly(x, v) */
251
static GEN
252
RgM_adj_from_char(GEN x, long v, GEN C)
253
{
254
if (varn(C) != v) /* problem with variable priorities */
255
{
256
C = gdiv(gsub(C, gsubst(C, v, gen_0)), pol_x(v));
257
if (odd(lg(x))) C = RgX_neg(C); /* even dimension */
258
return gsubst(C, v, x);
259
}
260
else
261
{
262
C = RgX_shift_shallow(C, -1);
263
if (odd(lg(x))) C = RgX_neg(C); /* even dimension */
264
return RgX_RgM_eval(C, x);
265
}
266
}
267
/* assume x square matrice */
268
static GEN
269
mattrace(GEN x)
270
{
271
long i, lx = lg(x);
272
GEN t;
273
if (lx < 3) return lx == 1? gen_0: gcopy(gcoeff(x,1,1));
274
t = gcoeff(x,1,1);
275
for (i = 2; i < lx; i++) t = gadd(t, gcoeff(x,i,i));
276
return t;
277
}
278
static int
279
bad_char(GEN q, long n)
280
{
281
forprime_t S;
282
ulong p;
283
if (!signe(q)) return 0;
284
(void)u_forprime_init(&S, 2, n);
285
while ((p = u_forprime_next(&S)))
286
if (!umodiu(q, p)) return 1;
287
return 0;
288
}
289
/* Using traces: return the characteristic polynomial of x (in variable v).
290
* If py != NULL, the adjoint matrix is put there. */
291
GEN
292
caradj(GEN x, long v, GEN *py)
293
{
294
pari_sp av, av0;
295
long i, k, n;
296
GEN T, y, t;
297
298
if ((T = easychar(x, v)))
299
{
300
if (py)
301
{
302
if (typ(x) != t_MAT) pari_err_TYPE("matadjoint",x);
303
*py = cgetg(1,t_MAT);
304
}
305
return T;
306
}
307
308
n = lg(x)-1; av0 = avma;
309
T = cgetg(n+3,t_POL); T[1] = evalsigne(1) | evalvarn(v);
310
gel(T,n+2) = gen_1;
311
if (!n) { if (py) *py = cgetg(1,t_MAT); return T; }
312
av = avma; t = gerepileupto(av, gneg(mattrace(x)));
313
gel(T,n+1) = t;
314
if (n == 1) {
315
T = fix_pol(av0, T);
316
if (py) *py = matid(1); return T;
317
}
318
if (n == 2) {
319
GEN a = gcoeff(x,1,1), b = gcoeff(x,1,2);
320
GEN c = gcoeff(x,2,1), d = gcoeff(x,2,2);
321
av = avma;
322
gel(T,2) = gerepileupto(av, gsub(gmul(a,d), gmul(b,c)));
323
T = fix_pol(av0, T);
324
if (py) {
325
y = cgetg(3, t_MAT);
326
gel(y,1) = mkcol2(gcopy(d), gneg(c));
327
gel(y,2) = mkcol2(gneg(b), gcopy(a));
328
*py = y;
329
}
330
return T;
331
}
332
/* l > 3 */
333
if (bad_char(residual_characteristic(x), n))
334
{ /* n! not invertible in base ring */
335
T = charpoly(x, v);
336
if (!py) return gerepileupto(av, T);
337
*py = RgM_adj_from_char(x, v, T);
338
gerepileall(av, 2, &T,py);
339
return T;
340
}
341
av = avma; y = RgM_shallowcopy(x);
342
for (i = 1; i <= n; i++) gcoeff(y,i,i) = gadd(gcoeff(y,i,i), t);
343
for (k = 2; k < n; k++)
344
{
345
GEN y0 = y;
346
y = RgM_mul(y, x);
347
t = gdivgs(mattrace(y), -k);
348
for (i = 1; i <= n; i++) gcoeff(y,i,i) = gadd(gcoeff(y,i,i), t);
349
y = gclone(y);
350
gel(T,n-k+2) = gerepilecopy(av, t); av = avma;
351
if (k > 2) gunclone(y0);
352
}
353
t = gmul(gcoeff(x,1,1),gcoeff(y,1,1));
354
for (i=2; i<=n; i++) t = gadd(t, gmul(gcoeff(x,1,i),gcoeff(y,i,1)));
355
gel(T,2) = gerepileupto(av, gneg(t));
356
T = fix_pol(av0, T);
357
if (py) *py = odd(n)? gcopy(y): RgM_neg(y);
358
gunclone(y); return T;
359
}
360
361
GEN
362
adj(GEN x)
363
{
364
GEN y;
365
(void)caradj(x, fetch_var(), &y);
366
(void)delete_var(); return y;
367
}
368
369
GEN
370
adjsafe(GEN x)
371
{
372
const long v = fetch_var();
373
pari_sp av = avma;
374
GEN C, A;
375
if (typ(x) != t_MAT) pari_err_TYPE("matadjoint",x);
376
if (lg(x) < 3) return gcopy(x);
377
C = charpoly(x,v);
378
A = RgM_adj_from_char(x, v, C);
379
(void)delete_var(); return gerepileupto(av, A);
380
}
381
382
GEN
383
matadjoint0(GEN x, long flag)
384
{
385
switch(flag)
386
{
387
case 0: return adj(x);
388
case 1: return adjsafe(x);
389
}
390
pari_err_FLAG("matadjoint");
391
return NULL; /* LCOV_EXCL_LINE */
392
}
393
394
/*******************************************************************/
395
/* */
396
/* Frobenius form */
397
/* */
398
/*******************************************************************/
399
400
/* The following section implement a mix of Ozello and Storjohann algorithms
401
402
P. Ozello, doctoral thesis (in French):
403
Calcul exact des formes de Jordan et de Frobenius d'une matrice, Chapitre 2
404
http://tel.archives-ouvertes.fr/tel-00323705
405
406
A. Storjohann, Diss. ETH No. 13922
407
Algorithms for Matrix Canonical Forms, Chapter 9
408
https://cs.uwaterloo.ca/~astorjoh/diss2up.pdf
409
410
We use Storjohann Lemma 9.14 (step1, step2, step3) Ozello theorem 4,
411
and Storjohann Lemma 9.18
412
*/
413
414
/* Elementary transforms */
415
416
/* M <- U^(-1) M U, U = E_{i,j}(k) => U^(-1) = E{i,j}(-k)
417
* P = U * P */
418
static void
419
transL(GEN M, GEN P, GEN k, long i, long j)
420
{
421
long l, n = lg(M)-1;
422
for(l=1; l<=n; l++) /* M[,j]-=k*M[,i] */
423
gcoeff(M,l,j) = gsub(gcoeff(M,l,j), gmul(gcoeff(M,l,i), k));
424
for(l=1; l<=n; l++) /* M[i,]+=k*M[j,] */
425
gcoeff(M,i,l) = gadd(gcoeff(M,i,l), gmul(gcoeff(M,j,l), k));
426
if (P)
427
for(l=1; l<=n; l++)
428
gcoeff(P,i,l) = gadd(gcoeff(P,i,l), gmul(gcoeff(P,j,l), k));
429
}
430
431
/* j = a or b */
432
static void
433
transD(GEN M, GEN P, long a, long b, long j)
434
{
435
long l, n;
436
GEN k = gcoeff(M,a,b), ki;
437
438
if (gequal1(k)) return;
439
ki = ginv(k); n = lg(M)-1;
440
for(l=1; l<=n; l++)
441
if (l!=j)
442
{
443
gcoeff(M,l,j) = gmul(gcoeff(M,l,j), k);
444
gcoeff(M,j,l) = (j==a && l==b)? gen_1: gmul(gcoeff(M,j,l), ki);
445
}
446
if (P)
447
for(l=1; l<=n; l++)
448
gcoeff(P,j,l) = gmul(gcoeff(P,j,l), ki);
449
}
450
451
static void
452
transS(GEN M, GEN P, long i, long j)
453
{
454
long l, n = lg(M)-1;
455
swap(gel(M,i), gel(M,j));
456
for (l=1; l<=n; l++)
457
swap(gcoeff(M,i,l), gcoeff(M,j,l));
458
if (P)
459
for (l=1; l<=n; l++)
460
swap(gcoeff(P,i,l), gcoeff(P,j,l));
461
}
462
463
/* Convert companion matrix to polynomial*/
464
static GEN
465
minpoly_polslice(GEN M, long i, long j, long v)
466
{
467
long k, d = j+1-i;
468
GEN P = cgetg(d+3,t_POL);
469
P[1] = evalsigne(1)|evalvarn(v);
470
for (k=0; k<d; k++)
471
gel(P,k+2) = gneg(gcoeff(M,i+k, j));
472
gel(P,d+2) = gen_1;
473
return P;
474
}
475
476
static GEN
477
minpoly_listpolslice(GEN M, GEN V, long v)
478
{
479
long i, n = lg(M)-1, nb = lg(V)-1;
480
GEN W = cgetg(nb+1, t_VEC);
481
for (i=1; i<=nb; i++)
482
gel(W,i) = minpoly_polslice(M, V[i], i < nb? V[i+1]-1: n, v);
483
return W;
484
}
485
486
static int
487
minpoly_dvdslice(GEN M, long i, long j, long k)
488
{
489
pari_sp av = avma;
490
long r = signe(RgX_rem(minpoly_polslice(M, i, j-1, 0),
491
minpoly_polslice(M, j, k, 0)));
492
return gc_bool(av, r == 0);
493
}
494
495
static void
496
RgM_replace(GEN M, GEN M2)
497
{
498
long n = lg(M)-1, m = nbrows(M), i, j;
499
for(i=1; i<=n; i++)
500
for(j=1; j<=m; j++)
501
gcoeff(M, i, j) = gcoeff(M2, i, j);
502
}
503
504
static void
505
gerepilemat2_inplace(pari_sp av, GEN M, GEN P)
506
{
507
GEN M2 = M, P2 = P;
508
gerepileall(av, P ? 2: 1, &M2, &P2);
509
RgM_replace(M, M2);
510
if (P) RgM_replace(P, P2);
511
}
512
513
/* Lemma 9.14 */
514
static long
515
weakfrobenius_step1(GEN M, GEN P, long j0)
516
{
517
pari_sp av = avma;
518
long n = lg(M)-1, k, j;
519
for (j = j0; j < n; ++j)
520
{
521
if (gequal0(gcoeff(M, j+1, j)))
522
{
523
for (k = j+2; k <= n; ++k)
524
if (!gequal0(gcoeff(M,k,j))) break;
525
if (k > n) return j;
526
transS(M, P, k, j+1);
527
}
528
transD(M, P, j+1, j, j+1);
529
/* Now M[j+1,j] = 1 */
530
for (k = 1; k <= n; ++k)
531
if (k != j+1 && !gequal0(gcoeff(M,k,j))) /* zero M[k,j] */
532
{
533
transL(M, P, gneg(gcoeff(M,k,j)), k, j+1);
534
gcoeff(M,k,j) = gen_0; /* avoid approximate 0 */
535
}
536
if (gc_needed(av,1))
537
{
538
if (DEBUGMEM > 1)
539
pari_warn(warnmem,"RgM_minpoly stage 1: j0=%ld, j=%ld", j0, j);
540
gerepilemat2_inplace(av, M, P);
541
}
542
}
543
return n;
544
}
545
546
static void
547
weakfrobenius_step2(GEN M, GEN P, long j)
548
{
549
pari_sp av = avma;
550
long i, k, n = lg(M)-1;
551
for(i=j; i>=2; i--)
552
{
553
for(k=j+1; k<=n; k++)
554
if (!gequal0(gcoeff(M,i,k)))
555
transL(M, P, gcoeff(M,i,k), i-1, k);
556
if (gc_needed(av,1))
557
{
558
if (DEBUGMEM > 1)
559
pari_warn(warnmem,"RgM_minpoly stage 2: j=%ld, i=%ld", j, i);
560
gerepilemat2_inplace(av, M, P);
561
}
562
}
563
}
564
565
static long
566
weakfrobenius_step3(GEN M, GEN P, long j0, long j)
567
{
568
long i, k, n = lg(M)-1;
569
if (j == n) return 0;
570
if (gequal0(gcoeff(M, j0, j+1)))
571
{
572
for (k=j+2; k<=n; k++)
573
if (!gequal0(gcoeff(M, j0, k))) break;
574
if (k > n) return 0;
575
transS(M, P, k, j+1);
576
}
577
transD(M, P, j0, j+1, j+1);
578
for (i=j+2; i<=n; i++)
579
if (!gequal0(gcoeff(M, j0, i)))
580
transL(M, P, gcoeff(M, j0, i),j+1, i);
581
return 1;
582
}
583
584
/* flag: 0 -> full Frobenius from , 1 -> weak Frobenius form */
585
GEN
586
RgM_Frobenius(GEN M, long flag, GEN *pt_P, GEN *pt_v)
587
{
588
pari_sp av = avma, av2, ltop;
589
long n = lg(M)-1, eps, j0 = 1, nb = 0;
590
GEN v, P;
591
v = cgetg(n+1, t_VECSMALL);
592
ltop = avma;
593
P = pt_P ? matid(n): NULL;
594
M = RgM_shallowcopy(M);
595
av2 = avma;
596
while (j0 <= n)
597
{
598
long j = weakfrobenius_step1(M, P, j0);
599
weakfrobenius_step2(M, P, j);
600
eps = weakfrobenius_step3(M, P, j0, j);
601
if (eps == 0)
602
{
603
v[++nb] = j0;
604
if (flag == 0 && nb > 1 && !minpoly_dvdslice(M, v[nb-1], j0, j))
605
{
606
j = j0; j0 = v[nb-1]; nb -= 2;
607
transL(M, P, gen_1, j, j0); /*lemma 9.18*/
608
} else
609
j0 = j+1;
610
}
611
else
612
transS(M, P, j0, j+1); /*theorem 4*/
613
if (gc_needed(av,1))
614
{
615
if (DEBUGMEM > 1)
616
pari_warn(warnmem,"weakfrobenius j0=%ld",j0);
617
gerepilemat2_inplace(av2, M, P);
618
}
619
}
620
fixlg(v, nb+1);
621
if (pt_v) *pt_v = v;
622
gerepileall(pt_v ? ltop: av, P? 2: 1, &M, &P);
623
if (pt_P) *pt_P = P;
624
return M;
625
}
626
627
static GEN
628
RgM_minpoly(GEN M, long v)
629
{
630
pari_sp av = avma;
631
GEN V, W;
632
if (lg(M) == 1) return pol_1(v);
633
M = RgM_Frobenius(M, 1, NULL, &V);
634
W = minpoly_listpolslice(M, V, v);
635
if (varncmp(v,gvar2(W)) >= 0)
636
pari_err_PRIORITY("matfrobenius", M, "<=", v);
637
return gerepileupto(av, RgX_normalize(glcm0(W, NULL)));
638
}
639
640
GEN
641
Frobeniusform(GEN V, long n)
642
{
643
long i, j, k;
644
GEN M = zeromatcopy(n,n);
645
for (k=1,i=1;i<lg(V);i++,k++)
646
{
647
GEN P = gel(V,i);
648
long d = degpol(P);
649
if (k+d-1 > n) pari_err_PREC("matfrobenius");
650
for (j=0; j<d-1; j++, k++) gcoeff(M,k+1,k) = gen_1;
651
for (j=0; j<d; j++) gcoeff(M,k-j,k) = gneg(gel(P, 1+d-j));
652
}
653
return M;
654
}
655
656
GEN
657
matfrobenius(GEN M, long flag, long v)
658
{
659
long n;
660
if (typ(M)!=t_MAT) pari_err_TYPE("matfrobenius",M);
661
if (v < 0) v = 0;
662
n = lg(M)-1;
663
if (n && lgcols(M)!=n+1) pari_err_DIM("matfrobenius");
664
if (flag > 2) pari_err_FLAG("matfrobenius");
665
switch (flag)
666
{
667
case 0:
668
return RgM_Frobenius(M, 0, NULL, NULL);
669
case 1:
670
{
671
pari_sp av = avma;
672
GEN V, W, F;
673
F = RgM_Frobenius(M, 0, NULL, &V);
674
W = minpoly_listpolslice(F, V, v);
675
if (varncmp(v, gvar2(W)) >= 0)
676
pari_err_PRIORITY("matfrobenius", M, "<=", v);
677
return gerepileupto(av, W);
678
}
679
case 2:
680
{
681
GEN P, F, R = cgetg(3, t_VEC);
682
F = RgM_Frobenius(M, 0, &P, NULL);
683
gel(R,1) = F; gel(R,2) = P;
684
return R;
685
}
686
default:
687
pari_err_FLAG("matfrobenius");
688
}
689
return NULL; /*LCOV_EXCL_LINE*/
690
}
691
692
/*******************************************************************/
693
/* */
694
/* MINIMAL POLYNOMIAL */
695
/* */
696
/*******************************************************************/
697
698
static GEN
699
RgXQ_minpoly_naive(GEN y, GEN P)
700
{
701
long n = lgpol(P);
702
GEN M = ker(RgXQ_matrix_pow(y,n,n,P));
703
return content(RgM_to_RgXV(M,varn(P)));
704
}
705
706
static GEN
707
RgXQ_minpoly_FpXQ(GEN x, GEN y, GEN p)
708
{
709
pari_sp av = avma;
710
GEN r;
711
if (lgefint(p) == 3)
712
{
713
ulong pp = uel(p, 2);
714
r = Flx_to_ZX_inplace(Flxq_minpoly(RgX_to_Flx(x, pp),
715
RgX_to_Flx(y, pp), pp));
716
}
717
else
718
r = FpXQ_minpoly(RgX_to_FpX(x, p), RgX_to_FpX(y, p), p);
719
return gerepileupto(av, FpX_to_mod(r, p));
720
}
721
722
static GEN
723
RgXQ_minpoly_FpXQXQ(GEN x, GEN S, GEN pol, GEN p)
724
{
725
pari_sp av = avma;
726
GEN r;
727
GEN T = RgX_to_FpX(pol, p);
728
if (signe(T)==0) pari_err_OP("minpoly",x,S);
729
if (lgefint(p) == 3)
730
{
731
ulong pp = uel(p, 2);
732
GEN Tp = ZX_to_Flx(T, pp);
733
r = FlxX_to_ZXX(FlxqXQ_minpoly(RgX_to_FlxqX(x, Tp, pp),
734
RgX_to_FlxqX(S, Tp, pp), Tp, pp));
735
}
736
else
737
r = FpXQXQ_minpoly(RgX_to_FpXQX(x, T, p), RgX_to_FpXQX(S, T, p), T, p);
738
return gerepileupto(av, FpXQX_to_mod(r, T, p));
739
}
740
741
#define code(t1,t2) ((t1 << 6) | t2)
742
743
static GEN
744
RgXQ_minpoly_fast(GEN x, GEN y)
745
{
746
GEN p, pol;
747
long pa;
748
long t = RgX_type2(x,y, &p,&pol,&pa);
749
switch(t)
750
{
751
case t_INTMOD: return RgXQ_minpoly_FpXQ(x, y, p);
752
case t_FFELT: return FFXQ_minpoly(x, y, pol);
753
case code(t_POLMOD, t_INTMOD):
754
return RgXQ_minpoly_FpXQXQ(x, y, pol, p);
755
default: return NULL;
756
}
757
}
758
759
#undef code
760
761
static GEN
762
easymin(GEN x, long v)
763
{
764
GEN G, R, dR;
765
long tx = typ(x);
766
if (tx == t_FFELT)
767
{
768
R = FpX_to_mod(FF_minpoly(x), FF_p_i(x));
769
setvarn(R,v); return R;
770
}
771
if (tx == t_POLMOD)
772
{
773
GEN a = gel(x,2), b = gel(x,1);
774
if (typ(a) != t_POL || varn(a) != varn(b))
775
{
776
if (varncmp(gvar(a), v) <= 0) pari_err_PRIORITY("minpoly", x, "<", v);
777
return deg1pol(gen_1, gneg_i(a), v);
778
}
779
R = RgXQ_minpoly_fast(a,b);
780
if (R) { setvarn(R, v); return R; }
781
if (!issquarefree(b))
782
{
783
R = RgXQ_minpoly_naive(a, b);
784
setvarn(R,v); return R;
785
}
786
}
787
R = easychar(x, v); if (!R) return NULL;
788
dR = RgX_deriv(R); if (!lgpol(dR)) return NULL;
789
G = RgX_normalize(RgX_gcd(R,dR));
790
return RgX_div(R,G);
791
}
792
GEN
793
minpoly(GEN x, long v)
794
{
795
pari_sp av = avma;
796
GEN P;
797
if (v < 0) v = 0;
798
P = easymin(x,v);
799
if (P) return gerepileupto(av,P);
800
/* typ(x) = t_MAT */
801
set_avma(av); return RgM_minpoly(x,v);
802
}
803
804
/*******************************************************************/
805
/* */
806
/* HESSENBERG FORM */
807
/* */
808
/*******************************************************************/
809
static int
810
relative0(GEN a, GEN a0, long bit)
811
{ return gequal0(a) || (bit && gexpo(a)-gexpo(a0) < bit); }
812
/* assume x a nonempty square t_MAT */
813
static GEN
814
RgM_hess(GEN x0, long prec)
815
{
816
pari_sp av;
817
long lx = lg(x0), bit = prec? 8 - bit_accuracy(prec): 0, m, i, j;
818
GEN x;
819
820
if (bit) x0 = RgM_shallowcopy(x0);
821
av = avma; x = RgM_shallowcopy(x0);
822
for (m=2; m<lx-1; m++)
823
{
824
GEN t = NULL;
825
for (i=m+1; i<lx; i++)
826
{
827
t = gcoeff(x,i,m-1);
828
if (!relative0(t, gcoeff(x0,i,m-1), bit)) break;
829
}
830
if (i == lx) continue;
831
for (j=m-1; j<lx; j++) swap(gcoeff(x,i,j), gcoeff(x,m,j));
832
swap(gel(x,i), gel(x,m));
833
if (bit)
834
{
835
for (j=m-1; j<lx; j++) swap(gcoeff(x0,i,j), gcoeff(x0,m,j));
836
swap(gel(x0,i), gel(x0,m));
837
}
838
t = ginv(t);
839
840
for (i=m+1; i<lx; i++)
841
{
842
GEN c = gcoeff(x,i,m-1);
843
if (gequal0(c)) continue;
844
845
c = gmul(c,t); gcoeff(x,i,m-1) = gen_0;
846
for (j=m; j<lx; j++)
847
gcoeff(x,i,j) = gsub(gcoeff(x,i,j), gmul(c,gcoeff(x,m,j)));
848
for (j=1; j<lx; j++)
849
gcoeff(x,j,m) = gadd(gcoeff(x,j,m), gmul(c,gcoeff(x,j,i)));
850
if (gc_needed(av,2))
851
{
852
if (DEBUGMEM>1) pari_warn(warnmem,"hess, m = %ld", m);
853
gerepileall(av,2, &x, &t);
854
}
855
}
856
}
857
return x;
858
}
859
860
GEN
861
hess(GEN x)
862
{
863
pari_sp av = avma;
864
GEN p = NULL, T = NULL;
865
long prec, lx = lg(x);
866
if (typ(x) != t_MAT) pari_err_TYPE("hess",x);
867
if (lx == 1) return cgetg(1,t_MAT);
868
if (lgcols(x) != lx) pari_err_DIM("hess");
869
switch(RgM_type(x, &p, &T, &prec))
870
{
871
case t_REAL:
872
case t_COMPLEX: break;
873
default: prec = 0;
874
}
875
return gerepilecopy(av, RgM_hess(x,prec));
876
}
877
878
GEN
879
Flm_hess(GEN x, ulong p)
880
{
881
long lx = lg(x), m, i, j;
882
if (lx == 1) return cgetg(1,t_MAT);
883
if (lgcols(x) != lx) pari_err_DIM("hess");
884
885
x = Flm_copy(x);
886
for (m=2; m<lx-1; m++)
887
{
888
ulong t = 0;
889
for (i=m+1; i<lx; i++) { t = ucoeff(x,i,m-1); if (t) break; }
890
if (i == lx) continue;
891
for (j=m-1; j<lx; j++) lswap(ucoeff(x,i,j), ucoeff(x,m,j));
892
swap(gel(x,i), gel(x,m)); t = Fl_inv(t, p);
893
894
for (i=m+1; i<lx; i++)
895
{
896
ulong c = ucoeff(x,i,m-1);
897
if (!c) continue;
898
899
c = Fl_mul(c,t,p); ucoeff(x,i,m-1) = 0;
900
for (j=m; j<lx; j++)
901
ucoeff(x,i,j) = Fl_sub(ucoeff(x,i,j), Fl_mul(c,ucoeff(x,m,j), p), p);
902
for (j=1; j<lx; j++)
903
ucoeff(x,j,m) = Fl_add(ucoeff(x,j,m), Fl_mul(c,ucoeff(x,j,i), p), p);
904
}
905
}
906
return x;
907
}
908
GEN
909
FpM_hess(GEN x, GEN p)
910
{
911
pari_sp av = avma;
912
long lx = lg(x), m, i, j;
913
if (lx == 1) return cgetg(1,t_MAT);
914
if (lgcols(x) != lx) pari_err_DIM("hess");
915
if (lgefint(p) == 3)
916
{
917
ulong pp = p[2];
918
x = Flm_hess(ZM_to_Flm(x, pp), pp);
919
return gerepileupto(av, Flm_to_ZM(x));
920
}
921
x = RgM_shallowcopy(x);
922
for (m=2; m<lx-1; m++)
923
{
924
GEN t = NULL;
925
for (i=m+1; i<lx; i++) { t = gcoeff(x,i,m-1); if (signe(t)) break; }
926
if (i == lx) continue;
927
for (j=m-1; j<lx; j++) swap(gcoeff(x,i,j), gcoeff(x,m,j));
928
swap(gel(x,i), gel(x,m)); t = Fp_inv(t, p);
929
930
for (i=m+1; i<lx; i++)
931
{
932
GEN c = gcoeff(x,i,m-1);
933
if (!signe(c)) continue;
934
935
c = Fp_mul(c,t, p); gcoeff(x,i,m-1) = gen_0;
936
for (j=m; j<lx; j++)
937
gcoeff(x,i,j) = Fp_sub(gcoeff(x,i,j), Fp_mul(c,gcoeff(x,m,j),p), p);
938
for (j=1; j<lx; j++)
939
gcoeff(x,j,m) = Fp_add(gcoeff(x,j,m), Fp_mul(c,gcoeff(x,j,i),p), p);
940
if (gc_needed(av,2))
941
{
942
if (DEBUGMEM>1) pari_warn(warnmem,"hess, m = %ld", m);
943
gerepileall(av,2, &x, &t);
944
}
945
}
946
}
947
return gerepilecopy(av,x);
948
}
949
/* H in Hessenberg form */
950
static GEN
951
RgM_hess_charpoly(GEN H, long v)
952
{
953
long n = lg(H), r, i;
954
GEN y = cgetg(n+1, t_VEC);
955
gel(y,1) = pol_1(v);
956
for (r = 1; r < n; r++)
957
{
958
pari_sp av2 = avma;
959
GEN z, a = gen_1, b = pol_0(v);
960
for (i = r-1; i; i--)
961
{
962
a = gmul(a, gcoeff(H,i+1,i));
963
if (gequal0(a)) break;
964
b = RgX_add(b, RgX_Rg_mul(gel(y,i), gmul(a,gcoeff(H,i,r))));
965
}
966
z = RgX_sub(RgX_shift_shallow(gel(y,r), 1),
967
RgX_Rg_mul(gel(y,r), gcoeff(H,r,r)));
968
gel(y,r+1) = gerepileupto(av2, RgX_sub(z, b)); /* (X - H[r,r])y[r] - b */
969
}
970
return gel(y,n);
971
}
972
GEN
973
carhess(GEN x, long v)
974
{
975
pari_sp av;
976
GEN y;
977
if ((y = easychar(x,v))) return y;
978
av = avma; y = RgM_hess_charpoly(hess(x), v);
979
return fix_pol(av, y);
980
}
981
982
/* Bound for sup norm of charpoly(M/dM), M integral: let B = |M|oo / |dM|,
983
* s = max_k binomial(n,k) (kB^2)^(k/2),
984
* return ceil(log2(s)) */
985
static long
986
charpoly_bound(GEN M, GEN dM, GEN N)
987
{
988
pari_sp av = avma;
989
GEN B = itor(N, LOWDEFAULTPREC);
990
GEN s = real_0(LOWDEFAULTPREC), bin, B2;
991
long n = lg(M)-1, k;
992
bin = gen_1;
993
if (dM) B = divri(B, dM);
994
B2 = sqrr(B);
995
for (k = n; k >= (n+1)>>1; k--)
996
{
997
GEN t = mulri(powruhalf(mulur(k, B2), k), bin);
998
if (abscmprr(t, s) > 0) s = t;
999
bin = diviuexact(muliu(bin, k), n-k+1);
1000
}
1001
return gc_long(av, ceil(dbllog2(s)));
1002
}
1003
1004
static GEN
1005
QM_charpoly_ZX_slice(GEN A, GEN dM, GEN P, GEN *mod)
1006
{
1007
pari_sp av = avma;
1008
long i, n = lg(P)-1;
1009
GEN H, T;
1010
if (n == 1)
1011
{
1012
ulong p = uel(P,1), dp = dM ? umodiu(dM, p): 1;
1013
GEN Hp, a = ZM_to_Flm(A, p);
1014
Hp = Flm_charpoly_i(a, p);
1015
if (dp != 1) Hp = Flx_rescale(Hp, Fl_inv(dp, p), p);
1016
Hp = gerepileupto(av, Flx_to_ZX(Hp));
1017
*mod = utoipos(p); return Hp;
1018
}
1019
T = ZV_producttree(P);
1020
A = ZM_nv_mod_tree(A, P, T);
1021
H = cgetg(n+1, t_VEC);
1022
for(i=1; i <= n; i++)
1023
{
1024
ulong p = uel(P,i), dp = dM ? umodiu(dM, p): 1;
1025
gel(H,i) = Flm_charpoly(gel(A, i), p);
1026
if (dp != 1) gel(H,i) = Flx_rescale(gel(H,i), Fl_inv(dp, p), p);
1027
}
1028
H = nxV_chinese_center_tree(H, P, T, ZV_chinesetree(P,T));
1029
*mod = gmael(T, lg(T)-1, 1);
1030
gerepileall(av, 2, &H, mod); return H;
1031
}
1032
1033
GEN
1034
QM_charpoly_ZX_worker(GEN P, GEN M, GEN dM)
1035
{
1036
GEN V = cgetg(3, t_VEC);
1037
gel(V,1) = QM_charpoly_ZX_slice(M, equali1(dM) ? NULL:dM, P, &gel(V,2));
1038
return V;
1039
}
1040
1041
/* Assume M a square ZM, dM integer. Return charpoly(M / dM) in Z[X] */
1042
static GEN
1043
QM_charpoly_ZX_i(GEN M, GEN dM, long bound)
1044
{
1045
long n = lg(M)-1;
1046
forprime_t S;
1047
GEN worker = snm_closure(is_entry("_QM_charpoly_ZX_worker"),
1048
mkvec2(M, dM? dM: gen_1));
1049
if (!n) return pol_1(0);
1050
if (bound < 0)
1051
{
1052
GEN N = ZM_supnorm(M);
1053
if (signe(N) == 0) return monomial(gen_1, n, 0);
1054
bound = charpoly_bound(M, dM, N) + 1;
1055
}
1056
if (DEBUGLEVEL>5) err_printf("ZM_charpoly: bound 2^%ld\n", bound);
1057
init_modular_big(&S);
1058
return gen_crt("QM_charpoly_ZX", worker, &S, dM, bound, 0, NULL,
1059
nxV_chinese_center, FpX_center);
1060
}
1061
1062
GEN
1063
QM_charpoly_ZX_bound(GEN M, long bit)
1064
{
1065
pari_sp av = avma;
1066
GEN dM; M = Q_remove_denom(M, &dM);
1067
return gerepilecopy(av, QM_charpoly_ZX_i(M, dM, bit));
1068
}
1069
GEN
1070
QM_charpoly_ZX(GEN M)
1071
{
1072
pari_sp av = avma;
1073
GEN dM; M = Q_remove_denom(M, &dM);
1074
return gerepilecopy(av, QM_charpoly_ZX_i(M, dM, -1));
1075
}
1076
GEN
1077
ZM_charpoly(GEN M)
1078
{
1079
pari_sp av = avma;
1080
return gerepilecopy(av, QM_charpoly_ZX_i(M, NULL, -1));
1081
}
1082
1083
/*******************************************************************/
1084
/* */
1085
/* CHARACTERISTIC POLYNOMIAL (BERKOWITZ'S ALGORITHM) */
1086
/* */
1087
/*******************************************************************/
1088
GEN
1089
carberkowitz(GEN x, long v)
1090
{
1091
long lx, i, j, k, r;
1092
GEN V, S, C, Q;
1093
pari_sp av0, av;
1094
if ((V = easychar(x,v))) return V;
1095
lx = lg(x); av0 = avma;
1096
V = cgetg(lx+1, t_VEC);
1097
S = cgetg(lx+1, t_VEC);
1098
C = cgetg(lx+1, t_VEC);
1099
Q = cgetg(lx+1, t_VEC);
1100
av = avma;
1101
gel(C,1) = gen_m1;
1102
gel(V,1) = gen_m1;
1103
for (i=2;i<=lx; i++) gel(C,i) = gel(Q,i) = gel(S,i) = gel(V,i) = gen_0;
1104
gel(V,2) = gcoeff(x,1,1);
1105
for (r = 2; r < lx; r++)
1106
{
1107
pari_sp av2;
1108
GEN t;
1109
1110
for (i = 1; i < r; i++) gel(S,i) = gcoeff(x,i,r);
1111
gel(C,2) = gcoeff(x,r,r);
1112
for (i = 1; i < r-1; i++)
1113
{
1114
av2 = avma; t = gmul(gcoeff(x,r,1), gel(S,1));
1115
for (j = 2; j < r; j++) t = gadd(t, gmul(gcoeff(x,r,j), gel(S,j)));
1116
gel(C,i+2) = gerepileupto(av2, t);
1117
for (j = 1; j < r; j++)
1118
{
1119
av2 = avma; t = gmul(gcoeff(x,j,1), gel(S,1));
1120
for (k = 2; k < r; k++) t = gadd(t, gmul(gcoeff(x,j,k), gel(S,k)));
1121
gel(Q,j) = gerepileupto(av2, t);
1122
}
1123
for (j = 1; j < r; j++) gel(S,j) = gel(Q,j);
1124
}
1125
av2 = avma; t = gmul(gcoeff(x,r,1), gel(S,1));
1126
for (j = 2; j < r; j++) t = gadd(t, gmul(gcoeff(x,r,j), gel(S,j)));
1127
gel(C,r+1) = gerepileupto(av2, t);
1128
if (gc_needed(av0,1))
1129
{
1130
if (DEBUGMEM>1) pari_warn(warnmem,"carberkowitz");
1131
gerepileall(av, 2, &C, &V);
1132
}
1133
for (i = 1; i <= r+1; i++)
1134
{
1135
av2 = avma; t = gmul(gel(C,i), gel(V,1));
1136
for (j = 2; j <= minss(r,i); j++)
1137
t = gadd(t, gmul(gel(C,i+1-j), gel(V,j)));
1138
gel(Q,i) = gerepileupto(av2, t);
1139
}
1140
for (i = 1; i <= r+1; i++) gel(V,i) = gel(Q,i);
1141
}
1142
V = RgV_to_RgX(vecreverse(V), v); /* not gtopoly: fail if v > gvar(V) */
1143
V = odd(lx)? gcopy(V): RgX_neg(V);
1144
return fix_pol(av0, V);
1145
}
1146
1147
/*******************************************************************/
1148
/* */
1149
/* NORMS */
1150
/* */
1151
/*******************************************************************/
1152
GEN
1153
gnorm(GEN x)
1154
{
1155
pari_sp av;
1156
long lx, i;
1157
GEN y;
1158
1159
switch(typ(x))
1160
{
1161
case t_INT: return sqri(x);
1162
case t_REAL: return sqrr(x);
1163
case t_FRAC: return sqrfrac(x);
1164
case t_COMPLEX: av = avma; return gerepileupto(av, cxnorm(x));
1165
case t_QUAD: av = avma; return gerepileupto(av, quadnorm(x));
1166
1167
case t_POL: case t_SER: case t_RFRAC: av = avma;
1168
return gerepileupto(av, greal(gmul(conj_i(x),x)));
1169
1170
case t_FFELT:
1171
y = cgetg(3, t_INTMOD);
1172
gel(y,1) = FF_p(x);
1173
gel(y,2) = FF_norm(x); return y;
1174
1175
case t_POLMOD:
1176
{
1177
GEN T = gel(x,1), a = gel(x,2);
1178
if (typ(a) != t_POL || varn(a) != varn(T)) return gpowgs(a, degpol(T));
1179
return RgXQ_norm(a, T);
1180
}
1181
case t_VEC: case t_COL: case t_MAT:
1182
y = cgetg_copy(x, &lx);
1183
for (i=1; i<lx; i++) gel(y,i) = gnorm(gel(x,i));
1184
return y;
1185
}
1186
pari_err_TYPE("gnorm",x);
1187
return NULL; /* LCOV_EXCL_LINE */
1188
}
1189
1190
/* return |q|^2, complex modulus */
1191
static GEN
1192
cxquadnorm(GEN q, long prec)
1193
{
1194
GEN X = gel(q,1), c = gel(X,2); /* (1-D)/4, -D/4 */
1195
if (signe(c) > 0) return quadnorm(q); /* imaginary */
1196
if (!prec) pari_err_TYPE("gnorml2", q);
1197
return sqrr(quadtofp(q, prec));
1198
}
1199
1200
static GEN
1201
gnorml2_i(GEN x, long prec)
1202
{
1203
pari_sp av;
1204
long i, lx;
1205
GEN s;
1206
1207
switch(typ(x))
1208
{
1209
case t_INT: return sqri(x);
1210
case t_REAL: return sqrr(x);
1211
case t_FRAC: return sqrfrac(x);
1212
case t_COMPLEX: av = avma; return gerepileupto(av, cxnorm(x));
1213
case t_QUAD: av = avma; return gerepileupto(av, cxquadnorm(x,prec));
1214
1215
case t_POL: lx = lg(x)-1; x++; break;
1216
1217
case t_VEC:
1218
case t_COL:
1219
case t_MAT: lx = lg(x); break;
1220
1221
default: pari_err_TYPE("gnorml2",x);
1222
return NULL; /* LCOV_EXCL_LINE */
1223
}
1224
if (lx == 1) return gen_0;
1225
av = avma;
1226
s = gnorml2(gel(x,1));
1227
for (i=2; i<lx; i++)
1228
{
1229
s = gadd(s, gnorml2(gel(x,i)));
1230
if (gc_needed(av,1))
1231
{
1232
if (DEBUGMEM>1) pari_warn(warnmem,"gnorml2");
1233
s = gerepileupto(av, s);
1234
}
1235
}
1236
return gerepileupto(av,s);
1237
}
1238
GEN
1239
gnorml2(GEN x) { return gnorml2_i(x, 0); }
1240
1241
static GEN pnormlp(GEN,GEN,long);
1242
static GEN
1243
pnormlpvec(long i0, GEN x, GEN p, long prec)
1244
{
1245
pari_sp av = avma;
1246
long i, lx = lg(x);
1247
GEN s = gen_0;
1248
for (i=i0; i<lx; i++)
1249
{
1250
s = gadd(s, pnormlp(gel(x,i),p,prec));
1251
if (gc_needed(av,1))
1252
{
1253
if (DEBUGMEM>1) pari_warn(warnmem,"gnormlp, i = %ld", i);
1254
s = gerepileupto(av, s);
1255
}
1256
}
1257
return s;
1258
}
1259
/* (||x||_p)^p */
1260
static GEN
1261
pnormlp(GEN x, GEN p, long prec)
1262
{
1263
switch(typ(x))
1264
{
1265
case t_INT: case t_REAL: x = mpabs(x); break;
1266
case t_FRAC: x = absfrac(x); break;
1267
case t_COMPLEX: case t_QUAD: x = gabs(x,prec); break;
1268
case t_POL: return pnormlpvec(2, x, p, prec);
1269
case t_VEC: case t_COL: case t_MAT: return pnormlpvec(1, x, p, prec);
1270
default: pari_err_TYPE("gnormlp",x);
1271
}
1272
return gpow(x, p, prec);
1273
}
1274
1275
GEN
1276
gnormlp(GEN x, GEN p, long prec)
1277
{
1278
pari_sp av = avma;
1279
if (!p || (typ(p) == t_INFINITY && inf_get_sign(p) > 0))
1280
return gsupnorm(x, prec);
1281
if (gsigne(p) <= 0) pari_err_DOMAIN("normlp", "p", "<=", gen_0, p);
1282
if (is_scalar_t(typ(x))) return gabs(x, prec);
1283
if (typ(p) == t_INT)
1284
{
1285
ulong pp = itou_or_0(p);
1286
switch(pp)
1287
{
1288
case 1: return gnorml1(x, prec);
1289
case 2: x = gnorml2_i(x, prec); break;
1290
default: x = pnormlp(x, p, prec); break;
1291
}
1292
if (pp && typ(x) == t_INT && Z_ispowerall(x, pp, &x))
1293
return gerepileuptoleaf(av, x);
1294
if (pp == 2) return gerepileupto(av, gsqrt(x, prec));
1295
}
1296
else
1297
x = pnormlp(x, p, prec);
1298
x = gpow(x, ginv(p), prec);
1299
return gerepileupto(av, x);
1300
}
1301
1302
GEN
1303
gnorml1(GEN x,long prec)
1304
{
1305
pari_sp av = avma;
1306
long lx,i;
1307
GEN s;
1308
switch(typ(x))
1309
{
1310
case t_INT: case t_REAL: return mpabs(x);
1311
case t_FRAC: return absfrac(x);
1312
1313
case t_COMPLEX: case t_QUAD:
1314
return gabs(x,prec);
1315
1316
case t_POL:
1317
lx = lg(x); s = gen_0;
1318
for (i=2; i<lx; i++) s = gadd(s, gnorml1(gel(x,i),prec));
1319
break;
1320
1321
case t_VEC: case t_COL: case t_MAT:
1322
lx = lg(x); s = gen_0;
1323
for (i=1; i<lx; i++) s = gadd(s, gnorml1(gel(x,i),prec));
1324
break;
1325
1326
default: pari_err_TYPE("gnorml1",x);
1327
return NULL; /* LCOV_EXCL_LINE */
1328
}
1329
return gerepileupto(av, s);
1330
}
1331
/* As gnorml1, except for t_QUAD and t_COMPLEX: |x + wy| := |x| + |y|
1332
* Still a norm of R-vector spaces, and can be cheaply computed without
1333
* square roots */
1334
GEN
1335
gnorml1_fake(GEN x)
1336
{
1337
pari_sp av = avma;
1338
long lx, i;
1339
GEN s;
1340
switch(typ(x))
1341
{
1342
case t_INT: case t_REAL: return mpabs(x);
1343
case t_FRAC: return absfrac(x);
1344
1345
case t_COMPLEX:
1346
s = gadd(gnorml1_fake(gel(x,1)), gnorml1_fake(gel(x,2)));
1347
break;
1348
case t_QUAD:
1349
s = gadd(gnorml1_fake(gel(x,2)), gnorml1_fake(gel(x,3)));
1350
break;
1351
1352
case t_POL:
1353
lx = lg(x); s = gen_0;
1354
for (i=2; i<lx; i++) s = gadd(s, gnorml1_fake(gel(x,i)));
1355
break;
1356
1357
case t_VEC: case t_COL: case t_MAT:
1358
lx = lg(x); s = gen_0;
1359
for (i=1; i<lx; i++) s = gadd(s, gnorml1_fake(gel(x,i)));
1360
break;
1361
1362
default: pari_err_TYPE("gnorml1_fake",x);
1363
return NULL; /* LCOV_EXCL_LINE */
1364
}
1365
return gerepileupto(av, s);
1366
}
1367
1368
static void
1369
store(GEN z, GEN *m) { if (!*m || gcmp(z, *m) > 0) *m = z; }
1370
/* compare |x| to *m or |x|^2 to *msq, whichever is easiest, and update
1371
* the pointed value if x is larger */
1372
void
1373
gsupnorm_aux(GEN x, GEN *m, GEN *msq, long prec)
1374
{
1375
long i, lx;
1376
GEN z;
1377
switch(typ(x))
1378
{
1379
case t_COMPLEX: z = cxnorm(x); store(z, msq); return;
1380
case t_QUAD: z = cxquadnorm(x,prec); store(z, msq); return;
1381
case t_INT: case t_REAL: z = mpabs(x); store(z,m); return;
1382
case t_FRAC: z = absfrac(x); store(z,m); return;
1383
1384
case t_POL: lx = lg(x)-1; x++; break;
1385
1386
case t_VEC:
1387
case t_COL:
1388
case t_MAT: lx = lg(x); break;
1389
1390
default: pari_err_TYPE("gsupnorm",x);
1391
return; /* LCOV_EXCL_LINE */
1392
}
1393
for (i=1; i<lx; i++) gsupnorm_aux(gel(x,i), m, msq, prec);
1394
}
1395
GEN
1396
gsupnorm(GEN x, long prec)
1397
{
1398
GEN m = NULL, msq = NULL;
1399
pari_sp av = avma;
1400
gsupnorm_aux(x, &m, &msq, prec);
1401
/* now set m = max (m, sqrt(msq)) */
1402
if (msq) {
1403
msq = gsqrt(msq, prec);
1404
if (!m || gcmp(m, msq) < 0) m = msq;
1405
} else if (!m) m = gen_0;
1406
return gerepilecopy(av, m);
1407
}
1408
1409
/*******************************************************************/
1410
/* */
1411
/* TRACES */
1412
/* */
1413
/*******************************************************************/
1414
GEN
1415
matcompanion(GEN x)
1416
{
1417
long n = degpol(x), j;
1418
GEN y, c;
1419
1420
if (typ(x)!=t_POL) pari_err_TYPE("matcompanion",x);
1421
if (!signe(x)) pari_err_DOMAIN("matcompanion","polynomial","=",gen_0,x);
1422
if (n == 0) return cgetg(1, t_MAT);
1423
1424
y = cgetg(n+1,t_MAT);
1425
for (j=1; j < n; j++) gel(y,j) = col_ei(n, j+1);
1426
c = cgetg(n+1,t_COL); gel(y,n) = c;
1427
if (gequal1(gel(x, n+2)))
1428
for (j=1; j<=n; j++) gel(c,j) = gneg(gel(x,j+1));
1429
else
1430
{ /* not monic. Hardly ever used */
1431
pari_sp av = avma;
1432
GEN d = gclone(gneg(gel(x,n+2)));
1433
set_avma(av);
1434
for (j=1; j<=n; j++) gel(c,j) = gdiv(gel(x,j+1), d);
1435
gunclone(d);
1436
}
1437
return y;
1438
}
1439
1440
GEN
1441
gtrace(GEN x)
1442
{
1443
pari_sp av;
1444
long i, lx, tx = typ(x);
1445
GEN y, z;
1446
1447
switch(tx)
1448
{
1449
case t_INT: case t_REAL: case t_FRAC:
1450
return gmul2n(x,1);
1451
1452
case t_COMPLEX:
1453
return gmul2n(gel(x,1),1);
1454
1455
case t_QUAD:
1456
y = gel(x,1);
1457
if (!gequal0(gel(y,3)))
1458
{ /* assume quad. polynomial is either x^2 + d or x^2 - x + d */
1459
av = avma;
1460
return gerepileupto(av, gadd(gel(x,3), gmul2n(gel(x,2),1)));
1461
}
1462
return gmul2n(gel(x,2),1);
1463
1464
case t_POL:
1465
y = cgetg_copy(x, &lx); y[1] = x[1];
1466
for (i=2; i<lx; i++) gel(y,i) = gtrace(gel(x,i));
1467
return normalizepol_lg(y, lx);
1468
1469
case t_SER:
1470
if (ser_isexactzero(x)) return gcopy(x);
1471
y = cgetg_copy(x, &lx); y[1] = x[1];
1472
for (i=2; i<lx; i++) gel(y,i) = gtrace(gel(x,i));
1473
return normalize(y);
1474
1475
case t_POLMOD:
1476
y = gel(x,1); z = gel(x,2);
1477
if (typ(z) != t_POL || varn(y) != varn(z)) return gmulsg(degpol(y), z);
1478
av = avma;
1479
return gerepileupto(av, RgXQ_trace(z, y));
1480
1481
case t_FFELT:
1482
y=cgetg(3, t_INTMOD);
1483
gel(y,1) = FF_p(x);
1484
gel(y,2) = FF_trace(x);
1485
return y;
1486
1487
case t_RFRAC:
1488
av = avma; return gerepileupto(av, gadd(x, conj_i(x)));
1489
1490
case t_VEC: case t_COL:
1491
y = cgetg_copy(x, &lx);
1492
for (i=1; i<lx; i++) gel(y,i) = gtrace(gel(x,i));
1493
return y;
1494
1495
case t_MAT:
1496
lx = lg(x); if (lx == 1) return gen_0;
1497
/*now lx >= 2*/
1498
if (lx != lgcols(x)) pari_err_DIM("gtrace");
1499
av = avma; return gerepileupto(av, mattrace(x));
1500
}
1501
pari_err_TYPE("gtrace",x);
1502
return NULL; /* LCOV_EXCL_LINE */
1503
}
1504
1505
/* Cholesky decomposition for positive definite matrix a
1506
* [GTM138, Algo 2.7.6, matrix Q]
1507
* If a is not positive definite return NULL. */
1508
GEN
1509
qfgaussred_positive(GEN a)
1510
{
1511
pari_sp av = avma;
1512
GEN b;
1513
long i,j,k, n = lg(a);
1514
1515
if (typ(a)!=t_MAT) pari_err_TYPE("qfgaussred_positive",a);
1516
if (n == 1) return cgetg(1, t_MAT);
1517
if (lgcols(a)!=n) pari_err_DIM("qfgaussred_positive");
1518
b = cgetg(n,t_MAT);
1519
for (j=1; j<n; j++)
1520
{
1521
GEN p1=cgetg(n,t_COL), p2=gel(a,j);
1522
1523
gel(b,j) = p1;
1524
for (i=1; i<=j; i++) gel(p1,i) = gel(p2,i);
1525
for ( ; i<n ; i++) gel(p1,i) = gen_0;
1526
}
1527
for (k=1; k<n; k++)
1528
{
1529
GEN bk, p = gcoeff(b,k,k), invp;
1530
if (gsigne(p)<=0) return gc_NULL(av); /* not positive definite */
1531
invp = ginv(p);
1532
bk = row(b, k);
1533
for (i=k+1; i<n; i++) gcoeff(b,k,i) = gmul(gel(bk,i), invp);
1534
for (i=k+1; i<n; i++)
1535
{
1536
GEN c = gel(bk, i);
1537
for (j=i; j<n; j++)
1538
gcoeff(b,i,j) = gsub(gcoeff(b,i,j), gmul(c,gcoeff(b,k,j)));
1539
}
1540
if (gc_needed(av,1))
1541
{
1542
if (DEBUGMEM>1) pari_warn(warnmem,"qfgaussred_positive");
1543
b=gerepilecopy(av,b);
1544
}
1545
}
1546
return gerepilecopy(av,b);
1547
}
1548
1549
/* Maximal pivot strategy: x is a suitable pivot if it is non zero and either
1550
* - an exact type, or
1551
* - it is maximal among remaining nonzero (t_REAL) pivots */
1552
static int
1553
suitable(GEN x, long k, GEN *pp, long *pi)
1554
{
1555
long t = typ(x);
1556
switch(t)
1557
{
1558
case t_INT: return signe(x) != 0;
1559
case t_FRAC: return 1;
1560
case t_REAL: {
1561
GEN p = *pp;
1562
if (signe(x) && (!p || abscmprr(p, x) < 0)) { *pp = x; *pi = k; }
1563
return 0;
1564
}
1565
default: return !gequal0(x);
1566
}
1567
}
1568
1569
/* Gauss reduction (arbitrary symetric matrix, only the part above the
1570
* diagonal is considered). If signature is nonzero, return only the
1571
* signature, in which case gsigne() should be defined for elements of a. */
1572
static GEN
1573
gaussred(GEN a, long signature)
1574
{
1575
GEN r, ak, al;
1576
pari_sp av = avma, av1;
1577
long n = lg(a), i, j, k, l, sp, sn, t;
1578
1579
if (typ(a) != t_MAT) pari_err_TYPE("gaussred",a);
1580
if (n == 1) return signature? mkvec2(gen_0, gen_0): cgetg(1, t_MAT);
1581
if (lgcols(a) != n) pari_err_DIM("gaussred");
1582
n--;
1583
r = const_vecsmall(n, 1); av1= avma;
1584
a = RgM_shallowcopy(a);
1585
t = n; sp = sn = 0;
1586
while (t)
1587
{
1588
long pind = 0;
1589
GEN invp, p = NULL;
1590
k=1; while (k<=n && (!r[k] || !suitable(gcoeff(a,k,k), k, &p, &pind))) k++;
1591
if (k > n && p) k = pind;
1592
if (k <= n)
1593
{
1594
p = gcoeff(a,k,k); invp = ginv(p); /* != 0 */
1595
if (signature) { /* skip if (!signature): gsigne may fail ! */
1596
if (gsigne(p) > 0) sp++; else sn++;
1597
}
1598
r[k] = 0; t--;
1599
ak = row(a, k);
1600
for (i=1; i<=n; i++)
1601
gcoeff(a,k,i) = r[i]? gmul(gcoeff(a,k,i), invp): gen_0;
1602
1603
for (i=1; i<=n; i++) if (r[i])
1604
{
1605
GEN c = gel(ak,i); /* - p * a[k,i] */
1606
if (gequal0(c)) continue;
1607
for (j=1; j<=n; j++) if (r[j])
1608
gcoeff(a,i,j) = gsub(gcoeff(a,i,j), gmul(c,gcoeff(a,k,j)));
1609
}
1610
gcoeff(a,k,k) = p;
1611
if (gc_needed(av1,1))
1612
{
1613
if(DEBUGMEM>1) pari_warn(warnmem,"gaussred (t = %ld)", t);
1614
a = gerepilecopy(av1, a);
1615
}
1616
}
1617
else
1618
{ /* all remaining diagonal coeffs are currently 0 */
1619
for (k=1; k<=n; k++) if (r[k])
1620
{
1621
l=k+1; while (l<=n && (!r[l] || !suitable(gcoeff(a,k,l), l, &p, &pind))) l++;
1622
if (l > n && p) l = pind;
1623
if (l > n) continue;
1624
1625
p = gcoeff(a,k,l); invp = ginv(p);
1626
sp++; sn++;
1627
r[k] = r[l] = 0; t -= 2;
1628
ak = row(a, k);
1629
al = row(a, l);
1630
for (i=1; i<=n; i++) if (r[i])
1631
{
1632
gcoeff(a,k,i) = gmul(gcoeff(a,k,i), invp);
1633
gcoeff(a,l,i) = gmul(gcoeff(a,l,i), invp);
1634
} else {
1635
gcoeff(a,k,i) = gen_0;
1636
gcoeff(a,l,i) = gen_0;
1637
}
1638
1639
for (i=1; i<=n; i++) if (r[i])
1640
{ /* c = a[k,i] * p, d = a[l,i] * p; */
1641
GEN c = gel(ak,i), d = gel(al,i);
1642
for (j=1; j<=n; j++) if (r[j])
1643
gcoeff(a,i,j) = gsub(gcoeff(a,i,j),
1644
gadd(gmul(gcoeff(a,l,j), c),
1645
gmul(gcoeff(a,k,j), d)));
1646
}
1647
for (i=1; i<=n; i++) if (r[i])
1648
{
1649
GEN c = gcoeff(a,k,i), d = gcoeff(a,l,i);
1650
gcoeff(a,k,i) = gadd(c, d);
1651
gcoeff(a,l,i) = gsub(c, d);
1652
}
1653
gcoeff(a,k,l) = gen_1;
1654
gcoeff(a,l,k) = gen_m1;
1655
gcoeff(a,k,k) = gmul2n(p,-1);
1656
gcoeff(a,l,l) = gneg(gcoeff(a,k,k));
1657
if (gc_needed(av1,1))
1658
{
1659
if(DEBUGMEM>1) pari_warn(warnmem,"gaussred");
1660
a = gerepilecopy(av1, a);
1661
}
1662
break;
1663
}
1664
if (k > n) break;
1665
}
1666
}
1667
if (!signature) return gerepilecopy(av, a);
1668
set_avma(av); return mkvec2s(sp, sn);
1669
}
1670
1671
GEN
1672
qfgaussred(GEN a) { return gaussred(a,0); }
1673
1674
GEN
1675
qfsign(GEN a) { return gaussred(a,1); }
1676
1677
/* x -= s(y+u*x) */
1678
/* y += s(x-u*y), simultaneously */
1679
static void
1680
rot(GEN x, GEN y, GEN s, GEN u) {
1681
GEN x1 = subrr(x, mulrr(s,addrr(y,mulrr(u,x))));
1682
GEN y1 = addrr(y, mulrr(s,subrr(x,mulrr(u,y))));
1683
affrr(x1,x);
1684
affrr(y1,y);
1685
}
1686
1687
/* Diagonalization of a REAL symetric matrix. Return a vector [L, r]:
1688
* L = vector of eigenvalues
1689
* r = matrix of eigenvectors */
1690
GEN
1691
jacobi(GEN a, long prec)
1692
{
1693
pari_sp av;
1694
long de, e, e1, e2, i, j, p, q, l = lg(a);
1695
GEN c, ja, L, r, L2, r2, unr;
1696
1697
if (typ(a) != t_MAT) pari_err_TYPE("jacobi",a);
1698
ja = cgetg(3,t_VEC);
1699
L = cgetg(l,t_COL); gel(ja,1) = L;
1700
r = cgetg(l,t_MAT); gel(ja,2) = r;
1701
if (l == 1) return ja;
1702
if (lgcols(a) != l) pari_err_DIM("jacobi");
1703
1704
e1 = HIGHEXPOBIT-1;
1705
for (j=1; j<l; j++)
1706
{
1707
GEN z = gtofp(gcoeff(a,j,j), prec);
1708
gel(L,j) = z;
1709
e = expo(z); if (e < e1) e1 = e;
1710
}
1711
for (j=1; j<l; j++)
1712
{
1713
gel(r,j) = cgetg(l,t_COL);
1714
for (i=1; i<l; i++) gcoeff(r,i,j) = utor(i==j? 1: 0, prec);
1715
}
1716
av = avma;
1717
1718
e2 = -(long)HIGHEXPOBIT; p = q = 1;
1719
c = cgetg(l,t_MAT);
1720
for (j=1; j<l; j++)
1721
{
1722
gel(c,j) = cgetg(j,t_COL);
1723
for (i=1; i<j; i++)
1724
{
1725
GEN z = gtofp(gcoeff(a,i,j), prec);
1726
gcoeff(c,i,j) = z;
1727
if (!signe(z)) continue;
1728
e = expo(z); if (e > e2) { e2 = e; p = i; q = j; }
1729
}
1730
}
1731
a = c; unr = real_1(prec);
1732
de = prec2nbits(prec);
1733
1734
/* e1 = min expo(a[i,i])
1735
* e2 = max expo(a[i,j]), i != j */
1736
while (e1-e2 < de)
1737
{
1738
pari_sp av2 = avma;
1739
GEN x, y, t, c, s, u;
1740
/* compute attached rotation in the plane formed by basis vectors number
1741
* p and q */
1742
x = subrr(gel(L,q),gel(L,p));
1743
if (signe(x))
1744
{
1745
x = divrr(x, shiftr(gcoeff(a,p,q),1));
1746
y = sqrtr(addrr(unr, sqrr(x)));
1747
t = invr((signe(x)>0)? addrr(x,y): subrr(x,y));
1748
}
1749
else
1750
y = t = unr;
1751
c = sqrtr(addrr(unr,sqrr(t)));
1752
s = divrr(t,c);
1753
u = divrr(t,addrr(unr,c));
1754
1755
/* compute successive transforms of a and the matrix of accumulated
1756
* rotations (r) */
1757
for (i=1; i<p; i++) rot(gcoeff(a,i,p), gcoeff(a,i,q), s,u);
1758
for (i=p+1; i<q; i++) rot(gcoeff(a,p,i), gcoeff(a,i,q), s,u);
1759
for (i=q+1; i<l; i++) rot(gcoeff(a,p,i), gcoeff(a,q,i), s,u);
1760
y = gcoeff(a,p,q);
1761
t = mulrr(t, y); shiftr_inplace(y, -de - 1);
1762
x = gel(L,p); subrrz(x,t, x);
1763
y = gel(L,q); addrrz(y,t, y);
1764
for (i=1; i<l; i++) rot(gcoeff(r,i,p), gcoeff(r,i,q), s,u);
1765
1766
e2 = -(long)HIGHEXPOBIT; p = q = 1;
1767
for (j=1; j<l; j++)
1768
{
1769
for (i=1; i<j; i++)
1770
{
1771
GEN z = gcoeff(a,i,j);
1772
if (!signe(z)) continue;
1773
e = expo(z); if (e > e2) { e2=e; p=i; q=j; }
1774
}
1775
for (i=j+1; i<l; i++)
1776
{
1777
GEN z = gcoeff(a,j,i);
1778
if (!signe(z)) continue;
1779
e = expo(z); if (e > e2) { e2=e; p=j; q=i; }
1780
}
1781
}
1782
set_avma(av2);
1783
}
1784
/* sort eigenvalues from smallest to largest */
1785
c = indexsort(L);
1786
r2 = vecpermute(r, c); for (i=1; i<l; i++) gel(r,i) = gel(r2,i);
1787
L2 = vecpermute(L, c); for (i=1; i<l; i++) gel(L,i) = gel(L2,i);
1788
set_avma(av); return ja;
1789
}
1790
1791
/*************************************************************************/
1792
/** **/
1793
/** Q-vector space -> Z-modules **/
1794
/** **/
1795
/*************************************************************************/
1796
1797
GEN
1798
matrixqz0(GEN x,GEN p)
1799
{
1800
if (typ(x) != t_MAT) pari_err_TYPE("matrixqz",x);
1801
if (!p) return QM_minors_coprime(x,NULL);
1802
if (typ(p) != t_INT) pari_err_TYPE("matrixqz",p);
1803
if (signe(p)>=0) return QM_minors_coprime(x,p);
1804
if (!RgM_is_QM(x)) pari_err_TYPE("matrixqz", x);
1805
if (absequaliu(p,1)) return QM_ImZ(x); /* p = -1 */
1806
if (absequaliu(p,2)) return QM_ImQ(x); /* p = -2 */
1807
pari_err_FLAG("QM_minors_coprime");
1808
return NULL; /* LCOV_EXCL_LINE */
1809
}
1810
1811
GEN
1812
QM_minors_coprime(GEN x, GEN D)
1813
{
1814
pari_sp av = avma, av1;
1815
long i, j, m, n, lP;
1816
GEN P, y;
1817
1818
n = lg(x)-1; if (!n) return gcopy(x);
1819
m = nbrows(x);
1820
if (n > m) pari_err_DOMAIN("QM_minors_coprime","n",">",strtoGENstr("m"),x);
1821
y = x; x = cgetg(n+1,t_MAT);
1822
for (j=1; j<=n; j++)
1823
{
1824
gel(x,j) = Q_primpart(gel(y,j));
1825
RgV_check_ZV(gel(x,j), "QM_minors_coprime");
1826
}
1827
/* x now a ZM */
1828
if (n==m)
1829
{
1830
if (gequal0(ZM_det(x)))
1831
pari_err_DOMAIN("QM_minors_coprime", "rank(A)", "<",stoi(n),x);
1832
set_avma(av); return matid(n);
1833
}
1834
/* m > n */
1835
if (!D || gequal0(D))
1836
{
1837
pari_sp av2 = avma;
1838
D = ZM_detmult(shallowtrans(x));
1839
if (is_pm1(D)) { set_avma(av2); return ZM_copy(x); }
1840
}
1841
P = gel(Z_factor(D), 1); lP = lg(P);
1842
av1 = avma;
1843
for (i=1; i < lP; i++)
1844
{
1845
GEN p = gel(P,i), pov2 = shifti(p, -1);
1846
for(;;)
1847
{
1848
GEN N, M = FpM_ker(x, p);
1849
long lM = lg(M);
1850
if (lM==1) break;
1851
1852
FpM_center_inplace(M, p, pov2);
1853
N = ZM_Z_divexact(ZM_mul(x,M), p);
1854
for (j=1; j<lM; j++)
1855
{
1856
long k = n; while (!signe(gcoeff(M,k,j))) k--;
1857
gel(x,k) = gel(N,j);
1858
}
1859
if (gc_needed(av1,1))
1860
{
1861
if (DEBUGMEM>1) pari_warn(warnmem,"QM_minors_coprime, p = %Ps", p);
1862
x = gerepilecopy(av1, x); pov2 = shifti(p, -1);
1863
}
1864
}
1865
}
1866
return gerepilecopy(av, x);
1867
}
1868
1869
static GEN
1870
QM_ImZ_all_i(GEN A, GEN *U, long remove, long hnf, long linindep)
1871
{
1872
GEN V = NULL, D;
1873
A = Q_remove_denom(A,&D);
1874
if (D)
1875
{
1876
long l, lA;
1877
V = matkermod(A,D,NULL);
1878
l = lg(V); lA = lg(A);
1879
if (l == 1) V = scalarmat_shallow(D, lA-1);
1880
else
1881
{
1882
if (l < lA) V = hnfmodid(V,D);
1883
A = ZM_Z_divexact(ZM_mul(A, V), D);
1884
}
1885
}
1886
if (!linindep && ZM_rank(A)==lg(A)-1) linindep = 1;
1887
if (hnf || !linindep) A = ZM_hnflll(A, U, remove);
1888
if (U && V)
1889
{
1890
if (hnf) *U = ZM_mul(V,*U);
1891
else *U = V;
1892
}
1893
return A;
1894
}
1895
GEN
1896
QM_ImZ_all(GEN x, GEN *U, long remove, long hnf)
1897
{
1898
pari_sp av = avma;
1899
x = QM_ImZ_all_i(x, U, remove, hnf, 0);
1900
if (U) gerepileall(av, 2, &x, &U); else x = gerepilecopy(av, x);
1901
return x;
1902
}
1903
GEN
1904
QM_ImZ_hnfall(GEN x, GEN *U, long remove) { return QM_ImZ_all(x, U, remove, 1); }
1905
GEN
1906
QM_ImZ_hnf(GEN x) { return QM_ImZ_hnfall(x, NULL, 1); }
1907
GEN
1908
QM_ImZ(GEN x) { return QM_ImZ_all(x, NULL, 1, 0); }
1909
1910
GEN
1911
QM_ImQ_all(GEN x, GEN *U, long remove, long hnf)
1912
{
1913
pari_sp av = avma;
1914
long i, n = lg(x), m;
1915
GEN ir, V, D, c, K = NULL;
1916
1917
if (U) *U = matid(n-1);
1918
if (n==1) return gcopy(x);
1919
m = lg(gel(x,1));
1920
1921
x = RgM_shallowcopy(x);
1922
for (i=1; i<n; i++)
1923
{
1924
gel(x,i) = Q_primitive_part(gel(x,i), &c);
1925
if (U && c && signe(c)) gcoeff(*U,i,i) = ginv(c);
1926
}
1927
1928
ir = ZM_indexrank(x);
1929
if (U)
1930
{
1931
*U = vecpermute(*U, gel(ir,2));
1932
if (remove < 2) K = ZM_ker(x);
1933
}
1934
x = vecpermute(x, gel(ir,2));
1935
1936
D = absi(ZM_det(rowpermute(x,gel(ir,1))));
1937
x = RgM_Rg_div(x, D);
1938
x = QM_ImZ_all_i(x, U? &V: NULL, remove, hnf, 1);
1939
1940
if (U)
1941
{
1942
*U = RgM_Rg_div(RgM_mul(*U,V),D);
1943
if (remove < 2) *U = shallowconcat(K,*U);
1944
if (!remove) x = shallowconcat(zeromatcopy(m-1,lg(K)-1), x);
1945
gerepileall(av, 2, &x, U);
1946
}
1947
else x = gerepilecopy(av,x);
1948
return x;
1949
}
1950
GEN
1951
QM_ImQ_hnfall(GEN x, GEN *U, long remove) { return QM_ImQ_all(x, U, remove, 1); }
1952
GEN
1953
QM_ImQ_hnf(GEN x) { return QM_ImQ_hnfall(x, NULL, 1); }
1954
GEN
1955
QM_ImQ(GEN x) { return QM_ImQ_all(x, NULL, 1, 0); }
1956
1957
GEN
1958
intersect(GEN x, GEN y)
1959
{
1960
long j, lx = lg(x);
1961
pari_sp av;
1962
GEN z;
1963
1964
if (typ(x)!=t_MAT) pari_err_TYPE("intersect",x);
1965
if (typ(y)!=t_MAT) pari_err_TYPE("intersect",y);
1966
if (lx==1 || lg(y)==1) return cgetg(1,t_MAT);
1967
1968
av = avma; z = ker(shallowconcat(x,y));
1969
for (j=lg(z)-1; j; j--) setlg(z[j], lx);
1970
return gerepileupto(av, image(RgM_mul(x,z)));
1971
}
1972
1973