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. */13#include "pari.h"14#include "paripriv.h"1516#define DEBUGLEVEL DEBUGLEVEL_mathnf1718#define dbg_printf(lvl) if (DEBUGLEVEL >= (lvl) + 3) err_printf1920/********************************************************************/21/** **/22/** BLACK BOX HERMITE RINGS AND HOWELL NORMAL FORM **/23/** contributed by Aurel Page (2017) **/24/** **/25/********************************************************************/2627/*28bb_hermite R:29- add(a,b): a+b30- neg(a): -a31- mul(a,b): a*b32- extgcd(a,b,&small): [d,U] with d in R and U in GL_2(R) such that [0;d] = [a;b]*U.33set small==1 to assert that U is a 'small' operation (no red needed).34- rann(a): b in R such that b*R = {x in R | a*x==0}35- lquo(a,b,&r): q in R such that r=a-b*q is a canonical representative36of the image of a in R/b*R. The canonical lift of 0 must be 0.37- unit(a): u unit in R^* such that a*u is a canonical generator of the ideal a*R38- equal0(a): a==0?39- equal1(a): a==1?40- s(n): image of the small integer n in R41- red(a): unique representative of a as an element of R4243op encoding of elementary operations:44- t_VECSMALL: the corresponding permutation (vecpermute)45- [Vecsmall([i,j])]: the transposition Ci <-> Cj46- [Vecsmall([i]),u], u in R^*: Ci <- Ci*u47- [Vecsmall([i,j]),a], a in R: Ci <- Ci + Cj*a48- [Vecsmall([i,j,0]),U], U in GL_2(R): (Ci|Cj) <- (Ci|Cj)*U49*/5051struct bb_hermite52{53GEN (*add)(void*, GEN, GEN);54GEN (*neg)(void*, GEN);55GEN (*mul)(void*, GEN, GEN);56GEN (*extgcd)(void*, GEN, GEN, int*);57GEN (*rann)(void*, GEN);58GEN (*lquo)(void*, GEN, GEN, GEN*);59GEN (*unit)(void*, GEN);60int (*equal0)(GEN);61int (*equal1)(GEN);62GEN (*s)(void*, long);63GEN (*red)(void*, GEN);64};6566static GEN67_Fp_add(void *data, GEN x, GEN y) { (void) data; return addii(x,y); }6869static GEN70_Fp_neg(void *data, GEN x) { (void) data; return negi(x); }7172static GEN73_Fp_mul(void *data, GEN x, GEN y) { (void) data; return mulii(x,y); }7475static GEN76_Fp_rann(void *data, GEN x)77{78GEN d, N = (GEN)data;79if (!signe(x)) return gen_1;80d = gcdii(x,N);81return modii(diviiexact(N,d),N);82}8384static GEN85_Fp_lquo(void *data, GEN x, GEN y, GEN* r) { (void) data; return truedvmdii(x,y,r); }8687/* D=MN, p|M => !p|a, p|N => p|a, return M */88static GEN89Z_split(GEN D, GEN a)90{91long i, n;92GEN N;93n = expi(D);94n = n<2 ? 1 : expu(n)+1;95for (i=1;i<=n;i++)96a = Fp_sqr(a,D);97N = gcdii(a,D);98return diviiexact(D,N);99}100101/* c s.t. gcd(a+cb,N) = gcd(a,b,N) without factoring */102static GEN103Z_stab(GEN a, GEN b, GEN N)104{105GEN g, a2, N2;106g = gcdii(a,b);107g = gcdii(g,N);108N2 = diviiexact(N,g);109a2 = diviiexact(a,g);110return Z_split(N2,a2);111}112113static GEN114_Fp_unit(void *data, GEN x)115{116GEN g,s,v,d,N=(GEN)data,N2;117long i;118if (!signe(x)) return NULL;119g = bezout(x,N,&s,&v);120if (equali1(g) || equali1(gcdii(s,N))) return mkvec2(g,s);121N2 = diviiexact(N,g);122for (i=0; i<5; i++)123{124s = addii(s,N2);125if (equali1(gcdii(s,N))) return mkvec2(g,s);126}127d = Z_stab(s,N2,N);128d = mulii(d,N2);129v = Fp_add(s,d,N);130if (equali1(v)) return NULL;131return mkvec2(g,v);132}133134static GEN135_Fp_extgcd(void *data, GEN x, GEN y, int* smallop)136{137GEN d,u,v,m;138if (equali1(y))139{140*smallop = 1;141return mkvec2(y,mkmat2(142mkcol2(gen_1,Fp_neg(x,(GEN)data)),143mkcol2(gen_0,gen_1)));144}145*smallop = 0;146d = bezout(x,y,&u,&v);147if (!signe(d)) return mkvec2(d,matid(2));148m = cgetg(3,t_MAT);149m = mkmat2(150mkcol2(diviiexact(y,d),negi(diviiexact(x,d))),151mkcol2(u,v));152return mkvec2(d,m);153}154155static int156_Fp_equal0(GEN x) { return !signe(x); }157158static int159_Fp_equal1(GEN x) { return equali1(x); }160161static GEN162_Fp_s(void *data, long x)163{164if (!x) return gen_0;165if (x==1) return gen_1;166return modsi(x,(GEN)data);167}168169static GEN170_Fp_red(void *data, GEN x) { return Fp_red(x, (GEN)data); }171172/* p not necessarily prime */173static const struct bb_hermite Fp_hermite=174{_Fp_add,_Fp_neg,_Fp_mul,_Fp_extgcd,_Fp_rann,_Fp_lquo,_Fp_unit,_Fp_equal0,_Fp_equal1,_Fp_s,_Fp_red};175176static const struct bb_hermite*177get_Fp_hermite(void **data, GEN p)178{179*data = (void*)p; return &Fp_hermite;180}181182static void183gen_redcol(GEN C, long lim, void* data, const struct bb_hermite *R)184{185long i;186for (i=1; i<=lim; i++)187if (!R->equal0(gel(C,i)))188gel(C,i) = R->red(data, gel(C,i));189}190191static GEN192/* return NULL if a==0 */193/* assume C*a is zero after lim */194gen_rightmulcol(GEN C, GEN a, long lim, int fillzeros, void* data, const struct bb_hermite *R)195{196GEN Ca,zero;197long i;198if (R->equal1(a)) return C;199if (R->equal0(a)) return NULL;200Ca = cgetg(lg(C),t_COL);201for (i=1; i<=lim; i++)202gel(Ca,i) = R->mul(data, gel(C,i), a);203if (fillzeros && lim+1 < lg(C))204{205zero = R->s(data,0);206for (i=lim+1; i<lg(C); i++)207gel(Ca,i) = zero;208}209return Ca;210}211212static void213/* C1 <- C1 + C2 */214/* assume C2[i]==0 for i>lim */215gen_addcol(GEN C1, GEN C2, long lim, void* data, const struct bb_hermite *R)216{217long i;218for (i=1; i<=lim; i++)219gel(C1,i) = R->add(data, gel(C1,i), gel(C2,i));220}221222static void223/* H[,i] <- H[,i] + C*a */224/* assume C is zero after lim */225gen_addrightmul(GEN H, GEN C, GEN a, long i, long lim, void* data, const struct bb_hermite *R)226{227GEN Ca;228if (R->equal0(a)) return;229Ca = gen_rightmulcol(C, a, lim, 0, data, R);230gen_addcol(gel(H,i), Ca, lim, data, R);231}232233static GEN234gen_zerocol(long n, void* data, const struct bb_hermite *R)235{236GEN C = cgetg(n+1,t_COL), zero = R->s(data, 0);237long i;238for (i=1; i<=n; i++) gel(C,i) = zero;239return C;240}241242static GEN243gen_zeromat(long m, long n, void* data, const struct bb_hermite *R)244{245GEN M = cgetg(n+1,t_MAT);246long i;247for (i=1; i<=n; i++) gel(M,i) = gen_zerocol(m, data, R);248return M;249}250251static GEN252gen_colei(long n, long i, void* data, const struct bb_hermite *R)253{254GEN C = cgetg(n+1,t_COL), zero = R->s(data, 0);255long j;256for (j=1; j<=n; j++)257if (i!=j) gel(C,j) = zero;258else gel(C,j) = R->s(data,1);259return C;260}261262static GEN263gen_matid_hermite(long n, void* data, const struct bb_hermite *R)264{265GEN M = cgetg(n+1,t_MAT);266long i;267for (i=1; i<=n; i++) gel(M,i) = gen_colei(n, i, data, R);268return M;269}270271static GEN272gen_matmul_hermite(GEN A, GEN B, void* data, const struct bb_hermite *R)273{274GEN M,sum,prod,zero = R->s(data,0);275long a,b,c,c2,i,j,k;276RgM_dimensions(A,&a,&c);277RgM_dimensions(B,&c2,&b);278if (c!=c2) pari_err_DIM("gen_matmul_hermite");279M = cgetg(b+1,t_MAT);280for (j=1; j<=b; j++)281{282gel(M,j) = cgetg(a+1,t_COL);283for (i=1; i<=a; i++)284{285sum = zero;286for (k=1; k<=c; k++)287{288prod = R->mul(data, gcoeff(A,i,k), gcoeff(B,k,j));289sum = R->add(data, sum, prod);290}291gcoeff(M,i,j) = sum;292}293gen_redcol(gel(M,j), a, data, R);294}295return M;296}297298static void299/* U = [u1,u2]~, C <- A*u1 + B*u2 */300/* assume both A, B and C are zero after lim */301gen_rightlincomb(GEN A, GEN B, GEN U, GEN *C, long lim, void* data, const struct bb_hermite *R)302{303GEN Au1, Bu2;304Au1 = gen_rightmulcol(A, gel(U,1), lim, 1, data, R);305Bu2 = gen_rightmulcol(B, gel(U,2), lim, 1, data, R);306if (!Au1 && !Bu2) { *C = gen_zerocol(lg(A)-1, data, R); return; }307if (!Au1) { *C = Bu2; return; }308if (!Bu2) { *C = Au1; return; }309gen_addcol(Au1, Bu2, lim, data, R);310*C = Au1;311}312313static void314/* (H[,i] | H[,j]) <- (H[,i] | H[,j]) * U */315/* assume both columns are zero after lim */316gen_elem(GEN H, GEN U, long i, long j, long lim, void* data, const struct bb_hermite *R)317{318GEN Hi, Hj;319Hi = shallowcopy(gel(H,i));320Hj = shallowcopy(gel(H,j));321gen_rightlincomb(Hi, Hj, gel(U,1), &gel(H,i), lim, data, R);322gen_rightlincomb(Hi, Hj, gel(U,2), &gel(H,j), lim, data, R);323}324325static int326/* assume C is zero after lim */327gen_is_zerocol(GEN C, long lim, void* data, const struct bb_hermite *R)328{329long i;330(void) data;331for (i=1; i<=lim; i++)332if (!R->equal0(gel(C,i))) return 0;333return 1;334}335336/* The mkop* functions return NULL if the corresponding operation is the identity */337338static GEN339/* Ci <- Ci + Cj*a */340mkoptransv(long i, long j, GEN a, void* data, const struct bb_hermite *R)341{342a = R->red(data,a);343if (R->equal0(a)) return NULL;344return mkvec2(mkvecsmall2(i,j),a);345}346347static GEN348/* (Ci|Cj) <- (Ci|Cj)*U */349mkopU(long i, long j, GEN U, void* data, const struct bb_hermite *R)350{351if (R->equal1(gcoeff(U,1,1)) && R->equal0(gcoeff(U,1,2))352&& R->equal1(gcoeff(U,2,2))) return mkoptransv(i,j,gcoeff(U,2,1),data,R);353return mkvec2(mkvecsmall3(i,j,0),U);354}355356static GEN357/* Ci <- Ci*u */358mkopmul(long i, GEN u, const struct bb_hermite *R)359{360if (R->equal1(u)) return NULL;361return mkvec2(mkvecsmall(i),u);362}363364static GEN365/* Ci <-> Cj */366mkopswap(long i, long j)367{368return mkvec(mkvecsmall2(i,j));369}370371/* M: t_MAT. Apply the operation op to M by right multiplication. */372static void373gen_rightapply(GEN M, GEN op, void* data, const struct bb_hermite *R)374{375GEN M2, ind, X;376long i, j, m = lg(gel(M,1))-1;377switch (typ(op))378{379case t_VECSMALL:380M2 = vecpermute(M,op);381for (i=1; i<lg(M); i++) gel(M,i) = gel(M2,i);382return;383case t_VEC:384ind = gel(op,1);385switch (lg(op))386{387case 2:388swap(gel(M,ind[1]),gel(M,ind[2]));389return;390case 3:391X = gel(op,2);392i = ind[1];393switch (lg(ind))394{395case 2:396gel(M,i) = gen_rightmulcol(gel(M,i), X, m, 0, data, R);397gen_redcol(gel(M,i), m, data, R);398return;399case 3:400gen_addrightmul(M, gel(M,ind[2]), X, i, m, data, R);401gen_redcol(gel(M,i), m, data, R);402return;403case 4:404j = ind[2];405gen_elem(M, X, i, j, m, data, R);406gen_redcol(gel(M,i), m, data, R);407gen_redcol(gel(M,j), m, data, R);408return;409}410}411}412}413414/* C: t_COL. Apply the operation op to C by left multiplication. */415static void416gen_leftapply(GEN C, GEN op, void* data, const struct bb_hermite *R)417{418GEN C2, ind, X;419long i, j;420switch (typ(op))421{422case t_VECSMALL:423C2 = vecpermute(C,perm_inv(op));424for (i=1; i<lg(C); i++) gel(C,i) = gel(C2,i);425return;426case t_VEC:427ind = gel(op,1);428switch (lg(op))429{430case 2:431swap(gel(C,ind[1]),gel(C,ind[2]));432return;433case 3:434X = gel(op,2);435i = ind[1];436switch (lg(ind))437{438case 2:439gel(C,i) = R->mul(data, X, gel(C,i));440gel(C,i) = R->red(data, gel(C,i));441return;442case 3:443j = ind[2];444if (R->equal0(gel(C,i))) return;445gel(C,j) = R->add(data, gel(C,j), R->mul(data, X, gel(C,i)));446return;447case 4:448j = ind[2];449C2 = gen_matmul_hermite(X, mkmat(mkcol2(gel(C,i),gel(C,j))), data, R);450gel(C,i) = gcoeff(C2,1,1);451gel(C,j) = gcoeff(C2,2,1);452return;453}454}455}456}457458/* \prod_i det ops[i]. Only makes sense if R is commutative. */459static GEN460gen_detops(GEN ops, void* data, const struct bb_hermite *R)461{462GEN d = R->s(data,1);463long i, l = lg(ops);464for (i = 1; i < l; i++)465{466GEN X, op = gel(ops,i);467switch (typ(op))468{469case t_VECSMALL:470if (perm_sign(op) < 0) d = R->neg(data,d);471break;472case t_VEC:473switch (lg(op))474{475case 2:476d = R->neg(data,d);477break;478case 3:479X = gel(op,2);480switch (lg(gel(op,1)))481{482case 2:483d = R->mul(data, d, X);484d = R->red(data, d);485break;486case 4:487{488GEN A = gcoeff(X,1,1), B = gcoeff(X,1,2);489GEN C = gcoeff(X,2,1), D = gcoeff(X,2,2);490GEN AD = R->mul(data,A,D);491GEN BC = R->mul(data,B,C);492d = R->mul(data, d, R->add(data, AD, R->neg(data,BC)));493d = R->red(data, d);494break;495}496}497break;498}499break;500}501}502return d;503}504505static int506gen_is_inv(GEN x, void* data, const struct bb_hermite *R)507{508GEN u = R->unit(data, x);509if (!u) return R->equal1(x);510return R->equal1(gel(u,1));511}512513static long514gen_last_inv_diago(GEN A, void* data, const struct bb_hermite *R)515{516long i,m,n,j;517RgM_dimensions(A,&m,&n);518for (i=1,j=n-m+1; i<=m; i++,j++)519if (!gen_is_inv(gcoeff(A,i,j),data,R)) return i-1;520return m;521}522523static GEN524/* remove_zerocols: 0 none, 1 until square, 2 all */525/* early abort: if not right-invertible, abort, return NULL, and set ops to the526* noninvertible pivot */527gen_howell_i(GEN A, long remove_zerocols, long permute_zerocols, long early_abort, long only_triangular, GEN* ops, void *data, const struct bb_hermite *R)528{529pari_sp av = avma, av1;530GEN H,U,piv=gen_0,u,q,a,perm,iszero,C,zero=R->s(data,0),d,g,r,op,one=R->s(data,1);531long m,n,i,j,s,si,i2,si2,nbz,lim,extra,maxop=0,nbop=0,lastinv=0;532int smallop;533534av1 = avma;535536RgM_dimensions(A,&m,&n);537if (early_abort && n<m)538{539if (ops) *ops = zero;540return NULL;541}542if (n<m+1)543{544extra = m+1-n;545H = shallowmatconcat(mkvec2(gen_zeromat(m,extra,data,R),A));546}547else548{549extra = 0;550H = RgM_shallowcopy(A);551}552RgM_dimensions(H,&m,&n);553s = n-m; /* shift */554555if(ops)556{557maxop = m*n + (m*(m+1))/2 + 1;558*ops = zerovec(maxop); /* filled with placeholders so gerepile can handle it */559}560561/* put in triangular form */562for (i=m,si=s+m; i>0 && si>extra; i--,si--) /* si = s+i */563{564if (R->red) gcoeff(H,i,si) = R->red(data, gcoeff(H,i,si));565/* bottom-right diagonal */566for (j = extra+1; j < si; j++)567{568if (R->red) gcoeff(H,i,j) = R->red(data, gcoeff(H,i,j));569if (R->equal0(gcoeff(H,i,j))) continue;570U = R->extgcd(data, gcoeff(H,i,j), gcoeff(H,i,si), &smallop);571d = gel(U,1);572U = gel(U,2);573if (n>10)574{575/* normalize diagonal coefficient -> faster reductions on this row */576u = R->unit(data, d);577if (u)578{579g = gel(u,1);580u = gel(u,2);581gcoeff(U,1,2) = R->mul(data, gcoeff(U,1,2), u);582gcoeff(U,2,2) = R->mul(data, gcoeff(U,2,2), u);583d = g;584}585}586gen_elem(H, U, j, si, i-1, data, R);587if (ops)588{589op = mkopU(j,si,U,data,R);590if (op) { nbop++; gel(*ops, nbop) = op; }591}592gcoeff(H,i,j) = zero;593gcoeff(H,i,si) = d;594if (R->red && !smallop)595{596gen_redcol(gel(H,si), i-1, data, R);597gen_redcol(gel(H,j), i-1, data, R);598}599}600601if (early_abort)602{603d = gcoeff(H,i,si);604u = R->unit(data, d);605if (u) d = gel(u,1);606if (!R->equal1(d))607{608if (ops) *ops = d;609return NULL;610}611}612613if (gc_needed(av,1))614{615if (DEBUGMEM>1) pari_warn(warnmem,"gen_howell[1]. i=%ld",i);616gerepileall(av1,ops?2:1,&H,ops);617}618}619620if (!ops)621lastinv = gen_last_inv_diago(H, data, R);622623/* put in reduced Howell form */624if (!only_triangular)625{626for (i=m,si=s+m; i>0; i--,si--) /* si = s+i */627{628/* normalize diagonal coefficient */629if (i<=lastinv) /* lastinv>0 => !ops */630gcoeff(H,i,si) = one;631else632{633u = R->unit(data,gcoeff(H,i,si));634if (u)635{636g = gel(u,1);637u = gel(u,2);638gel(H,si) = gen_rightmulcol(gel(H,si), u, i-1, 1, data, R);639gcoeff(H,i,si) = g;640if (R->red) gen_redcol(gel(H,si), i-1, data, R);641if (ops)642{643op = mkopmul(si,u,R);644if (op) { nbop++; gel(*ops,nbop) = op; }645}646}647else if (R->red) gcoeff(H,i,si) = R->red(data, gcoeff(H,i,si));648}649piv = gcoeff(H,i,si);650651/* reduce above diagonal */652if (!R->equal0(piv))653{654C = gel(H,si);655for (j=si+1; j<=n; j++)656{657if (i<=lastinv) /* lastinv>0 => !ops */658gcoeff(H,i,j) = zero;659else660{661gcoeff(H,i,j) = R->red(data, gcoeff(H,i,j));662if (R->equal1(piv)) { q = gcoeff(H,i,j); r = zero; }663else q = R->lquo(data, gcoeff(H,i,j), piv, &r);664q = R->neg(data,q);665gen_addrightmul(H, C, q, j, i-1, data, R);666if (ops)667{668op = mkoptransv(j,si,q,data,R);669if (op) { nbop++; gel(*ops,nbop) = op; }670}671gcoeff(H,i,j) = r;672}673}674}675676/* ensure Howell property */677if (i>1)678{679a = R->rann(data, piv);680if (!R->equal0(a))681{682gel(H,1) = gen_rightmulcol(gel(H,si), a, i-1, 1, data, R);683if (gel(H,1) == gel(H,si)) gel(H,1) = shallowcopy(gel(H,1)); /* in case rightmulcol cheated */684if (ops)685{686op = mkoptransv(1,si,a,data,R);687if (op) { nbop++; gel(*ops,nbop) = op; }688}689for (i2=i-1,si2=s+i2; i2>0; i2--,si2--)690{691if (R->red) gcoeff(H,i2,1) = R->red(data, gcoeff(H,i2,1));692if (R->equal0(gcoeff(H,i2,1))) continue;693if (R->red) gcoeff(H,i2,si2) = R->red(data, gcoeff(H,i2,si2));694if (R->equal0(gcoeff(H,i2,si2)))695{696swap(gel(H,1), gel(H,si2));697if (ops) { nbop++; gel(*ops,nbop) = mkopswap(1,si2); }698continue;699}700U = R->extgcd(data, gcoeff(H,i2,1), gcoeff(H,i2,si2), &smallop);701d = gel(U,1);702U = gel(U,2);703gen_elem(H, U, 1, si2, i2-1, data, R);704if (ops)705{706op = mkopU(1,si2,U,data,R);707if (op) { nbop++; gel(*ops,nbop) = op; }708}709gcoeff(H,i2,1) = zero;710gcoeff(H,i2,si2) = d;711if (R->red && !smallop)712{713gen_redcol(gel(H,si2), i2, data, R);714gen_redcol(gel(H,1), i2-1, data, R);715}716}717}718}719720if (gc_needed(av,1))721{722if (DEBUGMEM>1) pari_warn(warnmem,"gen_howell[2]. i=%ld",i);723gerepileall(av1,ops?3:2,&H,&piv,ops);724}725}726}727728if (R->red)729for (j=1; j<=n; j++)730{731lim = maxss(0,m-n+j);732gen_redcol(gel(H,j), lim, data, R);733for (i=lim+1; i<=m; i++) gcoeff(H,i,j) = zero;734}735736/* put zero columns first */737iszero = cgetg(n+1,t_VECSMALL);738739nbz = 0;740for (i=1; i<=n; i++)741{742iszero[i] = gen_is_zerocol(gel(H,i), maxss(0,m-n+i), data, R);743if (iszero[i]) nbz++;744}745746j = 1;747if (permute_zerocols)748{749perm = cgetg(n+1, t_VECSMALL);750for (i=1; i<=n; i++)751if (iszero[i])752{753perm[j] = i;754j++;755}756}757else perm = cgetg(n-nbz+1, t_VECSMALL);758for (i=1; i<=n; i++)759if (!iszero[i])760{761perm[j] = i;762j++;763}764765if (permute_zerocols || remove_zerocols==2) H = vecpermute(H, perm);766if (permute_zerocols && remove_zerocols==2) H = vecslice(H, nbz+1, n);767if (remove_zerocols==1) H = vecslice(H, s+1, n);768if (permute_zerocols && ops) { nbop++; gel(*ops,nbop) = perm; }769770if (ops) { setlg(*ops, nbop+1); } /* should have nbop <= maxop */771772return H;773}774775static GEN776gen_howell(GEN A, long remove_zerocols, long permute_zerocols, long early_abort, long only_triangular, GEN* ops, void *data, const struct bb_hermite *R)777{778pari_sp av = avma;779GEN H = gen_howell_i(A, remove_zerocols, permute_zerocols, early_abort, only_triangular, ops, data, R);780gerepileall(av, ops?2:1, &H, ops);781return H;782}783784static GEN785gen_matimage(GEN A, GEN* U, void *data, const struct bb_hermite *R)786{787GEN ops, H;788if (U)789{790pari_sp av = avma, av1;791long m, n, i, r, n2, pergc;792RgM_dimensions(A,&m,&n);793H = gen_howell_i(A, 2, 1, 0, 0, &ops, data, R);794av1 = avma;795r = lg(H)-1;796*U = shallowmatconcat(mkvec2(gen_zeromat(n, maxss(0,m-n+1), data, R), gen_matid_hermite(n, data, R)));797n2 = lg(*U)-1;798pergc = maxss(m,n);799for (i=1; i<lg(ops); i++)800{801gen_rightapply(*U, gel(ops,i), data, R);802if (!(i%pergc) && gc_needed(av,1))803{804if (DEBUGMEM>1) pari_warn(warnmem,"gen_matimage. i=%ld",i);805gerepileall(av1,1,U);806}807}808if (r<n2) *U = vecslice(*U, n2-r+1, n2);809gerepileall(av, 2, &H, U);810return H;811}812else return gen_howell(A, 2, 0, 0, 0, NULL, data, R);813}814815static GEN816/* H in true Howell form: no zero columns */817gen_kernel_howell(GEN H, void *data, const struct bb_hermite *R)818{819GEN K, piv, FK;820long m, n, j, j2, i;821RgM_dimensions(H,&m,&n);822K = gen_zeromat(n, n, data, R);823for (j=n,i=m; j>0; j--)824{825while (R->equal0(gcoeff(H,i,j))) i--;826piv = gcoeff(H,i,j);827if (R->equal0(piv)) continue;828gcoeff(K,j,j) = R->rann(data, piv);829if (j<n)830{831FK = gen_matmul_hermite(matslice(H,i,i,j+1,n), matslice(K, j+1, n, j+1, n), data, R);832for (j2=j+1; j2<=n; j2++)833gcoeff(K,j,j2) = R->neg(data, R->lquo(data, gcoeff(FK,1,j2-j), piv, NULL));834/* remainder has to be zero */835}836}837return K;838}839840static GEN841/* (H,ops) Howell form of A, n = number of columns of A, return a kernel of A */842gen_kernel_from_howell(GEN H, GEN ops, long n, void *data, const struct bb_hermite *R)843{844pari_sp av;845GEN K, KH, zC;846long m, r, n2, nbz, i, o, extra, j;847RgM_dimensions(H,&m,&r);848if (!r) return gen_matid_hermite(n, data, R); /* zerology: what if 0==1 in R? */849n2 = maxss(n,m+1);850extra = n2-n;851nbz = n2-r;852/* compute kernel of augmented matrix */853KH = gen_kernel_howell(H, data, R);854zC = gen_zerocol(nbz, data, R);855K = cgetg(nbz+r+1, t_MAT);856for (i=1; i<=nbz; i++)857gel(K,i) = gen_colei(nbz+r, i, data, R);858for (i=1; i<=r; i++)859gel(K,nbz+i) = shallowconcat(zC, gel(KH,i));860for (i=1; i<lg(K); i++)861{862av = avma;863for (o=lg(ops)-1; o>0; o--)864gen_leftapply(gel(K,i), gel(ops,o), data, R);865gen_redcol(gel(K,i), nbz+r, data, R);866gerepileall(av, 1, &gel(K,i));867}868/* deduce kernel of original matrix */869K = rowpermute(K, cyclic_perm(n2,extra));870K = gen_howell_i(K, 2, 0, 0, 0, NULL, data, R);871for (j=lg(K)-1, i=n2; j>0; j--)872{873while (R->equal0(gcoeff(K,i,j))) i--;874if (i<=n) return matslice(K, 1, n, 1, j);875}876return cgetg(1,t_MAT);877}878879/* not GC-clean */880static GEN881gen_kernel(GEN A, GEN* im, void *data, const struct bb_hermite *R)882{883pari_sp av = avma;884long n = lg(A)-1;885GEN H, ops, K;886H = gen_howell_i(A, 2, 1, 0, 0, &ops, data, R);887gerepileall(av,2,&H,&ops);888K = gen_kernel_from_howell(H, ops, n, data, R);889if (im) *im = H;890return K;891}892893/* right inverse, not GC-clean */894static GEN895gen_inv(GEN A, void* data, const struct bb_hermite *R)896{897pari_sp av;898GEN ops, H, U, un=R->s(data,1);899long m,n,j,o,n2;900RgM_dimensions(A,&m,&n);901av = avma;902H = gen_howell_i(A, 0, 0, 1, 0, &ops, data, R);903if (!H) pari_err_INV("gen_inv", ops);904n2 = lg(H)-1;905ops = gerepilecopy(av,ops); /* get rid of H */906U = gen_zeromat(n2, m, data, R);907for (j=1; j<=m; j++)908gcoeff(U,j+n2-m,j) = un;909for (j=1; j<=m; j++)910{911av = avma;912for (o=lg(ops)-1; o>0; o--)913gen_leftapply(gel(U,j), gel(ops,o), data, R);914gen_redcol(gel(U,j), n2, data, R);915gerepileall(av, 1, &gel(U,j));916}917if (n2>n) U = rowslice(U, n2-n+1, n2);918return U;919}920921/*922H true Howell form (no zero columns).923Compute Z = Y - HX canonical representative of R^m mod H(R^n)924*/925static GEN926gen_reduce_mod_howell(GEN H, GEN Y, GEN *X, void* data, const struct bb_hermite *R)927{928pari_sp av = avma;929long m,n,i,j;930GEN Z, q, r, C;931RgM_dimensions(H,&m,&n);932if (X) *X = gen_zerocol(n,data,R);933Z = shallowcopy(Y);934i = m;935for (j=n; j>0; j--)936{937while (R->equal0(gcoeff(H,i,j))) i--;938q = R->lquo(data, gel(Z,i), gcoeff(H,i,j), &r);939gel(Z,i) = r;940C = gen_rightmulcol(gel(H,j), R->neg(data,q), i, 0, data, R);941if (C) gen_addcol(Z, C, i-1, data, R);942if (X) gel(*X,j) = q;943}944gen_redcol(Z, lg(Z)-1, data, R);945if (X)946{947gerepileall(av, 2, &Z, X);948return Z;949}950return gerepilecopy(av, Z);951}952953/* H: Howell form of A */954/* (m,n): dimensions of original matrix A */955/* return NULL if no solution */956static GEN957gen_solve_from_howell(GEN H, GEN ops, long m, long n, GEN Y, void* data, const struct bb_hermite *R)958{959pari_sp av = avma;960GEN Z, X;961long n2, mH, nH, i;962RgM_dimensions(H,&mH,&nH);963n2 = maxss(n,m+1);964Z = gen_reduce_mod_howell(H, Y, &X, data, R);965if (!gen_is_zerocol(Z,m,data,R)) { set_avma(av); return NULL; }966967X = shallowconcat(zerocol(n2-nH),X);968for (i=lg(ops)-1; i>0; i--) gen_leftapply(X, gel(ops,i), data, R);969X = vecslice(X, n2-n+1, n2);970gen_redcol(X, n, data, R);971return gerepilecopy(av, X);972}973974/* return NULL if no solution, not GC-clean */975static GEN976gen_solve(GEN A, GEN Y, GEN* K, void* data, const struct bb_hermite *R)977{978GEN H, ops, X;979long m,n;980981RgM_dimensions(A,&m,&n);982if (!n) m = lg(Y)-1;983H = gen_howell_i(A, 2, 1, 0, 0, &ops, data, R);984X = gen_solve_from_howell(H, ops, m, n, Y, data, R);985if (!X) return NULL;986if (K) *K = gen_kernel_from_howell(H, ops, n, data, R);987return X;988}989990GEN991matimagemod(GEN A, GEN d, GEN* U)992{993void *data;994const struct bb_hermite* R;995if (typ(A)!=t_MAT || !RgM_is_ZM(A)) pari_err_TYPE("matimagemod", A);996if (typ(d)!=t_INT) pari_err_TYPE("matimagemod", d);997if (signe(d)<=0) pari_err_DOMAIN("matimagemod", "d", "<=", gen_0, d);998if (equali1(d)) return cgetg(1,t_MAT);999R = get_Fp_hermite(&data, d);1000return gen_matimage(A, U, data, R);1001}10021003/* for testing purpose */1004/*1005GEN1006ZM_hnfmodid2(GEN A, GEN d)1007{1008pari_sp av = avma;1009void *data;1010long i;1011const struct bb_hermite* R = get_Fp_hermite(&data, d);1012GEN H;1013if (typ(A)!=t_MAT || !RgM_is_ZM(A)) pari_err_TYPE("ZM_hnfmodid2", A);1014if (typ(d)!=t_INT) pari_err_TYPE("ZM_hnfmodid2", d);1015H = gen_howell_i(A, 1, 0, 0, 0, NULL, data, R);1016for (i=1; i<lg(H); i++)1017if (!signe(gcoeff(H,i,i))) gcoeff(H,i,i) = d;1018return gerepilecopy(av,H);1019}1020*/10211022GEN1023matdetmod(GEN A, GEN d)1024{1025pari_sp av = avma;1026void *data;1027const struct bb_hermite* R;1028long n = lg(A)-1, i;1029GEN D, H, ops;1030if (typ(A)!=t_MAT || !RgM_is_ZM(A)) pari_err_TYPE("matdetmod", A);1031if (typ(d)!=t_INT) pari_err_TYPE("matdetmod", d);1032if (signe(d)<=0) pari_err_DOMAIN("matdetmod", "d", "<=", gen_0, d);1033if (!n) return equali1(d) ? gen_0 : gen_1;1034if (n != nbrows(A)) pari_err_DIM("matdetmod");1035if (equali1(d)) return gen_0;1036R = get_Fp_hermite(&data, d);1037H = gen_howell_i(A, 1, 0, 0, 1, &ops, data, R);1038D = gen_detops(ops, data, R);1039D = Fp_inv(D, d);1040for (i = 1; i <= n; i++) D = Fp_mul(D, gcoeff(H,i,i), d);1041return gerepileuptoint(av, D);1042}10431044GEN1045matkermod(GEN A, GEN d, GEN* im)1046{1047pari_sp av = avma;1048void *data;1049const struct bb_hermite* R;1050long m,n;1051GEN K;1052if (typ(A)!=t_MAT || !RgM_is_ZM(A)) pari_err_TYPE("matkermod", A);1053if (typ(d)!=t_INT) pari_err_TYPE("matkermod", d);1054if (signe(d)<=0) pari_err_DOMAIN("makermod", "d", "<=", gen_0, d);1055if (equali1(d)) return cgetg(1,t_MAT);1056R = get_Fp_hermite(&data, d);1057RgM_dimensions(A,&m,&n);1058if (!im && m>2*n) /* TODO tune treshold */1059A = shallowtrans(matimagemod(shallowtrans(A),d,NULL));1060K = gen_kernel(A, im, data, R);1061gerepileall(av,im?2:1,&K,im);1062return K;1063}10641065/* left inverse */1066GEN1067matinvmod(GEN A, GEN d)1068{1069pari_sp av = avma;1070void *data;1071const struct bb_hermite* R;1072GEN U;1073if (typ(A)!=t_MAT || !RgM_is_ZM(A)) pari_err_TYPE("matinvmod", A);1074if (typ(d)!=t_INT) pari_err_TYPE("matinvmod", d);1075if (signe(d)<=0) pari_err_DOMAIN("matinvmod", "d", "<=", gen_0, d);1076if (equali1(d)) {1077long m,n;1078RgM_dimensions(A,&m,&n);1079if (m<n) pari_err_INV("matinvmod",A);1080return zeromatcopy(n,m);1081}1082R = get_Fp_hermite(&data, d);1083U = gen_inv(shallowtrans(A), data, R);1084return gerepilecopy(av, shallowtrans(U));1085}10861087/* assume all D[i]>0, not GC-clean */1088static GEN1089matsolvemod_finite(GEN M, GEN D, GEN Y, long flag)1090{1091void *data;1092const struct bb_hermite* R;1093GEN X, K, d;1094long m, n;10951096RgM_dimensions(M,&m,&n);1097if (typ(D)==t_COL)1098{ /* create extra variables for the D[i] */1099GEN MD;1100long i, i2, extra = 0;1101d = gen_1;1102for (i=1; i<lg(D); i++) d = lcmii(d,gel(D,i));1103for (i=1; i<lg(D); i++) if (!equalii(gel(D,i),d)) extra++;1104MD = cgetg(extra+1,t_MAT);1105i2 = 1;1106for (i=1; i<lg(D); i++)1107if (!equalii(gel(D,i),d))1108{1109gel(MD,i2) = Rg_col_ei(gel(D,i),m,i);1110i2++;1111}1112M = shallowconcat(M,MD);1113}1114else d = D;11151116R = get_Fp_hermite(&data, d);1117X = gen_solve(M, Y, flag? &K: NULL, data, R);1118if (!X) return gen_0;1119X = vecslice(X,1,n);11201121if (flag)1122{1123K = rowslice(K,1,n);1124K = hnfmodid(shallowconcat(zerocol(n),K),d);1125X = mkvec2(X,K);1126}1127return X;1128}11291130/* Return a solution of congruence system sum M[i,j] X_j = Y_i mod D_i1131* If pU != NULL, put in *pU a Z-basis of the homogeneous system.1132* Works for all non negative D_i but inefficient compared to1133* matsolvemod_finite; to be used only when one D_i is 0 */1134static GEN1135gaussmoduloall(GEN M, GEN D, GEN Y, GEN *pU)1136{1137pari_sp av = avma;1138long n, m, j, l, lM = lg(M);1139GEN delta, H, U, u1, u2, x;11401141if (lM == 1)1142{1143long lY = 0;1144switch(typ(Y))1145{1146case t_INT: break;1147case t_COL: lY = lg(Y); break;1148default: pari_err_TYPE("gaussmodulo",Y);1149}1150switch(typ(D))1151{1152case t_INT: break;1153case t_COL: if (lY && lY != lg(D)) pari_err_DIM("gaussmodulo");1154break;1155default: pari_err_TYPE("gaussmodulo",D);1156}1157if (pU) *pU = cgetg(1, t_MAT);1158return cgetg(1,t_COL);1159}1160n = nbrows(M);1161switch(typ(D))1162{1163case t_COL:1164if (lg(D)-1!=n) pari_err_DIM("gaussmodulo");1165delta = diagonal_shallow(D); break;1166case t_INT: delta = scalarmat_shallow(D,n); break;1167default: pari_err_TYPE("gaussmodulo",D);1168return NULL; /* LCOV_EXCL_LINE */1169}1170switch(typ(Y))1171{1172case t_INT: Y = const_col(n, Y); break;1173case t_COL:1174if (lg(Y)-1!=n) pari_err_DIM("gaussmodulo");1175break;1176default: pari_err_TYPE("gaussmodulo",Y);1177return NULL; /* LCOV_EXCL_LINE */1178}1179H = ZM_hnfall_i(shallowconcat(M,delta), &U, 1);1180Y = hnf_solve(H,Y); if (!Y) return gen_0;1181l = lg(H); /* may be smaller than lM if some moduli are 0 */1182n = l-1;1183m = lg(U)-1 - n;1184u1 = cgetg(m+1,t_MAT);1185u2 = cgetg(n+1,t_MAT);1186for (j=1; j<=m; j++) { GEN c = gel(U,j); setlg(c,lM); gel(u1,j) = c; }1187U += m;1188for (j=1; j<=n; j++) { GEN c = gel(U,j); setlg(c,lM); gel(u2,j) = c; }1189/* (u1 u2)1190* (M D) (* * ) = (0 H) */1191u1 = ZM_lll(u1, 0.75, LLL_INPLACE);1192Y = ZM_ZC_mul(u2,Y);1193x = ZC_reducemodmatrix(Y, u1);1194if (!pU) x = gerepileupto(av, x);1195else1196{1197gerepileall(av, 2, &x, &u1);1198*pU = u1;1199}1200return x;1201}1202/* to be used when one D_i is 0 */1203static GEN1204msolvemod0(GEN M, GEN D, GEN Y, long flag)1205{1206pari_sp av = avma;1207GEN y, x, U;1208if (!flag) return gaussmoduloall(M,D,Y,NULL);1209y = cgetg(3,t_VEC);1210x = gaussmoduloall(M,D,Y,&U);1211if (x == gen_0) { set_avma(av); return gen_0; }1212gel(y,1) = x;1213gel(y,2) = U; return y;12141215}1216GEN1217matsolvemod(GEN M, GEN D, GEN Y, long flag)1218{1219pari_sp av = avma;1220long m, n, i, char0 = 0;1221if (typ(M)!=t_MAT || !RgM_is_ZM(M)) pari_err_TYPE("matsolvemod (M)",M);1222RgM_dimensions(M,&m,&n);1223if (typ(D)!=t_INT && (typ(D)!=t_COL || !RgV_is_ZV(D)))1224pari_err_TYPE("matsolvemod (D)",D);1225if (n)1226{ if (typ(D)==t_COL && lg(D)!=m+1) pari_err_DIM("matsolvemod [1]"); }1227else1228{ if (typ(D)==t_COL) m = lg(D)-1; }1229if (typ(Y)==t_INT)1230Y = const_col(m,Y);1231else if (typ(Y)!=t_COL || !RgV_is_ZV(Y)) pari_err_TYPE("matsolvemod (Y)",Y);1232if (!n && !m) m = lg(Y)-1;1233else if (m != lg(Y)-1) pari_err_DIM("matsolvemod [2]");1234if (typ(D)==t_INT)1235{1236if (signe(D)<0) pari_err_DOMAIN("matsolvemod","D","<",gen_0,D);1237if (!signe(D)) char0 = 1;1238}1239else /*typ(D)==t_COL*/1240for (i=1; i<=m; i++)1241{1242if (signe(gel(D,i))<0)1243pari_err_DOMAIN("matsolvemod","D[i]","<",gen_0,gel(D,i));1244if (!signe(gel(D,i))) char0 = 1;1245}1246return gerepilecopy(av, char0? msolvemod0(M,D,Y,flag)1247: matsolvemod_finite(M,D,Y,flag));1248}1249GEN1250gaussmodulo2(GEN M, GEN D, GEN Y) { return matsolvemod(M,D,Y, 1); }1251GEN1252gaussmodulo(GEN M, GEN D, GEN Y) { return matsolvemod(M,D,Y, 0); }125312541255