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