Testing latest pari + WASM + node.js... and it works?! Wow.
License: GPL3
ubuntu2004
/* Copyright (C) 2000 The PARI group.12This file is part of the PARI/GP package.34PARI/GP is free software; you can redistribute it and/or modify it under the5terms of the GNU General Public License as published by the Free Software6Foundation; either version 2 of the License, or (at your option) any later7version. It is distributed in the hope that it will be useful, but WITHOUT8ANY WARRANTY WHATSOEVER.910Check the License for details. You should have received a copy of it, along11with the package; see the file 'COPYING'. If not, write to the Free Software12Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */1314/********************************************************************/15/** **/16/** LINEAR ALGEBRA **/17/** (second part) **/18/** **/19/********************************************************************/20#include "pari.h"21#include "paripriv.h"2223#define DEBUGLEVEL DEBUGLEVEL_mat2425/*******************************************************************/26/* */27/* CHARACTERISTIC POLYNOMIAL */28/* */29/*******************************************************************/3031static GEN32Flm_charpoly_i(GEN x, ulong p)33{34long lx = lg(x), r, i;35GEN H, y = cgetg(lx+1, t_VEC);36gel(y,1) = pol1_Flx(0); H = Flm_hess(x, p);37for (r = 1; r < lx; r++)38{39pari_sp av2 = avma;40ulong a = 1;41GEN z, b = zero_Flx(0);42for (i = r-1; i; i--)43{44a = Fl_mul(a, ucoeff(H,i+1,i), p);45if (!a) break;46b = Flx_add(b, Flx_Fl_mul(gel(y,i), Fl_mul(a,ucoeff(H,i,r),p), p), p);47}48z = Flx_sub(Flx_shift(gel(y,r), 1),49Flx_Fl_mul(gel(y,r), ucoeff(H,r,r), p), p);50/* (X - H[r,r])y[r] - b */51gel(y,r+1) = gerepileuptoleaf(av2, Flx_sub(z, b, p));52}53return gel(y,lx);54}5556GEN57Flm_charpoly(GEN x, ulong p)58{59pari_sp av = avma;60return gerepileuptoleaf(av, Flm_charpoly_i(x,p));61}6263GEN64FpM_charpoly(GEN x, GEN p)65{66pari_sp av = avma;67long lx, r, i;68GEN y, H;6970if (lgefint(p) == 3)71{72ulong pp = p[2];73y = Flx_to_ZX(Flm_charpoly_i(ZM_to_Flm(x,pp), pp));74return gerepileupto(av, y);75}76lx = lg(x); y = cgetg(lx+1, t_VEC);77gel(y,1) = pol_1(0); H = FpM_hess(x, p);78for (r = 1; r < lx; r++)79{80pari_sp av2 = avma;81GEN z, a = gen_1, b = pol_0(0);82for (i = r-1; i; i--)83{84a = Fp_mul(a, gcoeff(H,i+1,i), p);85if (!signe(a)) break;86b = ZX_add(b, ZX_Z_mul(gel(y,i), Fp_mul(a,gcoeff(H,i,r),p)));87}88b = FpX_red(b, p);89z = FpX_sub(RgX_shift_shallow(gel(y,r), 1),90FpX_Fp_mul(gel(y,r), gcoeff(H,r,r), p), p);91z = FpX_sub(z,b,p);92if (r+1 == lx) { gel(y,lx) = z; break; }93gel(y,r+1) = gerepileupto(av2, z); /* (X - H[r,r])y[r] - b */94}95return gerepileupto(av, gel(y,lx));96}9798GEN99charpoly0(GEN x, long v, long flag)100{101if (v<0) v = 0;102switch(flag)103{104case 0: return caradj(x,v,NULL);105case 1: return caract(x,v);106case 2: return carhess(x,v);107case 3: return carberkowitz(x,v);108case 4:109if (typ(x) != t_MAT) pari_err_TYPE("charpoly",x);110RgM_check_ZM(x, "charpoly");111x = ZM_charpoly(x); setvarn(x, v); return x;112case 5:113return charpoly(x, v);114}115pari_err_FLAG("charpoly");116return NULL; /* LCOV_EXCL_LINE */117}118119/* characteristic pol. Easy cases. Return NULL in case it's not so easy. */120static GEN121easychar(GEN x, long v)122{123pari_sp av;124long lx;125GEN p1;126127switch(typ(x))128{129case t_INT: case t_REAL: case t_INTMOD:130case t_FRAC: case t_PADIC:131p1=cgetg(4,t_POL);132p1[1]=evalsigne(1) | evalvarn(v);133gel(p1,2) = gneg(x); gel(p1,3) = gen_1;134return p1;135136case t_COMPLEX: case t_QUAD:137p1 = cgetg(5,t_POL);138p1[1] = evalsigne(1) | evalvarn(v);139gel(p1,2) = gnorm(x); av = avma;140gel(p1,3) = gerepileupto(av, gneg(gtrace(x)));141gel(p1,4) = gen_1; return p1;142143case t_FFELT: {144pari_sp ltop=avma;145p1 = FpX_to_mod(FF_charpoly(x), FF_p_i(x));146setvarn(p1,v); return gerepileupto(ltop,p1);147}148149case t_POLMOD:150{151GEN A = gel(x,2), T = gel(x,1);152if (typ(A)==t_POL && RgX_is_QX(A) && RgX_is_ZX(T))153return QXQ_charpoly(A, T, v);154else155return RgXQ_charpoly(A, T, v);156}157case t_MAT:158lx=lg(x);159if (lx==1) return pol_1(v);160if (lgcols(x) != lx) break;161return NULL;162}163pari_err_TYPE("easychar",x);164return NULL; /* LCOV_EXCL_LINE */165}166/* compute charpoly by mapping to Fp first, return lift to Z */167static GEN168RgM_Fp_charpoly(GEN x, GEN p, long v)169{170GEN T;171if (lgefint(p) == 3)172{173ulong pp = itou(p);174T = Flm_charpoly_i(RgM_to_Flm(x, pp), pp);175T = Flx_to_ZX(T);176}177else178T = FpM_charpoly(RgM_to_FpM(x, p), p);179setvarn(T, v); return T;180}181GEN182charpoly(GEN x, long v)183{184GEN T, p = NULL;185long prec;186if ((T = easychar(x,v))) return T;187switch(RgM_type(x, &p,&T,&prec))188{189case t_INT:190T = ZM_charpoly(x); setvarn(T, v); break;191case t_INTMOD:192if (!BPSW_psp(p)) T = carberkowitz(x, v);193else194{195pari_sp av = avma;196T = RgM_Fp_charpoly(x,p,v);197T = gerepileupto(av, FpX_to_mod(T,p));198}199break;200case t_REAL:201case t_COMPLEX:202case t_PADIC:203T = carhess(x, v);204break;205default:206T = carberkowitz(x, v);207}208return T;209}210211/* We possibly worked with an "invalid" polynomial p, satisfying212* varn(p) > gvar2(p). Fix this. */213static GEN214fix_pol(pari_sp av, GEN p)215{216long w = gvar2(p), v = varn(p);217if (w == v) pari_err_PRIORITY("charpoly", p, "=", w);218if (varncmp(w,v) < 0) p = gerepileupto(av, poleval(p, pol_x(v)));219return p;220}221GEN222caract(GEN x, long v)223{224pari_sp av = avma;225GEN T, C, x_k, Q;226long k, n;227228if ((T = easychar(x,v))) return T;229230n = lg(x)-1;231if (n == 1) return fix_pol(av, deg1pol(gen_1, gneg(gcoeff(x,1,1)), v));232233x_k = pol_x(v); /* to be modified in place */234T = scalarpol(det(x), v); C = utoineg(n); Q = pol_x(v);235for (k=1; k<=n; k++)236{237GEN mk = utoineg(k), d;238gel(x_k,2) = mk;239d = det(RgM_Rg_add_shallow(x, mk));240T = RgX_add(RgX_mul(T, x_k), RgX_Rg_mul(Q, gmul(C, d)));241if (k == n) break;242243Q = RgX_mul(Q, x_k);244C = diviuexact(mulsi(k-n,C), k+1); /* (-1)^k binomial(n,k) */245}246return fix_pol(av, RgX_Rg_div(T, mpfact(n)));247}248249/* C = charpoly(x, v) */250static GEN251RgM_adj_from_char(GEN x, long v, GEN C)252{253if (varn(C) != v) /* problem with variable priorities */254{255C = gdiv(gsub(C, gsubst(C, v, gen_0)), pol_x(v));256if (odd(lg(x))) C = RgX_neg(C); /* even dimension */257return gsubst(C, v, x);258}259else260{261C = RgX_shift_shallow(C, -1);262if (odd(lg(x))) C = RgX_neg(C); /* even dimension */263return RgX_RgM_eval(C, x);264}265}266/* assume x square matrice */267static GEN268mattrace(GEN x)269{270long i, lx = lg(x);271GEN t;272if (lx < 3) return lx == 1? gen_0: gcopy(gcoeff(x,1,1));273t = gcoeff(x,1,1);274for (i = 2; i < lx; i++) t = gadd(t, gcoeff(x,i,i));275return t;276}277static int278bad_char(GEN q, long n)279{280forprime_t S;281ulong p;282if (!signe(q)) return 0;283(void)u_forprime_init(&S, 2, n);284while ((p = u_forprime_next(&S)))285if (!umodiu(q, p)) return 1;286return 0;287}288/* Using traces: return the characteristic polynomial of x (in variable v).289* If py != NULL, the adjoint matrix is put there. */290GEN291caradj(GEN x, long v, GEN *py)292{293pari_sp av, av0;294long i, k, n;295GEN T, y, t;296297if ((T = easychar(x, v)))298{299if (py)300{301if (typ(x) != t_MAT) pari_err_TYPE("matadjoint",x);302*py = cgetg(1,t_MAT);303}304return T;305}306307n = lg(x)-1; av0 = avma;308T = cgetg(n+3,t_POL); T[1] = evalsigne(1) | evalvarn(v);309gel(T,n+2) = gen_1;310if (!n) { if (py) *py = cgetg(1,t_MAT); return T; }311av = avma; t = gerepileupto(av, gneg(mattrace(x)));312gel(T,n+1) = t;313if (n == 1) {314T = fix_pol(av0, T);315if (py) *py = matid(1); return T;316}317if (n == 2) {318GEN a = gcoeff(x,1,1), b = gcoeff(x,1,2);319GEN c = gcoeff(x,2,1), d = gcoeff(x,2,2);320av = avma;321gel(T,2) = gerepileupto(av, gsub(gmul(a,d), gmul(b,c)));322T = fix_pol(av0, T);323if (py) {324y = cgetg(3, t_MAT);325gel(y,1) = mkcol2(gcopy(d), gneg(c));326gel(y,2) = mkcol2(gneg(b), gcopy(a));327*py = y;328}329return T;330}331/* l > 3 */332if (bad_char(residual_characteristic(x), n))333{ /* n! not invertible in base ring */334T = charpoly(x, v);335if (!py) return gerepileupto(av, T);336*py = RgM_adj_from_char(x, v, T);337gerepileall(av, 2, &T,py);338return T;339}340av = avma; y = RgM_shallowcopy(x);341for (i = 1; i <= n; i++) gcoeff(y,i,i) = gadd(gcoeff(y,i,i), t);342for (k = 2; k < n; k++)343{344GEN y0 = y;345y = RgM_mul(y, x);346t = gdivgs(mattrace(y), -k);347for (i = 1; i <= n; i++) gcoeff(y,i,i) = gadd(gcoeff(y,i,i), t);348y = gclone(y);349gel(T,n-k+2) = gerepilecopy(av, t); av = avma;350if (k > 2) gunclone(y0);351}352t = gmul(gcoeff(x,1,1),gcoeff(y,1,1));353for (i=2; i<=n; i++) t = gadd(t, gmul(gcoeff(x,1,i),gcoeff(y,i,1)));354gel(T,2) = gerepileupto(av, gneg(t));355T = fix_pol(av0, T);356if (py) *py = odd(n)? gcopy(y): RgM_neg(y);357gunclone(y); return T;358}359360GEN361adj(GEN x)362{363GEN y;364(void)caradj(x, fetch_var(), &y);365(void)delete_var(); return y;366}367368GEN369adjsafe(GEN x)370{371const long v = fetch_var();372pari_sp av = avma;373GEN C, A;374if (typ(x) != t_MAT) pari_err_TYPE("matadjoint",x);375if (lg(x) < 3) return gcopy(x);376C = charpoly(x,v);377A = RgM_adj_from_char(x, v, C);378(void)delete_var(); return gerepileupto(av, A);379}380381GEN382matadjoint0(GEN x, long flag)383{384switch(flag)385{386case 0: return adj(x);387case 1: return adjsafe(x);388}389pari_err_FLAG("matadjoint");390return NULL; /* LCOV_EXCL_LINE */391}392393/*******************************************************************/394/* */395/* Frobenius form */396/* */397/*******************************************************************/398399/* The following section implement a mix of Ozello and Storjohann algorithms400401P. Ozello, doctoral thesis (in French):402Calcul exact des formes de Jordan et de Frobenius d'une matrice, Chapitre 2403http://tel.archives-ouvertes.fr/tel-00323705404405A. Storjohann, Diss. ETH No. 13922406Algorithms for Matrix Canonical Forms, Chapter 9407https://cs.uwaterloo.ca/~astorjoh/diss2up.pdf408409We use Storjohann Lemma 9.14 (step1, step2, step3) Ozello theorem 4,410and Storjohann Lemma 9.18411*/412413/* Elementary transforms */414415/* M <- U^(-1) M U, U = E_{i,j}(k) => U^(-1) = E{i,j}(-k)416* P = U * P */417static void418transL(GEN M, GEN P, GEN k, long i, long j)419{420long l, n = lg(M)-1;421for(l=1; l<=n; l++) /* M[,j]-=k*M[,i] */422gcoeff(M,l,j) = gsub(gcoeff(M,l,j), gmul(gcoeff(M,l,i), k));423for(l=1; l<=n; l++) /* M[i,]+=k*M[j,] */424gcoeff(M,i,l) = gadd(gcoeff(M,i,l), gmul(gcoeff(M,j,l), k));425if (P)426for(l=1; l<=n; l++)427gcoeff(P,i,l) = gadd(gcoeff(P,i,l), gmul(gcoeff(P,j,l), k));428}429430/* j = a or b */431static void432transD(GEN M, GEN P, long a, long b, long j)433{434long l, n;435GEN k = gcoeff(M,a,b), ki;436437if (gequal1(k)) return;438ki = ginv(k); n = lg(M)-1;439for(l=1; l<=n; l++)440if (l!=j)441{442gcoeff(M,l,j) = gmul(gcoeff(M,l,j), k);443gcoeff(M,j,l) = (j==a && l==b)? gen_1: gmul(gcoeff(M,j,l), ki);444}445if (P)446for(l=1; l<=n; l++)447gcoeff(P,j,l) = gmul(gcoeff(P,j,l), ki);448}449450static void451transS(GEN M, GEN P, long i, long j)452{453long l, n = lg(M)-1;454swap(gel(M,i), gel(M,j));455for (l=1; l<=n; l++)456swap(gcoeff(M,i,l), gcoeff(M,j,l));457if (P)458for (l=1; l<=n; l++)459swap(gcoeff(P,i,l), gcoeff(P,j,l));460}461462/* Convert companion matrix to polynomial*/463static GEN464minpoly_polslice(GEN M, long i, long j, long v)465{466long k, d = j+1-i;467GEN P = cgetg(d+3,t_POL);468P[1] = evalsigne(1)|evalvarn(v);469for (k=0; k<d; k++)470gel(P,k+2) = gneg(gcoeff(M,i+k, j));471gel(P,d+2) = gen_1;472return P;473}474475static GEN476minpoly_listpolslice(GEN M, GEN V, long v)477{478long i, n = lg(M)-1, nb = lg(V)-1;479GEN W = cgetg(nb+1, t_VEC);480for (i=1; i<=nb; i++)481gel(W,i) = minpoly_polslice(M, V[i], i < nb? V[i+1]-1: n, v);482return W;483}484485static int486minpoly_dvdslice(GEN M, long i, long j, long k)487{488pari_sp av = avma;489long r = signe(RgX_rem(minpoly_polslice(M, i, j-1, 0),490minpoly_polslice(M, j, k, 0)));491return gc_bool(av, r == 0);492}493494static void495RgM_replace(GEN M, GEN M2)496{497long n = lg(M)-1, m = nbrows(M), i, j;498for(i=1; i<=n; i++)499for(j=1; j<=m; j++)500gcoeff(M, i, j) = gcoeff(M2, i, j);501}502503static void504gerepilemat2_inplace(pari_sp av, GEN M, GEN P)505{506GEN M2 = M, P2 = P;507gerepileall(av, P ? 2: 1, &M2, &P2);508RgM_replace(M, M2);509if (P) RgM_replace(P, P2);510}511512/* Lemma 9.14 */513static long514weakfrobenius_step1(GEN M, GEN P, long j0)515{516pari_sp av = avma;517long n = lg(M)-1, k, j;518for (j = j0; j < n; ++j)519{520if (gequal0(gcoeff(M, j+1, j)))521{522for (k = j+2; k <= n; ++k)523if (!gequal0(gcoeff(M,k,j))) break;524if (k > n) return j;525transS(M, P, k, j+1);526}527transD(M, P, j+1, j, j+1);528/* Now M[j+1,j] = 1 */529for (k = 1; k <= n; ++k)530if (k != j+1 && !gequal0(gcoeff(M,k,j))) /* zero M[k,j] */531{532transL(M, P, gneg(gcoeff(M,k,j)), k, j+1);533gcoeff(M,k,j) = gen_0; /* avoid approximate 0 */534}535if (gc_needed(av,1))536{537if (DEBUGMEM > 1)538pari_warn(warnmem,"RgM_minpoly stage 1: j0=%ld, j=%ld", j0, j);539gerepilemat2_inplace(av, M, P);540}541}542return n;543}544545static void546weakfrobenius_step2(GEN M, GEN P, long j)547{548pari_sp av = avma;549long i, k, n = lg(M)-1;550for(i=j; i>=2; i--)551{552for(k=j+1; k<=n; k++)553if (!gequal0(gcoeff(M,i,k)))554transL(M, P, gcoeff(M,i,k), i-1, k);555if (gc_needed(av,1))556{557if (DEBUGMEM > 1)558pari_warn(warnmem,"RgM_minpoly stage 2: j=%ld, i=%ld", j, i);559gerepilemat2_inplace(av, M, P);560}561}562}563564static long565weakfrobenius_step3(GEN M, GEN P, long j0, long j)566{567long i, k, n = lg(M)-1;568if (j == n) return 0;569if (gequal0(gcoeff(M, j0, j+1)))570{571for (k=j+2; k<=n; k++)572if (!gequal0(gcoeff(M, j0, k))) break;573if (k > n) return 0;574transS(M, P, k, j+1);575}576transD(M, P, j0, j+1, j+1);577for (i=j+2; i<=n; i++)578if (!gequal0(gcoeff(M, j0, i)))579transL(M, P, gcoeff(M, j0, i),j+1, i);580return 1;581}582583/* flag: 0 -> full Frobenius from , 1 -> weak Frobenius form */584GEN585RgM_Frobenius(GEN M, long flag, GEN *pt_P, GEN *pt_v)586{587pari_sp av = avma, av2, ltop;588long n = lg(M)-1, eps, j0 = 1, nb = 0;589GEN v, P;590v = cgetg(n+1, t_VECSMALL);591ltop = avma;592P = pt_P ? matid(n): NULL;593M = RgM_shallowcopy(M);594av2 = avma;595while (j0 <= n)596{597long j = weakfrobenius_step1(M, P, j0);598weakfrobenius_step2(M, P, j);599eps = weakfrobenius_step3(M, P, j0, j);600if (eps == 0)601{602v[++nb] = j0;603if (flag == 0 && nb > 1 && !minpoly_dvdslice(M, v[nb-1], j0, j))604{605j = j0; j0 = v[nb-1]; nb -= 2;606transL(M, P, gen_1, j, j0); /*lemma 9.18*/607} else608j0 = j+1;609}610else611transS(M, P, j0, j+1); /*theorem 4*/612if (gc_needed(av,1))613{614if (DEBUGMEM > 1)615pari_warn(warnmem,"weakfrobenius j0=%ld",j0);616gerepilemat2_inplace(av2, M, P);617}618}619fixlg(v, nb+1);620if (pt_v) *pt_v = v;621gerepileall(pt_v ? ltop: av, P? 2: 1, &M, &P);622if (pt_P) *pt_P = P;623return M;624}625626static GEN627RgM_minpoly(GEN M, long v)628{629pari_sp av = avma;630GEN V, W;631if (lg(M) == 1) return pol_1(v);632M = RgM_Frobenius(M, 1, NULL, &V);633W = minpoly_listpolslice(M, V, v);634if (varncmp(v,gvar2(W)) >= 0)635pari_err_PRIORITY("matfrobenius", M, "<=", v);636return gerepileupto(av, RgX_normalize(glcm0(W, NULL)));637}638639GEN640Frobeniusform(GEN V, long n)641{642long i, j, k;643GEN M = zeromatcopy(n,n);644for (k=1,i=1;i<lg(V);i++,k++)645{646GEN P = gel(V,i);647long d = degpol(P);648if (k+d-1 > n) pari_err_PREC("matfrobenius");649for (j=0; j<d-1; j++, k++) gcoeff(M,k+1,k) = gen_1;650for (j=0; j<d; j++) gcoeff(M,k-j,k) = gneg(gel(P, 1+d-j));651}652return M;653}654655GEN656matfrobenius(GEN M, long flag, long v)657{658long n;659if (typ(M)!=t_MAT) pari_err_TYPE("matfrobenius",M);660if (v < 0) v = 0;661n = lg(M)-1;662if (n && lgcols(M)!=n+1) pari_err_DIM("matfrobenius");663if (flag > 2) pari_err_FLAG("matfrobenius");664switch (flag)665{666case 0:667return RgM_Frobenius(M, 0, NULL, NULL);668case 1:669{670pari_sp av = avma;671GEN V, W, F;672F = RgM_Frobenius(M, 0, NULL, &V);673W = minpoly_listpolslice(F, V, v);674if (varncmp(v, gvar2(W)) >= 0)675pari_err_PRIORITY("matfrobenius", M, "<=", v);676return gerepileupto(av, W);677}678case 2:679{680GEN P, F, R = cgetg(3, t_VEC);681F = RgM_Frobenius(M, 0, &P, NULL);682gel(R,1) = F; gel(R,2) = P;683return R;684}685default:686pari_err_FLAG("matfrobenius");687}688return NULL; /*LCOV_EXCL_LINE*/689}690691/*******************************************************************/692/* */693/* MINIMAL POLYNOMIAL */694/* */695/*******************************************************************/696697static GEN698RgXQ_minpoly_naive(GEN y, GEN P)699{700long n = lgpol(P);701GEN M = ker(RgXQ_matrix_pow(y,n,n,P));702return content(RgM_to_RgXV(M,varn(P)));703}704705static GEN706RgXQ_minpoly_FpXQ(GEN x, GEN y, GEN p)707{708pari_sp av = avma;709GEN r;710if (lgefint(p) == 3)711{712ulong pp = uel(p, 2);713r = Flx_to_ZX_inplace(Flxq_minpoly(RgX_to_Flx(x, pp),714RgX_to_Flx(y, pp), pp));715}716else717r = FpXQ_minpoly(RgX_to_FpX(x, p), RgX_to_FpX(y, p), p);718return gerepileupto(av, FpX_to_mod(r, p));719}720721static GEN722RgXQ_minpoly_FpXQXQ(GEN x, GEN S, GEN pol, GEN p)723{724pari_sp av = avma;725GEN r;726GEN T = RgX_to_FpX(pol, p);727if (signe(T)==0) pari_err_OP("minpoly",x,S);728if (lgefint(p) == 3)729{730ulong pp = uel(p, 2);731GEN Tp = ZX_to_Flx(T, pp);732r = FlxX_to_ZXX(FlxqXQ_minpoly(RgX_to_FlxqX(x, Tp, pp),733RgX_to_FlxqX(S, Tp, pp), Tp, pp));734}735else736r = FpXQXQ_minpoly(RgX_to_FpXQX(x, T, p), RgX_to_FpXQX(S, T, p), T, p);737return gerepileupto(av, FpXQX_to_mod(r, T, p));738}739740#define code(t1,t2) ((t1 << 6) | t2)741742static GEN743RgXQ_minpoly_fast(GEN x, GEN y)744{745GEN p, pol;746long pa;747long t = RgX_type2(x,y, &p,&pol,&pa);748switch(t)749{750case t_INTMOD: return RgXQ_minpoly_FpXQ(x, y, p);751case t_FFELT: return FFXQ_minpoly(x, y, pol);752case code(t_POLMOD, t_INTMOD):753return RgXQ_minpoly_FpXQXQ(x, y, pol, p);754default: return NULL;755}756}757758#undef code759760static GEN761easymin(GEN x, long v)762{763GEN G, R, dR;764long tx = typ(x);765if (tx == t_FFELT)766{767R = FpX_to_mod(FF_minpoly(x), FF_p_i(x));768setvarn(R,v); return R;769}770if (tx == t_POLMOD)771{772GEN a = gel(x,2), b = gel(x,1);773if (typ(a) != t_POL || varn(a) != varn(b))774{775if (varncmp(gvar(a), v) <= 0) pari_err_PRIORITY("minpoly", x, "<", v);776return deg1pol(gen_1, gneg_i(a), v);777}778R = RgXQ_minpoly_fast(a,b);779if (R) { setvarn(R, v); return R; }780if (!issquarefree(b))781{782R = RgXQ_minpoly_naive(a, b);783setvarn(R,v); return R;784}785}786R = easychar(x, v); if (!R) return NULL;787dR = RgX_deriv(R); if (!lgpol(dR)) return NULL;788G = RgX_normalize(RgX_gcd(R,dR));789return RgX_div(R,G);790}791GEN792minpoly(GEN x, long v)793{794pari_sp av = avma;795GEN P;796if (v < 0) v = 0;797P = easymin(x,v);798if (P) return gerepileupto(av,P);799/* typ(x) = t_MAT */800set_avma(av); return RgM_minpoly(x,v);801}802803/*******************************************************************/804/* */805/* HESSENBERG FORM */806/* */807/*******************************************************************/808static int809relative0(GEN a, GEN a0, long bit)810{ return gequal0(a) || (bit && gexpo(a)-gexpo(a0) < bit); }811/* assume x a nonempty square t_MAT */812static GEN813RgM_hess(GEN x0, long prec)814{815pari_sp av;816long lx = lg(x0), bit = prec? 8 - bit_accuracy(prec): 0, m, i, j;817GEN x;818819if (bit) x0 = RgM_shallowcopy(x0);820av = avma; x = RgM_shallowcopy(x0);821for (m=2; m<lx-1; m++)822{823GEN t = NULL;824for (i=m+1; i<lx; i++)825{826t = gcoeff(x,i,m-1);827if (!relative0(t, gcoeff(x0,i,m-1), bit)) break;828}829if (i == lx) continue;830for (j=m-1; j<lx; j++) swap(gcoeff(x,i,j), gcoeff(x,m,j));831swap(gel(x,i), gel(x,m));832if (bit)833{834for (j=m-1; j<lx; j++) swap(gcoeff(x0,i,j), gcoeff(x0,m,j));835swap(gel(x0,i), gel(x0,m));836}837t = ginv(t);838839for (i=m+1; i<lx; i++)840{841GEN c = gcoeff(x,i,m-1);842if (gequal0(c)) continue;843844c = gmul(c,t); gcoeff(x,i,m-1) = gen_0;845for (j=m; j<lx; j++)846gcoeff(x,i,j) = gsub(gcoeff(x,i,j), gmul(c,gcoeff(x,m,j)));847for (j=1; j<lx; j++)848gcoeff(x,j,m) = gadd(gcoeff(x,j,m), gmul(c,gcoeff(x,j,i)));849if (gc_needed(av,2))850{851if (DEBUGMEM>1) pari_warn(warnmem,"hess, m = %ld", m);852gerepileall(av,2, &x, &t);853}854}855}856return x;857}858859GEN860hess(GEN x)861{862pari_sp av = avma;863GEN p = NULL, T = NULL;864long prec, lx = lg(x);865if (typ(x) != t_MAT) pari_err_TYPE("hess",x);866if (lx == 1) return cgetg(1,t_MAT);867if (lgcols(x) != lx) pari_err_DIM("hess");868switch(RgM_type(x, &p, &T, &prec))869{870case t_REAL:871case t_COMPLEX: break;872default: prec = 0;873}874return gerepilecopy(av, RgM_hess(x,prec));875}876877GEN878Flm_hess(GEN x, ulong p)879{880long lx = lg(x), m, i, j;881if (lx == 1) return cgetg(1,t_MAT);882if (lgcols(x) != lx) pari_err_DIM("hess");883884x = Flm_copy(x);885for (m=2; m<lx-1; m++)886{887ulong t = 0;888for (i=m+1; i<lx; i++) { t = ucoeff(x,i,m-1); if (t) break; }889if (i == lx) continue;890for (j=m-1; j<lx; j++) lswap(ucoeff(x,i,j), ucoeff(x,m,j));891swap(gel(x,i), gel(x,m)); t = Fl_inv(t, p);892893for (i=m+1; i<lx; i++)894{895ulong c = ucoeff(x,i,m-1);896if (!c) continue;897898c = Fl_mul(c,t,p); ucoeff(x,i,m-1) = 0;899for (j=m; j<lx; j++)900ucoeff(x,i,j) = Fl_sub(ucoeff(x,i,j), Fl_mul(c,ucoeff(x,m,j), p), p);901for (j=1; j<lx; j++)902ucoeff(x,j,m) = Fl_add(ucoeff(x,j,m), Fl_mul(c,ucoeff(x,j,i), p), p);903}904}905return x;906}907GEN908FpM_hess(GEN x, GEN p)909{910pari_sp av = avma;911long lx = lg(x), m, i, j;912if (lx == 1) return cgetg(1,t_MAT);913if (lgcols(x) != lx) pari_err_DIM("hess");914if (lgefint(p) == 3)915{916ulong pp = p[2];917x = Flm_hess(ZM_to_Flm(x, pp), pp);918return gerepileupto(av, Flm_to_ZM(x));919}920x = RgM_shallowcopy(x);921for (m=2; m<lx-1; m++)922{923GEN t = NULL;924for (i=m+1; i<lx; i++) { t = gcoeff(x,i,m-1); if (signe(t)) break; }925if (i == lx) continue;926for (j=m-1; j<lx; j++) swap(gcoeff(x,i,j), gcoeff(x,m,j));927swap(gel(x,i), gel(x,m)); t = Fp_inv(t, p);928929for (i=m+1; i<lx; i++)930{931GEN c = gcoeff(x,i,m-1);932if (!signe(c)) continue;933934c = Fp_mul(c,t, p); gcoeff(x,i,m-1) = gen_0;935for (j=m; j<lx; j++)936gcoeff(x,i,j) = Fp_sub(gcoeff(x,i,j), Fp_mul(c,gcoeff(x,m,j),p), p);937for (j=1; j<lx; j++)938gcoeff(x,j,m) = Fp_add(gcoeff(x,j,m), Fp_mul(c,gcoeff(x,j,i),p), p);939if (gc_needed(av,2))940{941if (DEBUGMEM>1) pari_warn(warnmem,"hess, m = %ld", m);942gerepileall(av,2, &x, &t);943}944}945}946return gerepilecopy(av,x);947}948/* H in Hessenberg form */949static GEN950RgM_hess_charpoly(GEN H, long v)951{952long n = lg(H), r, i;953GEN y = cgetg(n+1, t_VEC);954gel(y,1) = pol_1(v);955for (r = 1; r < n; r++)956{957pari_sp av2 = avma;958GEN z, a = gen_1, b = pol_0(v);959for (i = r-1; i; i--)960{961a = gmul(a, gcoeff(H,i+1,i));962if (gequal0(a)) break;963b = RgX_add(b, RgX_Rg_mul(gel(y,i), gmul(a,gcoeff(H,i,r))));964}965z = RgX_sub(RgX_shift_shallow(gel(y,r), 1),966RgX_Rg_mul(gel(y,r), gcoeff(H,r,r)));967gel(y,r+1) = gerepileupto(av2, RgX_sub(z, b)); /* (X - H[r,r])y[r] - b */968}969return gel(y,n);970}971GEN972carhess(GEN x, long v)973{974pari_sp av;975GEN y;976if ((y = easychar(x,v))) return y;977av = avma; y = RgM_hess_charpoly(hess(x), v);978return fix_pol(av, y);979}980981/* Bound for sup norm of charpoly(M/dM), M integral: let B = |M|oo / |dM|,982* s = max_k binomial(n,k) (kB^2)^(k/2),983* return ceil(log2(s)) */984static long985charpoly_bound(GEN M, GEN dM, GEN N)986{987pari_sp av = avma;988GEN B = itor(N, LOWDEFAULTPREC);989GEN s = real_0(LOWDEFAULTPREC), bin, B2;990long n = lg(M)-1, k;991bin = gen_1;992if (dM) B = divri(B, dM);993B2 = sqrr(B);994for (k = n; k >= (n+1)>>1; k--)995{996GEN t = mulri(powruhalf(mulur(k, B2), k), bin);997if (abscmprr(t, s) > 0) s = t;998bin = diviuexact(muliu(bin, k), n-k+1);999}1000return gc_long(av, ceil(dbllog2(s)));1001}10021003static GEN1004QM_charpoly_ZX_slice(GEN A, GEN dM, GEN P, GEN *mod)1005{1006pari_sp av = avma;1007long i, n = lg(P)-1;1008GEN H, T;1009if (n == 1)1010{1011ulong p = uel(P,1), dp = dM ? umodiu(dM, p): 1;1012GEN Hp, a = ZM_to_Flm(A, p);1013Hp = Flm_charpoly_i(a, p);1014if (dp != 1) Hp = Flx_rescale(Hp, Fl_inv(dp, p), p);1015Hp = gerepileupto(av, Flx_to_ZX(Hp));1016*mod = utoipos(p); return Hp;1017}1018T = ZV_producttree(P);1019A = ZM_nv_mod_tree(A, P, T);1020H = cgetg(n+1, t_VEC);1021for(i=1; i <= n; i++)1022{1023ulong p = uel(P,i), dp = dM ? umodiu(dM, p): 1;1024gel(H,i) = Flm_charpoly(gel(A, i), p);1025if (dp != 1) gel(H,i) = Flx_rescale(gel(H,i), Fl_inv(dp, p), p);1026}1027H = nxV_chinese_center_tree(H, P, T, ZV_chinesetree(P,T));1028*mod = gmael(T, lg(T)-1, 1);1029gerepileall(av, 2, &H, mod); return H;1030}10311032GEN1033QM_charpoly_ZX_worker(GEN P, GEN M, GEN dM)1034{1035GEN V = cgetg(3, t_VEC);1036gel(V,1) = QM_charpoly_ZX_slice(M, equali1(dM) ? NULL:dM, P, &gel(V,2));1037return V;1038}10391040/* Assume M a square ZM, dM integer. Return charpoly(M / dM) in Z[X] */1041static GEN1042QM_charpoly_ZX_i(GEN M, GEN dM, long bound)1043{1044long n = lg(M)-1;1045forprime_t S;1046GEN worker = snm_closure(is_entry("_QM_charpoly_ZX_worker"),1047mkvec2(M, dM? dM: gen_1));1048if (!n) return pol_1(0);1049if (bound < 0)1050{1051GEN N = ZM_supnorm(M);1052if (signe(N) == 0) return monomial(gen_1, n, 0);1053bound = charpoly_bound(M, dM, N) + 1;1054}1055if (DEBUGLEVEL>5) err_printf("ZM_charpoly: bound 2^%ld\n", bound);1056init_modular_big(&S);1057return gen_crt("QM_charpoly_ZX", worker, &S, dM, bound, 0, NULL,1058nxV_chinese_center, FpX_center);1059}10601061GEN1062QM_charpoly_ZX_bound(GEN M, long bit)1063{1064pari_sp av = avma;1065GEN dM; M = Q_remove_denom(M, &dM);1066return gerepilecopy(av, QM_charpoly_ZX_i(M, dM, bit));1067}1068GEN1069QM_charpoly_ZX(GEN M)1070{1071pari_sp av = avma;1072GEN dM; M = Q_remove_denom(M, &dM);1073return gerepilecopy(av, QM_charpoly_ZX_i(M, dM, -1));1074}1075GEN1076ZM_charpoly(GEN M)1077{1078pari_sp av = avma;1079return gerepilecopy(av, QM_charpoly_ZX_i(M, NULL, -1));1080}10811082/*******************************************************************/1083/* */1084/* CHARACTERISTIC POLYNOMIAL (BERKOWITZ'S ALGORITHM) */1085/* */1086/*******************************************************************/1087GEN1088carberkowitz(GEN x, long v)1089{1090long lx, i, j, k, r;1091GEN V, S, C, Q;1092pari_sp av0, av;1093if ((V = easychar(x,v))) return V;1094lx = lg(x); av0 = avma;1095V = cgetg(lx+1, t_VEC);1096S = cgetg(lx+1, t_VEC);1097C = cgetg(lx+1, t_VEC);1098Q = cgetg(lx+1, t_VEC);1099av = avma;1100gel(C,1) = gen_m1;1101gel(V,1) = gen_m1;1102for (i=2;i<=lx; i++) gel(C,i) = gel(Q,i) = gel(S,i) = gel(V,i) = gen_0;1103gel(V,2) = gcoeff(x,1,1);1104for (r = 2; r < lx; r++)1105{1106pari_sp av2;1107GEN t;11081109for (i = 1; i < r; i++) gel(S,i) = gcoeff(x,i,r);1110gel(C,2) = gcoeff(x,r,r);1111for (i = 1; i < r-1; i++)1112{1113av2 = avma; t = gmul(gcoeff(x,r,1), gel(S,1));1114for (j = 2; j < r; j++) t = gadd(t, gmul(gcoeff(x,r,j), gel(S,j)));1115gel(C,i+2) = gerepileupto(av2, t);1116for (j = 1; j < r; j++)1117{1118av2 = avma; t = gmul(gcoeff(x,j,1), gel(S,1));1119for (k = 2; k < r; k++) t = gadd(t, gmul(gcoeff(x,j,k), gel(S,k)));1120gel(Q,j) = gerepileupto(av2, t);1121}1122for (j = 1; j < r; j++) gel(S,j) = gel(Q,j);1123}1124av2 = avma; t = gmul(gcoeff(x,r,1), gel(S,1));1125for (j = 2; j < r; j++) t = gadd(t, gmul(gcoeff(x,r,j), gel(S,j)));1126gel(C,r+1) = gerepileupto(av2, t);1127if (gc_needed(av0,1))1128{1129if (DEBUGMEM>1) pari_warn(warnmem,"carberkowitz");1130gerepileall(av, 2, &C, &V);1131}1132for (i = 1; i <= r+1; i++)1133{1134av2 = avma; t = gmul(gel(C,i), gel(V,1));1135for (j = 2; j <= minss(r,i); j++)1136t = gadd(t, gmul(gel(C,i+1-j), gel(V,j)));1137gel(Q,i) = gerepileupto(av2, t);1138}1139for (i = 1; i <= r+1; i++) gel(V,i) = gel(Q,i);1140}1141V = RgV_to_RgX(vecreverse(V), v); /* not gtopoly: fail if v > gvar(V) */1142V = odd(lx)? gcopy(V): RgX_neg(V);1143return fix_pol(av0, V);1144}11451146/*******************************************************************/1147/* */1148/* NORMS */1149/* */1150/*******************************************************************/1151GEN1152gnorm(GEN x)1153{1154pari_sp av;1155long lx, i;1156GEN y;11571158switch(typ(x))1159{1160case t_INT: return sqri(x);1161case t_REAL: return sqrr(x);1162case t_FRAC: return sqrfrac(x);1163case t_COMPLEX: av = avma; return gerepileupto(av, cxnorm(x));1164case t_QUAD: av = avma; return gerepileupto(av, quadnorm(x));11651166case t_POL: case t_SER: case t_RFRAC: av = avma;1167return gerepileupto(av, greal(gmul(conj_i(x),x)));11681169case t_FFELT:1170y = cgetg(3, t_INTMOD);1171gel(y,1) = FF_p(x);1172gel(y,2) = FF_norm(x); return y;11731174case t_POLMOD:1175{1176GEN T = gel(x,1), a = gel(x,2);1177if (typ(a) != t_POL || varn(a) != varn(T)) return gpowgs(a, degpol(T));1178return RgXQ_norm(a, T);1179}1180case t_VEC: case t_COL: case t_MAT:1181y = cgetg_copy(x, &lx);1182for (i=1; i<lx; i++) gel(y,i) = gnorm(gel(x,i));1183return y;1184}1185pari_err_TYPE("gnorm",x);1186return NULL; /* LCOV_EXCL_LINE */1187}11881189/* return |q|^2, complex modulus */1190static GEN1191cxquadnorm(GEN q, long prec)1192{1193GEN X = gel(q,1), c = gel(X,2); /* (1-D)/4, -D/4 */1194if (signe(c) > 0) return quadnorm(q); /* imaginary */1195if (!prec) pari_err_TYPE("gnorml2", q);1196return sqrr(quadtofp(q, prec));1197}11981199static GEN1200gnorml2_i(GEN x, long prec)1201{1202pari_sp av;1203long i, lx;1204GEN s;12051206switch(typ(x))1207{1208case t_INT: return sqri(x);1209case t_REAL: return sqrr(x);1210case t_FRAC: return sqrfrac(x);1211case t_COMPLEX: av = avma; return gerepileupto(av, cxnorm(x));1212case t_QUAD: av = avma; return gerepileupto(av, cxquadnorm(x,prec));12131214case t_POL: lx = lg(x)-1; x++; break;12151216case t_VEC:1217case t_COL:1218case t_MAT: lx = lg(x); break;12191220default: pari_err_TYPE("gnorml2",x);1221return NULL; /* LCOV_EXCL_LINE */1222}1223if (lx == 1) return gen_0;1224av = avma;1225s = gnorml2(gel(x,1));1226for (i=2; i<lx; i++)1227{1228s = gadd(s, gnorml2(gel(x,i)));1229if (gc_needed(av,1))1230{1231if (DEBUGMEM>1) pari_warn(warnmem,"gnorml2");1232s = gerepileupto(av, s);1233}1234}1235return gerepileupto(av,s);1236}1237GEN1238gnorml2(GEN x) { return gnorml2_i(x, 0); }12391240static GEN pnormlp(GEN,GEN,long);1241static GEN1242pnormlpvec(long i0, GEN x, GEN p, long prec)1243{1244pari_sp av = avma;1245long i, lx = lg(x);1246GEN s = gen_0;1247for (i=i0; i<lx; i++)1248{1249s = gadd(s, pnormlp(gel(x,i),p,prec));1250if (gc_needed(av,1))1251{1252if (DEBUGMEM>1) pari_warn(warnmem,"gnormlp, i = %ld", i);1253s = gerepileupto(av, s);1254}1255}1256return s;1257}1258/* (||x||_p)^p */1259static GEN1260pnormlp(GEN x, GEN p, long prec)1261{1262switch(typ(x))1263{1264case t_INT: case t_REAL: x = mpabs(x); break;1265case t_FRAC: x = absfrac(x); break;1266case t_COMPLEX: case t_QUAD: x = gabs(x,prec); break;1267case t_POL: return pnormlpvec(2, x, p, prec);1268case t_VEC: case t_COL: case t_MAT: return pnormlpvec(1, x, p, prec);1269default: pari_err_TYPE("gnormlp",x);1270}1271return gpow(x, p, prec);1272}12731274GEN1275gnormlp(GEN x, GEN p, long prec)1276{1277pari_sp av = avma;1278if (!p || (typ(p) == t_INFINITY && inf_get_sign(p) > 0))1279return gsupnorm(x, prec);1280if (gsigne(p) <= 0) pari_err_DOMAIN("normlp", "p", "<=", gen_0, p);1281if (is_scalar_t(typ(x))) return gabs(x, prec);1282if (typ(p) == t_INT)1283{1284ulong pp = itou_or_0(p);1285switch(pp)1286{1287case 1: return gnorml1(x, prec);1288case 2: x = gnorml2_i(x, prec); break;1289default: x = pnormlp(x, p, prec); break;1290}1291if (pp && typ(x) == t_INT && Z_ispowerall(x, pp, &x))1292return gerepileuptoleaf(av, x);1293if (pp == 2) return gerepileupto(av, gsqrt(x, prec));1294}1295else1296x = pnormlp(x, p, prec);1297x = gpow(x, ginv(p), prec);1298return gerepileupto(av, x);1299}13001301GEN1302gnorml1(GEN x,long prec)1303{1304pari_sp av = avma;1305long lx,i;1306GEN s;1307switch(typ(x))1308{1309case t_INT: case t_REAL: return mpabs(x);1310case t_FRAC: return absfrac(x);13111312case t_COMPLEX: case t_QUAD:1313return gabs(x,prec);13141315case t_POL:1316lx = lg(x); s = gen_0;1317for (i=2; i<lx; i++) s = gadd(s, gnorml1(gel(x,i),prec));1318break;13191320case t_VEC: case t_COL: case t_MAT:1321lx = lg(x); s = gen_0;1322for (i=1; i<lx; i++) s = gadd(s, gnorml1(gel(x,i),prec));1323break;13241325default: pari_err_TYPE("gnorml1",x);1326return NULL; /* LCOV_EXCL_LINE */1327}1328return gerepileupto(av, s);1329}1330/* As gnorml1, except for t_QUAD and t_COMPLEX: |x + wy| := |x| + |y|1331* Still a norm of R-vector spaces, and can be cheaply computed without1332* square roots */1333GEN1334gnorml1_fake(GEN x)1335{1336pari_sp av = avma;1337long lx, i;1338GEN s;1339switch(typ(x))1340{1341case t_INT: case t_REAL: return mpabs(x);1342case t_FRAC: return absfrac(x);13431344case t_COMPLEX:1345s = gadd(gnorml1_fake(gel(x,1)), gnorml1_fake(gel(x,2)));1346break;1347case t_QUAD:1348s = gadd(gnorml1_fake(gel(x,2)), gnorml1_fake(gel(x,3)));1349break;13501351case t_POL:1352lx = lg(x); s = gen_0;1353for (i=2; i<lx; i++) s = gadd(s, gnorml1_fake(gel(x,i)));1354break;13551356case t_VEC: case t_COL: case t_MAT:1357lx = lg(x); s = gen_0;1358for (i=1; i<lx; i++) s = gadd(s, gnorml1_fake(gel(x,i)));1359break;13601361default: pari_err_TYPE("gnorml1_fake",x);1362return NULL; /* LCOV_EXCL_LINE */1363}1364return gerepileupto(av, s);1365}13661367static void1368store(GEN z, GEN *m) { if (!*m || gcmp(z, *m) > 0) *m = z; }1369/* compare |x| to *m or |x|^2 to *msq, whichever is easiest, and update1370* the pointed value if x is larger */1371void1372gsupnorm_aux(GEN x, GEN *m, GEN *msq, long prec)1373{1374long i, lx;1375GEN z;1376switch(typ(x))1377{1378case t_COMPLEX: z = cxnorm(x); store(z, msq); return;1379case t_QUAD: z = cxquadnorm(x,prec); store(z, msq); return;1380case t_INT: case t_REAL: z = mpabs(x); store(z,m); return;1381case t_FRAC: z = absfrac(x); store(z,m); return;13821383case t_POL: lx = lg(x)-1; x++; break;13841385case t_VEC:1386case t_COL:1387case t_MAT: lx = lg(x); break;13881389default: pari_err_TYPE("gsupnorm",x);1390return; /* LCOV_EXCL_LINE */1391}1392for (i=1; i<lx; i++) gsupnorm_aux(gel(x,i), m, msq, prec);1393}1394GEN1395gsupnorm(GEN x, long prec)1396{1397GEN m = NULL, msq = NULL;1398pari_sp av = avma;1399gsupnorm_aux(x, &m, &msq, prec);1400/* now set m = max (m, sqrt(msq)) */1401if (msq) {1402msq = gsqrt(msq, prec);1403if (!m || gcmp(m, msq) < 0) m = msq;1404} else if (!m) m = gen_0;1405return gerepilecopy(av, m);1406}14071408/*******************************************************************/1409/* */1410/* TRACES */1411/* */1412/*******************************************************************/1413GEN1414matcompanion(GEN x)1415{1416long n = degpol(x), j;1417GEN y, c;14181419if (typ(x)!=t_POL) pari_err_TYPE("matcompanion",x);1420if (!signe(x)) pari_err_DOMAIN("matcompanion","polynomial","=",gen_0,x);1421if (n == 0) return cgetg(1, t_MAT);14221423y = cgetg(n+1,t_MAT);1424for (j=1; j < n; j++) gel(y,j) = col_ei(n, j+1);1425c = cgetg(n+1,t_COL); gel(y,n) = c;1426if (gequal1(gel(x, n+2)))1427for (j=1; j<=n; j++) gel(c,j) = gneg(gel(x,j+1));1428else1429{ /* not monic. Hardly ever used */1430pari_sp av = avma;1431GEN d = gclone(gneg(gel(x,n+2)));1432set_avma(av);1433for (j=1; j<=n; j++) gel(c,j) = gdiv(gel(x,j+1), d);1434gunclone(d);1435}1436return y;1437}14381439GEN1440gtrace(GEN x)1441{1442pari_sp av;1443long i, lx, tx = typ(x);1444GEN y, z;14451446switch(tx)1447{1448case t_INT: case t_REAL: case t_FRAC:1449return gmul2n(x,1);14501451case t_COMPLEX:1452return gmul2n(gel(x,1),1);14531454case t_QUAD:1455y = gel(x,1);1456if (!gequal0(gel(y,3)))1457{ /* assume quad. polynomial is either x^2 + d or x^2 - x + d */1458av = avma;1459return gerepileupto(av, gadd(gel(x,3), gmul2n(gel(x,2),1)));1460}1461return gmul2n(gel(x,2),1);14621463case t_POL:1464y = cgetg_copy(x, &lx); y[1] = x[1];1465for (i=2; i<lx; i++) gel(y,i) = gtrace(gel(x,i));1466return normalizepol_lg(y, lx);14671468case t_SER:1469if (ser_isexactzero(x)) return gcopy(x);1470y = cgetg_copy(x, &lx); y[1] = x[1];1471for (i=2; i<lx; i++) gel(y,i) = gtrace(gel(x,i));1472return normalize(y);14731474case t_POLMOD:1475y = gel(x,1); z = gel(x,2);1476if (typ(z) != t_POL || varn(y) != varn(z)) return gmulsg(degpol(y), z);1477av = avma;1478return gerepileupto(av, RgXQ_trace(z, y));14791480case t_FFELT:1481y=cgetg(3, t_INTMOD);1482gel(y,1) = FF_p(x);1483gel(y,2) = FF_trace(x);1484return y;14851486case t_RFRAC:1487av = avma; return gerepileupto(av, gadd(x, conj_i(x)));14881489case t_VEC: case t_COL:1490y = cgetg_copy(x, &lx);1491for (i=1; i<lx; i++) gel(y,i) = gtrace(gel(x,i));1492return y;14931494case t_MAT:1495lx = lg(x); if (lx == 1) return gen_0;1496/*now lx >= 2*/1497if (lx != lgcols(x)) pari_err_DIM("gtrace");1498av = avma; return gerepileupto(av, mattrace(x));1499}1500pari_err_TYPE("gtrace",x);1501return NULL; /* LCOV_EXCL_LINE */1502}15031504/* Cholesky decomposition for positive definite matrix a1505* [GTM138, Algo 2.7.6, matrix Q]1506* If a is not positive definite return NULL. */1507GEN1508qfgaussred_positive(GEN a)1509{1510pari_sp av = avma;1511GEN b;1512long i,j,k, n = lg(a);15131514if (typ(a)!=t_MAT) pari_err_TYPE("qfgaussred_positive",a);1515if (n == 1) return cgetg(1, t_MAT);1516if (lgcols(a)!=n) pari_err_DIM("qfgaussred_positive");1517b = cgetg(n,t_MAT);1518for (j=1; j<n; j++)1519{1520GEN p1=cgetg(n,t_COL), p2=gel(a,j);15211522gel(b,j) = p1;1523for (i=1; i<=j; i++) gel(p1,i) = gel(p2,i);1524for ( ; i<n ; i++) gel(p1,i) = gen_0;1525}1526for (k=1; k<n; k++)1527{1528GEN bk, p = gcoeff(b,k,k), invp;1529if (gsigne(p)<=0) return gc_NULL(av); /* not positive definite */1530invp = ginv(p);1531bk = row(b, k);1532for (i=k+1; i<n; i++) gcoeff(b,k,i) = gmul(gel(bk,i), invp);1533for (i=k+1; i<n; i++)1534{1535GEN c = gel(bk, i);1536for (j=i; j<n; j++)1537gcoeff(b,i,j) = gsub(gcoeff(b,i,j), gmul(c,gcoeff(b,k,j)));1538}1539if (gc_needed(av,1))1540{1541if (DEBUGMEM>1) pari_warn(warnmem,"qfgaussred_positive");1542b=gerepilecopy(av,b);1543}1544}1545return gerepilecopy(av,b);1546}15471548/* Maximal pivot strategy: x is a suitable pivot if it is non zero and either1549* - an exact type, or1550* - it is maximal among remaining nonzero (t_REAL) pivots */1551static int1552suitable(GEN x, long k, GEN *pp, long *pi)1553{1554long t = typ(x);1555switch(t)1556{1557case t_INT: return signe(x) != 0;1558case t_FRAC: return 1;1559case t_REAL: {1560GEN p = *pp;1561if (signe(x) && (!p || abscmprr(p, x) < 0)) { *pp = x; *pi = k; }1562return 0;1563}1564default: return !gequal0(x);1565}1566}15671568/* Gauss reduction (arbitrary symetric matrix, only the part above the1569* diagonal is considered). If signature is nonzero, return only the1570* signature, in which case gsigne() should be defined for elements of a. */1571static GEN1572gaussred(GEN a, long signature)1573{1574GEN r, ak, al;1575pari_sp av = avma, av1;1576long n = lg(a), i, j, k, l, sp, sn, t;15771578if (typ(a) != t_MAT) pari_err_TYPE("gaussred",a);1579if (n == 1) return signature? mkvec2(gen_0, gen_0): cgetg(1, t_MAT);1580if (lgcols(a) != n) pari_err_DIM("gaussred");1581n--;1582r = const_vecsmall(n, 1); av1= avma;1583a = RgM_shallowcopy(a);1584t = n; sp = sn = 0;1585while (t)1586{1587long pind = 0;1588GEN invp, p = NULL;1589k=1; while (k<=n && (!r[k] || !suitable(gcoeff(a,k,k), k, &p, &pind))) k++;1590if (k > n && p) k = pind;1591if (k <= n)1592{1593p = gcoeff(a,k,k); invp = ginv(p); /* != 0 */1594if (signature) { /* skip if (!signature): gsigne may fail ! */1595if (gsigne(p) > 0) sp++; else sn++;1596}1597r[k] = 0; t--;1598ak = row(a, k);1599for (i=1; i<=n; i++)1600gcoeff(a,k,i) = r[i]? gmul(gcoeff(a,k,i), invp): gen_0;16011602for (i=1; i<=n; i++) if (r[i])1603{1604GEN c = gel(ak,i); /* - p * a[k,i] */1605if (gequal0(c)) continue;1606for (j=1; j<=n; j++) if (r[j])1607gcoeff(a,i,j) = gsub(gcoeff(a,i,j), gmul(c,gcoeff(a,k,j)));1608}1609gcoeff(a,k,k) = p;1610if (gc_needed(av1,1))1611{1612if(DEBUGMEM>1) pari_warn(warnmem,"gaussred (t = %ld)", t);1613a = gerepilecopy(av1, a);1614}1615}1616else1617{ /* all remaining diagonal coeffs are currently 0 */1618for (k=1; k<=n; k++) if (r[k])1619{1620l=k+1; while (l<=n && (!r[l] || !suitable(gcoeff(a,k,l), l, &p, &pind))) l++;1621if (l > n && p) l = pind;1622if (l > n) continue;16231624p = gcoeff(a,k,l); invp = ginv(p);1625sp++; sn++;1626r[k] = r[l] = 0; t -= 2;1627ak = row(a, k);1628al = row(a, l);1629for (i=1; i<=n; i++) if (r[i])1630{1631gcoeff(a,k,i) = gmul(gcoeff(a,k,i), invp);1632gcoeff(a,l,i) = gmul(gcoeff(a,l,i), invp);1633} else {1634gcoeff(a,k,i) = gen_0;1635gcoeff(a,l,i) = gen_0;1636}16371638for (i=1; i<=n; i++) if (r[i])1639{ /* c = a[k,i] * p, d = a[l,i] * p; */1640GEN c = gel(ak,i), d = gel(al,i);1641for (j=1; j<=n; j++) if (r[j])1642gcoeff(a,i,j) = gsub(gcoeff(a,i,j),1643gadd(gmul(gcoeff(a,l,j), c),1644gmul(gcoeff(a,k,j), d)));1645}1646for (i=1; i<=n; i++) if (r[i])1647{1648GEN c = gcoeff(a,k,i), d = gcoeff(a,l,i);1649gcoeff(a,k,i) = gadd(c, d);1650gcoeff(a,l,i) = gsub(c, d);1651}1652gcoeff(a,k,l) = gen_1;1653gcoeff(a,l,k) = gen_m1;1654gcoeff(a,k,k) = gmul2n(p,-1);1655gcoeff(a,l,l) = gneg(gcoeff(a,k,k));1656if (gc_needed(av1,1))1657{1658if(DEBUGMEM>1) pari_warn(warnmem,"gaussred");1659a = gerepilecopy(av1, a);1660}1661break;1662}1663if (k > n) break;1664}1665}1666if (!signature) return gerepilecopy(av, a);1667set_avma(av); return mkvec2s(sp, sn);1668}16691670GEN1671qfgaussred(GEN a) { return gaussred(a,0); }16721673GEN1674qfsign(GEN a) { return gaussred(a,1); }16751676/* x -= s(y+u*x) */1677/* y += s(x-u*y), simultaneously */1678static void1679rot(GEN x, GEN y, GEN s, GEN u) {1680GEN x1 = subrr(x, mulrr(s,addrr(y,mulrr(u,x))));1681GEN y1 = addrr(y, mulrr(s,subrr(x,mulrr(u,y))));1682affrr(x1,x);1683affrr(y1,y);1684}16851686/* Diagonalization of a REAL symetric matrix. Return a vector [L, r]:1687* L = vector of eigenvalues1688* r = matrix of eigenvectors */1689GEN1690jacobi(GEN a, long prec)1691{1692pari_sp av;1693long de, e, e1, e2, i, j, p, q, l = lg(a);1694GEN c, ja, L, r, L2, r2, unr;16951696if (typ(a) != t_MAT) pari_err_TYPE("jacobi",a);1697ja = cgetg(3,t_VEC);1698L = cgetg(l,t_COL); gel(ja,1) = L;1699r = cgetg(l,t_MAT); gel(ja,2) = r;1700if (l == 1) return ja;1701if (lgcols(a) != l) pari_err_DIM("jacobi");17021703e1 = HIGHEXPOBIT-1;1704for (j=1; j<l; j++)1705{1706GEN z = gtofp(gcoeff(a,j,j), prec);1707gel(L,j) = z;1708e = expo(z); if (e < e1) e1 = e;1709}1710for (j=1; j<l; j++)1711{1712gel(r,j) = cgetg(l,t_COL);1713for (i=1; i<l; i++) gcoeff(r,i,j) = utor(i==j? 1: 0, prec);1714}1715av = avma;17161717e2 = -(long)HIGHEXPOBIT; p = q = 1;1718c = cgetg(l,t_MAT);1719for (j=1; j<l; j++)1720{1721gel(c,j) = cgetg(j,t_COL);1722for (i=1; i<j; i++)1723{1724GEN z = gtofp(gcoeff(a,i,j), prec);1725gcoeff(c,i,j) = z;1726if (!signe(z)) continue;1727e = expo(z); if (e > e2) { e2 = e; p = i; q = j; }1728}1729}1730a = c; unr = real_1(prec);1731de = prec2nbits(prec);17321733/* e1 = min expo(a[i,i])1734* e2 = max expo(a[i,j]), i != j */1735while (e1-e2 < de)1736{1737pari_sp av2 = avma;1738GEN x, y, t, c, s, u;1739/* compute attached rotation in the plane formed by basis vectors number1740* p and q */1741x = subrr(gel(L,q),gel(L,p));1742if (signe(x))1743{1744x = divrr(x, shiftr(gcoeff(a,p,q),1));1745y = sqrtr(addrr(unr, sqrr(x)));1746t = invr((signe(x)>0)? addrr(x,y): subrr(x,y));1747}1748else1749y = t = unr;1750c = sqrtr(addrr(unr,sqrr(t)));1751s = divrr(t,c);1752u = divrr(t,addrr(unr,c));17531754/* compute successive transforms of a and the matrix of accumulated1755* rotations (r) */1756for (i=1; i<p; i++) rot(gcoeff(a,i,p), gcoeff(a,i,q), s,u);1757for (i=p+1; i<q; i++) rot(gcoeff(a,p,i), gcoeff(a,i,q), s,u);1758for (i=q+1; i<l; i++) rot(gcoeff(a,p,i), gcoeff(a,q,i), s,u);1759y = gcoeff(a,p,q);1760t = mulrr(t, y); shiftr_inplace(y, -de - 1);1761x = gel(L,p); subrrz(x,t, x);1762y = gel(L,q); addrrz(y,t, y);1763for (i=1; i<l; i++) rot(gcoeff(r,i,p), gcoeff(r,i,q), s,u);17641765e2 = -(long)HIGHEXPOBIT; p = q = 1;1766for (j=1; j<l; j++)1767{1768for (i=1; i<j; i++)1769{1770GEN z = gcoeff(a,i,j);1771if (!signe(z)) continue;1772e = expo(z); if (e > e2) { e2=e; p=i; q=j; }1773}1774for (i=j+1; i<l; i++)1775{1776GEN z = gcoeff(a,j,i);1777if (!signe(z)) continue;1778e = expo(z); if (e > e2) { e2=e; p=j; q=i; }1779}1780}1781set_avma(av2);1782}1783/* sort eigenvalues from smallest to largest */1784c = indexsort(L);1785r2 = vecpermute(r, c); for (i=1; i<l; i++) gel(r,i) = gel(r2,i);1786L2 = vecpermute(L, c); for (i=1; i<l; i++) gel(L,i) = gel(L2,i);1787set_avma(av); return ja;1788}17891790/*************************************************************************/1791/** **/1792/** Q-vector space -> Z-modules **/1793/** **/1794/*************************************************************************/17951796GEN1797matrixqz0(GEN x,GEN p)1798{1799if (typ(x) != t_MAT) pari_err_TYPE("matrixqz",x);1800if (!p) return QM_minors_coprime(x,NULL);1801if (typ(p) != t_INT) pari_err_TYPE("matrixqz",p);1802if (signe(p)>=0) return QM_minors_coprime(x,p);1803if (!RgM_is_QM(x)) pari_err_TYPE("matrixqz", x);1804if (absequaliu(p,1)) return QM_ImZ(x); /* p = -1 */1805if (absequaliu(p,2)) return QM_ImQ(x); /* p = -2 */1806pari_err_FLAG("QM_minors_coprime");1807return NULL; /* LCOV_EXCL_LINE */1808}18091810GEN1811QM_minors_coprime(GEN x, GEN D)1812{1813pari_sp av = avma, av1;1814long i, j, m, n, lP;1815GEN P, y;18161817n = lg(x)-1; if (!n) return gcopy(x);1818m = nbrows(x);1819if (n > m) pari_err_DOMAIN("QM_minors_coprime","n",">",strtoGENstr("m"),x);1820y = x; x = cgetg(n+1,t_MAT);1821for (j=1; j<=n; j++)1822{1823gel(x,j) = Q_primpart(gel(y,j));1824RgV_check_ZV(gel(x,j), "QM_minors_coprime");1825}1826/* x now a ZM */1827if (n==m)1828{1829if (gequal0(ZM_det(x)))1830pari_err_DOMAIN("QM_minors_coprime", "rank(A)", "<",stoi(n),x);1831set_avma(av); return matid(n);1832}1833/* m > n */1834if (!D || gequal0(D))1835{1836pari_sp av2 = avma;1837D = ZM_detmult(shallowtrans(x));1838if (is_pm1(D)) { set_avma(av2); return ZM_copy(x); }1839}1840P = gel(Z_factor(D), 1); lP = lg(P);1841av1 = avma;1842for (i=1; i < lP; i++)1843{1844GEN p = gel(P,i), pov2 = shifti(p, -1);1845for(;;)1846{1847GEN N, M = FpM_ker(x, p);1848long lM = lg(M);1849if (lM==1) break;18501851FpM_center_inplace(M, p, pov2);1852N = ZM_Z_divexact(ZM_mul(x,M), p);1853for (j=1; j<lM; j++)1854{1855long k = n; while (!signe(gcoeff(M,k,j))) k--;1856gel(x,k) = gel(N,j);1857}1858if (gc_needed(av1,1))1859{1860if (DEBUGMEM>1) pari_warn(warnmem,"QM_minors_coprime, p = %Ps", p);1861x = gerepilecopy(av1, x); pov2 = shifti(p, -1);1862}1863}1864}1865return gerepilecopy(av, x);1866}18671868static GEN1869QM_ImZ_all_i(GEN A, GEN *U, long remove, long hnf, long linindep)1870{1871GEN V = NULL, D;1872A = Q_remove_denom(A,&D);1873if (D)1874{1875long l, lA;1876V = matkermod(A,D,NULL);1877l = lg(V); lA = lg(A);1878if (l == 1) V = scalarmat_shallow(D, lA-1);1879else1880{1881if (l < lA) V = hnfmodid(V,D);1882A = ZM_Z_divexact(ZM_mul(A, V), D);1883}1884}1885if (!linindep && ZM_rank(A)==lg(A)-1) linindep = 1;1886if (hnf || !linindep) A = ZM_hnflll(A, U, remove);1887if (U && V)1888{1889if (hnf) *U = ZM_mul(V,*U);1890else *U = V;1891}1892return A;1893}1894GEN1895QM_ImZ_all(GEN x, GEN *U, long remove, long hnf)1896{1897pari_sp av = avma;1898x = QM_ImZ_all_i(x, U, remove, hnf, 0);1899if (U) gerepileall(av, 2, &x, &U); else x = gerepilecopy(av, x);1900return x;1901}1902GEN1903QM_ImZ_hnfall(GEN x, GEN *U, long remove) { return QM_ImZ_all(x, U, remove, 1); }1904GEN1905QM_ImZ_hnf(GEN x) { return QM_ImZ_hnfall(x, NULL, 1); }1906GEN1907QM_ImZ(GEN x) { return QM_ImZ_all(x, NULL, 1, 0); }19081909GEN1910QM_ImQ_all(GEN x, GEN *U, long remove, long hnf)1911{1912pari_sp av = avma;1913long i, n = lg(x), m;1914GEN ir, V, D, c, K = NULL;19151916if (U) *U = matid(n-1);1917if (n==1) return gcopy(x);1918m = lg(gel(x,1));19191920x = RgM_shallowcopy(x);1921for (i=1; i<n; i++)1922{1923gel(x,i) = Q_primitive_part(gel(x,i), &c);1924if (U && c && signe(c)) gcoeff(*U,i,i) = ginv(c);1925}19261927ir = ZM_indexrank(x);1928if (U)1929{1930*U = vecpermute(*U, gel(ir,2));1931if (remove < 2) K = ZM_ker(x);1932}1933x = vecpermute(x, gel(ir,2));19341935D = absi(ZM_det(rowpermute(x,gel(ir,1))));1936x = RgM_Rg_div(x, D);1937x = QM_ImZ_all_i(x, U? &V: NULL, remove, hnf, 1);19381939if (U)1940{1941*U = RgM_Rg_div(RgM_mul(*U,V),D);1942if (remove < 2) *U = shallowconcat(K,*U);1943if (!remove) x = shallowconcat(zeromatcopy(m-1,lg(K)-1), x);1944gerepileall(av, 2, &x, U);1945}1946else x = gerepilecopy(av,x);1947return x;1948}1949GEN1950QM_ImQ_hnfall(GEN x, GEN *U, long remove) { return QM_ImQ_all(x, U, remove, 1); }1951GEN1952QM_ImQ_hnf(GEN x) { return QM_ImQ_hnfall(x, NULL, 1); }1953GEN1954QM_ImQ(GEN x) { return QM_ImQ_all(x, NULL, 1, 0); }19551956GEN1957intersect(GEN x, GEN y)1958{1959long j, lx = lg(x);1960pari_sp av;1961GEN z;19621963if (typ(x)!=t_MAT) pari_err_TYPE("intersect",x);1964if (typ(y)!=t_MAT) pari_err_TYPE("intersect",y);1965if (lx==1 || lg(y)==1) return cgetg(1,t_MAT);19661967av = avma; z = ker(shallowconcat(x,y));1968for (j=lg(z)-1; j; j--) setlg(z[j], lx);1969return gerepileupto(av, image(RgM_mul(x,z)));1970}197119721973