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