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#include "pari.h"15#include "paripriv.h"1617#define DEBUGLEVEL DEBUGLEVEL_mat1819/********************************************************************/20/** **/21/** REDUCTION **/22/** **/23/********************************************************************/24/* z in Z^n, return lift(Col(z) * Mod(1,p)) */25GEN26FpC_red(GEN x, GEN p)27{ pari_APPLY_type(t_COL, modii(gel(x,i), p)) }2829/* z in Z^n, return lift(Vec(z) * Mod(1,p)) */30GEN31FpV_red(GEN x, GEN p)32{ pari_APPLY_type(t_VEC, modii(gel(x,i), p)) }3334GEN35FpC_center(GEN x, GEN p, GEN pov2)36{ pari_APPLY_type(t_COL, Fp_center(gel(x,i), p, pov2)) }3738/* assume 0 <= u < p and ps2 = p>>1 */39INLINE void40Fp_center_inplace(GEN u, GEN p, GEN ps2)41{ if (abscmpii(u,ps2) > 0) subiiz(u,p,u); }4243void44FpC_center_inplace(GEN z, GEN p, GEN ps2)45{46long i,l = lg(z);47for (i=1; i<l; i++)48Fp_center_inplace(gel(z,i), p, ps2);49}5051GEN52Flv_center(GEN x, ulong p, ulong ps2)53{ pari_APPLY_ulong(Fl_center(uel(x,i),p,ps2)) }5455/* z in Mat m,n(Z), return lift(z * Mod(1,p)) */56GEN57FpM_red(GEN x, GEN p)58{ pari_APPLY_same(FpC_red(gel(x,i), p)) }5960GEN61FpM_center(GEN x, GEN p, GEN pov2)62{ pari_APPLY_same(FpC_center(gel(x,i), p,pov2)) }6364void65FpM_center_inplace(GEN z, GEN p, GEN pov2)66{67long i, l = lg(z);68for (i=1; i<l; i++) FpC_center_inplace(gel(z,i), p, pov2);69}70GEN71Flm_center(GEN x, ulong p, ulong ps2)72{ pari_APPLY_same(Flv_center(gel(x,i), p,ps2)) }7374GEN75random_FpV(long d, GEN p)76{77long i;78GEN y = cgetg(d+1,t_VEC);79for (i=1; i<=d; i++) gel(y,i) = randomi(p);80return y;81}8283GEN84random_FpC(long d, GEN p)85{86long i;87GEN y = cgetg(d+1,t_COL);88for (i=1; i<=d; i++) gel(y,i) = randomi(p);89return y;90}9192GEN93random_Flv(long n, ulong p)94{95GEN y = cgetg(n+1, t_VECSMALL);96long i;97for (i=1; i<=n; i++) uel(y,i) = random_Fl(p);98return y;99}100101/********************************************************************/102/** **/103/** ADD, SUB **/104/** **/105/********************************************************************/106GEN107FpC_add(GEN x, GEN y, GEN p)108{ pari_APPLY_type(t_COL, Fp_add(gel(x,i), gel(y,i), p)) }109110GEN111FpV_add(GEN x, GEN y, GEN p)112{ pari_APPLY_type(t_VEC, Fp_add(gel(x,i), gel(y,i), p)) }113114GEN115FpM_add(GEN x, GEN y, GEN p)116{ pari_APPLY_same(FpC_add(gel(x,i), gel(y,i), p)) }117118GEN119Flv_add(GEN x, GEN y, ulong p)120{121if (p==2)122pari_APPLY_ulong(x[i]^y[i])123else124pari_APPLY_ulong(Fl_add(x[i], y[i], p))125}126127void128Flv_add_inplace(GEN x, GEN y, ulong p)129{130long i, l = lg(x);131if (p==2)132for (i = 1; i < l; i++) x[i] ^= y[i];133else134for (i = 1; i < l; i++) x[i] = Fl_add(x[i], y[i], p);135}136137ulong138Flv_sum(GEN x, ulong p)139{140long i, l = lg(x);141ulong s = 0;142if (p==2)143for (i = 1; i < l; i++) s ^= x[i];144else145for (i = 1; i < l; i++) s = Fl_add(s, x[i], p);146return s;147}148149GEN150FpC_sub(GEN x, GEN y, GEN p)151{ pari_APPLY_type(t_COL, Fp_sub(gel(x,i), gel(y,i), p)) }152153GEN154FpV_sub(GEN x, GEN y, GEN p)155{ pari_APPLY_type(t_VEC, Fp_sub(gel(x,i), gel(y,i), p)) }156157GEN158FpM_sub(GEN x, GEN y, GEN p)159{ pari_APPLY_same(FpC_sub(gel(x,i), gel(y,i), p)) }160161GEN162Flv_sub(GEN x, GEN y, ulong p)163{ pari_APPLY_ulong(Fl_sub(uel(x,i), uel(y,i), p)) }164165void166Flv_sub_inplace(GEN x, GEN y, ulong p)167{168long i, l = lg(x);169for (i = 1; i < l; i++) x[i] = Fl_sub(x[i], y[i], p);170}171172GEN173Flm_Fl_add(GEN x, ulong y, ulong p)174{175long l = lg(x), i, j;176GEN z = cgetg(l,t_MAT);177178if (l==1) return z;179if (l != lgcols(x)) pari_err_OP( "+", x, utoi(y));180for (i=1; i<l; i++)181{182GEN zi = cgetg(l,t_VECSMALL), xi = gel(x,i);183gel(z,i) = zi;184for (j=1; j<l; j++) zi[j] = xi[j];185zi[i] = Fl_add(zi[i], y, p);186}187return z;188}189GEN190Flm_Fl_sub(GEN x, ulong y, ulong p) { return Flm_Fl_add(x, Fl_neg(y,p), p); }191192GEN193Flm_add(GEN x, GEN y, ulong p)194{ pari_APPLY_same(Flv_add(gel(x,i), gel(y,i), p)) }195196GEN197Flm_sub(GEN x, GEN y, ulong p)198{ pari_APPLY_same(Flv_sub(gel(x,i), gel(y,i), p)) }199200/********************************************************************/201/** **/202/** MULTIPLICATION **/203/** **/204/********************************************************************/205GEN206FpC_Fp_mul(GEN x, GEN y, GEN p)207{ pari_APPLY_type(t_COL, Fp_mul(gel(x,i), y, p)) }208209GEN210Flv_Fl_mul(GEN x, ulong y, ulong p)211{ pari_APPLY_ulong(Fl_mul(x[i], y, p)) }212213GEN214Flv_Fl_div(GEN x, ulong y, ulong p)215{216return Flv_Fl_mul(x, Fl_inv(y, p), p);217}218void219Flv_Fl_div_inplace(GEN x, ulong y, ulong p)220{221Flv_Fl_mul_inplace(x, Fl_inv(y, p), p);222}223GEN224FpM_Fp_mul(GEN X, GEN c, GEN p) {225long i, j, h, l = lg(X);226GEN A = cgetg(l, t_MAT);227if (l == 1) return A;228h = lgcols(X);229for (j=1; j<l; j++)230{231GEN a = cgetg(h, t_COL), x = gel(X, j);232for (i = 1; i < h; i++) gel(a,i) = Fp_mul(gel(x,i), c, p);233gel(A,j) = a;234}235return A;236}237238/* x *= y */239void240Flv_Fl_mul_part_inplace(GEN x, ulong y, ulong p, long l)241{242long i;243if (HIGHWORD(y | p))244for(i=1; i<=l; i++) x[i] = Fl_mul(x[i], y, p);245else246for(i=1; i<=l; i++) x[i] = (x[i] * y) % p;247}248void249Flv_Fl_mul_inplace(GEN x, ulong y, ulong p)250{ Flv_Fl_mul_part_inplace(x, y, p, lg(x)-1); }251252/* set y *= x */253void254Flm_Fl_mul_inplace(GEN y, ulong x, ulong p)255{256long i, j, m, l = lg(y);257if (l == 1) return;258m = lgcols(y);259if (HIGHWORD(x | p))260for(j=1; j<l; j++)261for(i=1; i<m; i++) ucoeff(y,i,j) = Fl_mul(ucoeff(y,i,j), x, p);262else263for(j=1; j<l; j++)264for(i=1; i<m; i++) ucoeff(y,i,j) = (ucoeff(y,i,j) * x) % p;265}266267/* return x * y */268static GEN269Flm_Fl_mul_pre_i(GEN y, ulong x, ulong p, ulong pi)270{271long i, j, m, l = lg(y);272GEN z = cgetg(l, t_MAT);273if (l == 1) return z;274m = lgcols(y);275for(j=1; j<l; j++) {276GEN c = cgetg(m, t_VECSMALL); gel(z,j) = c;277for(i=1; i<m; i++) uel(c,i) = Fl_mul_pre(ucoeff(y,i,j), x, p, pi);278}279return z;280}281282/* return x * y */283static GEN284Flm_Fl_mul_OK(GEN y, ulong x, ulong p)285{286long i, j, m, l = lg(y);287GEN z = cgetg(l, t_MAT);288if (l == 1) return z;289m = lgcols(y);290for(j=1; j<l; j++) {291GEN c = cgetg(m, t_VECSMALL); gel(z,j) = c;292for(i=1; i<m; i++) uel(c,i) = (ucoeff(y,i,j) * x) % p;293}294return z;295}296297/* return x * y */298GEN299Flm_Fl_mul_pre(GEN y, ulong x, ulong p, ulong pi)300{301if (HIGHWORD(p))302return Flm_Fl_mul_pre_i(y, x, p, pi);303else304return Flm_Fl_mul_OK(y, x, p);305}306307/* return x * y */308GEN309Flm_Fl_mul(GEN y, ulong x, ulong p)310{311if (HIGHWORD(p))312return Flm_Fl_mul_pre_i(y, x, p, get_Fl_red(p));313else314return Flm_Fl_mul_OK(y, x, p);315}316317GEN318Flv_neg(GEN x, ulong p)319{ pari_APPLY_ulong(Fl_neg(uel(x,i), p)) }320321void322Flv_neg_inplace(GEN v, ulong p)323{324long i;325for (i = 1; i < lg(v); ++i)326v[i] = Fl_neg(v[i], p);327}328329GEN330Flm_neg(GEN x, ulong p)331{ pari_APPLY_same(Flv_neg(gel(x,i), p)) }332333/* x[i,]*y. Assume lx > 1 and 0 < i < lgcols(x) */334static GEN335ZMrow_ZC_mul_i(GEN x, GEN y, long lx, long i)336{337GEN c = mulii(gcoeff(x,i,1), gel(y,1));338long k;339for (k = 2; k < lx; k++)340{341GEN t = mulii(gcoeff(x,i,k), gel(y,k));342if (signe(t)) c = addii(c, t);343}344return c;345}346347static long348zmrow_zc_mul(GEN x, GEN y, long lx, long i)349{350long k;351long c = coeff(x,i,1) * y[1];352for (k = 2; k < lx; k++)353c += coeff(x,i,k) * y[k];354return c;355}356357GEN358zm_zc_mul(GEN x, GEN y)359{360long lx = lg(x), l, i;361GEN z;362if (lx == 1) return cgetg(1, t_VECSMALL);363l = lg(gel(x,1));364z = cgetg(l,t_VECSMALL);365for (i=1; i<l; i++) z[i] = zmrow_zc_mul(x, y, lx, i);366return z;367}368369GEN370zm_mul(GEN x, GEN y)371{372long i,j,lx=lg(x), ly=lg(y);373GEN z;374if (ly==1) return cgetg(1,t_MAT);375z = cgetg(ly,t_MAT);376if (lx==1)377{378for (i=1; i<ly; i++) gel(z,i) = cgetg(1,t_VECSMALL);379return z;380}381for (j=1; j<ly; j++)382gel(z,j) = zm_zc_mul(x, gel(y,j));383return z;384}385386static ulong387Flmrow_Flc_mul_SMALL(GEN x, GEN y, ulong p, long lx, long i)388{389ulong c = ucoeff(x,i,1) * uel(y,1);390long k;391for (k = 2; k < lx; k++) {392c += ucoeff(x,i,k) * uel(y,k);393if (c & HIGHBIT) c %= p;394}395return c % p;396}397398static ulong399Flmrow_Flc_mul_i(GEN x, GEN y, ulong p, ulong pi, long lx, long i)400{401ulong l0, l1, v1, h0, h1;402long k = 1;403LOCAL_OVERFLOW;404LOCAL_HIREMAINDER;405l1 = mulll(ucoeff(x,i,k), uel(y,k)); h1 = hiremainder; v1 = 0;406while (++k < lx) {407l0 = mulll(ucoeff(x,i,k), uel(y,k)); h0 = hiremainder;408l1 = addll(l0, l1); h1 = addllx(h0, h1); v1 += overflow;409}410if (v1 == 0) return remll_pre(h1, l1, p, pi);411else return remlll_pre(v1, h1, l1, p, pi);412}413414static GEN415Flm_Flc_mul_i_2(GEN x, GEN y, long lx, long l)416{417long i,j;418GEN z = NULL;419420for (j=1; j<lx; j++)421{422if (!y[j]) continue;423if (!z) z = Flv_copy(gel(x,j));424else for (i = 1; i < l; i++) z[i] ^= coeff(x,i,j);425}426if (!z) z = zero_zv(l-1);427return z;428}429430static GEN431FpM_FpC_mul_i(GEN x, GEN y, long lx, long l, GEN p)432{433GEN z = cgetg(l,t_COL);434long i;435for (i = 1; i < l; i++)436{437pari_sp av = avma;438GEN c = ZMrow_ZC_mul_i(x, y, lx, i);439gel(z,i) = gerepileuptoint(av, modii(c,p));440}441return z;442}443444static void445__Flm_Flc_mul_i_SMALL(GEN z, GEN x, GEN y, long lx, long l, ulong p)446{447long i;448for (i = 1; i < l; i++) z[i] = Flmrow_Flc_mul_SMALL(x, y, p, lx, i);449}450static GEN451Flm_Flc_mul_i_SMALL(GEN x, GEN y, long lx, long l, ulong p)452{453GEN z = cgetg(l,t_VECSMALL);454__Flm_Flc_mul_i_SMALL(z, x, y, lx, l, p);455return z;456}457458static void459__Flm_Flc_mul_i(GEN z, GEN x, GEN y, long lx, long l, ulong p, ulong pi)460{461long i;462for (i = 1; i < l; i++) z[i] = Flmrow_Flc_mul_i(x, y, p, pi, lx, i);463}464static GEN465Flm_Flc_mul_i(GEN x, GEN y, long lx, long l, ulong p, ulong pi)466{467GEN z = cgetg(l,t_VECSMALL);468__Flm_Flc_mul_i(z, x, y, lx, l, p, pi);469return z;470}471472GEN473FpM_mul(GEN x, GEN y, GEN p)474{475long lx=lg(x), ly=lg(y);476GEN z;477pari_sp av = avma;478if (ly==1) return cgetg(1,t_MAT);479if (lx==1) return zeromat(0, ly-1);480if (lgefint(p) == 3)481{482ulong pp = uel(p,2);483if (pp == 2)484{485x = ZM_to_F2m(x);486y = ZM_to_F2m(y);487z = F2m_to_ZM(F2m_mul(x,y));488}489else490{491x = ZM_to_Flm(x, pp);492y = ZM_to_Flm(y, pp);493z = Flm_to_ZM(Flm_mul(x,y, pp));494}495} else496z = FpM_red(ZM_mul(x, y), p);497return gerepileupto(av, z);498}499500static GEN501Flm_mul_classical(GEN x, GEN y, long l, long lx, long ly, ulong p, ulong pi)502{503long j;504GEN z = cgetg(ly, t_MAT);505if (SMALL_ULONG(p))506for (j=1; j<ly; j++)507gel(z,j) = Flm_Flc_mul_i_SMALL(x, gel(y,j), lx, l, p);508else509for (j=1; j<ly; j++)510gel(z,j) = Flm_Flc_mul_i(x, gel(y,j), lx, l, p, pi);511return z;512}513514/*515Return A[ma+1..ma+da, na+1..na+ea] - B[mb+1..mb+db, nb+1..nb+eb]516as an (m x n)-matrix, padding the input with zeroes as necessary.517*/518static void519add_slices_ip(long m, long n,520GEN A, long ma, long da, long na, long ea,521GEN B, long mb, long db, long nb, long eb,522GEN M, long dy, long dx, ulong p)523{524long min_d = minss(da, db), min_e = minss(ea, eb), i, j;525GEN C;526527for (j = 1; j <= min_e; j++) {528C = gel(M, j + dx) + dy;529for (i = 1; i <= min_d; i++)530uel(C, i) = Fl_add(ucoeff(A, ma + i, na + j),531ucoeff(B, mb + i, nb + j), p);532for (; i <= da; i++)533uel(C, i) = ucoeff(A, ma + i, na + j);534for (; i <= db; i++)535uel(C, i) = ucoeff(B, mb + i, nb + j);536for (; i <= m; i++)537uel(C, i) = 0;538}539for (; j <= ea; j++) {540C = gel(M, j + dx) + dy;541for (i = 1; i <= da; i++)542uel(C, i) = ucoeff(A, ma + i, na + j);543for (; i <= m; i++)544uel(C, i) = 0;545}546for (; j <= eb; j++) {547C = gel(M, j + dx) + dy;548for (i = 1; i <= db; i++)549uel(C, i) = ucoeff(B, mb + i, nb + j);550for (; i <= m; i++)551uel(C, i) = 0;552}553for (; j <= n; j++) {554C = gel(M, j + dx) + dy;555for (i = 1; i <= m; i++)556uel(C, i) = 0;557}558}559560static GEN561add_slices(long m, long n,562GEN A, long ma, long da, long na, long ea,563GEN B, long mb, long db, long nb, long eb, ulong p)564{565GEN M;566long j;567M = cgetg(n + 1, t_MAT);568for (j = 1; j <= n; j++)569gel(M, j) = cgetg(m + 1, t_VECSMALL);570add_slices_ip(m, n, A, ma, da, na, ea, B, mb, db, nb, eb, M, 0, 0, p);571return M;572}573574/*575Return A[ma+1..ma+da, na+1..na+ea] - B[mb+1..mb+db, nb+1..nb+eb]576as an (m x n)-matrix, padding the input with zeroes as necessary.577*/578static GEN579subtract_slices(long m, long n,580GEN A, long ma, long da, long na, long ea,581GEN B, long mb, long db, long nb, long eb, ulong p)582{583long min_d = minss(da, db), min_e = minss(ea, eb), i, j;584GEN M = cgetg(n + 1, t_MAT), C;585586for (j = 1; j <= min_e; j++) {587gel(M, j) = C = cgetg(m + 1, t_VECSMALL);588for (i = 1; i <= min_d; i++)589uel(C, i) = Fl_sub(ucoeff(A, ma + i, na + j),590ucoeff(B, mb + i, nb + j), p);591for (; i <= da; i++)592uel(C, i) = ucoeff(A, ma + i, na + j);593for (; i <= db; i++)594uel(C, i) = Fl_neg(ucoeff(B, mb + i, nb + j), p);595for (; i <= m; i++)596uel(C, i) = 0;597}598for (; j <= ea; j++) {599gel(M, j) = C = cgetg(m + 1, t_VECSMALL);600for (i = 1; i <= da; i++)601uel(C, i) = ucoeff(A, ma + i, na + j);602for (; i <= m; i++)603uel(C, i) = 0;604}605for (; j <= eb; j++) {606gel(M, j) = C = cgetg(m + 1, t_VECSMALL);607for (i = 1; i <= db; i++)608uel(C, i) = Fl_neg(ucoeff(B, mb + i, nb + j), p);609for (; i <= m; i++)610uel(C, i) = 0;611}612for (; j <= n; j++)613gel(M, j) = zero_Flv(m);614return M;615}616617static GEN Flm_mul_i(GEN x, GEN y, long l, long lx, long ly, ulong p, ulong pi);618619/* Strassen-Winograd matrix product A (m x n) * B (n x p) */620static GEN621Flm_mul_sw(GEN A, GEN B, long m, long n, long p, ulong l, ulong li)622{623pari_sp av;624long m1 = (m + 1)/2, m2 = m/2,625n1 = (n + 1)/2, n2 = n/2,626p1 = (p + 1)/2, p2 = p/2;627GEN A11, A12, A22, B11, B21, B22,628S1, S2, S3, S4, T1, T2, T3, T4,629M1, M2, M3, M4, M5, M6, M7,630V1, V2, V3, C;631long j;632C = cgetg(p + 1, t_MAT);633for (j = 1; j <= p; j++)634gel(C, j) = cgetg(m + 1, t_VECSMALL);635av = avma;636T2 = subtract_slices(n1, p2, B, 0, n1, p1, p2, B, n1, n2, p1, p2, l);637S1 = subtract_slices(m2, n1, A, m1, m2, 0, n1, A, 0, m2, 0, n1, l);638M2 = Flm_mul_i(S1, T2, m2 + 1, n1 + 1, p2 + 1, l, li);639if (gc_needed(av, 1))640gerepileall(av, 2, &T2, &M2); /* destroy S1 */641T3 = subtract_slices(n1, p1, T2, 0, n1, 0, p2, B, 0, n1, 0, p1, l);642if (gc_needed(av, 1))643gerepileall(av, 2, &M2, &T3); /* destroy T2 */644S2 = add_slices(m2, n1, A, m1, m2, 0, n1, A, m1, m2, n1, n2, l);645T1 = subtract_slices(n1, p1, B, 0, n1, p1, p2, B, 0, n1, 0, p2, l);646M3 = Flm_mul_i(S2, T1, m2 + 1, n1 + 1, p2 + 1, l, li);647if (gc_needed(av, 1))648gerepileall(av, 4, &M2, &T3, &S2, &M3); /* destroy T1 */649S3 = subtract_slices(m1, n1, S2, 0, m2, 0, n1, A, 0, m1, 0, n1, l);650if (gc_needed(av, 1))651gerepileall(av, 4, &M2, &T3, &M3, &S3); /* destroy S2 */652A11 = matslice(A, 1, m1, 1, n1);653B11 = matslice(B, 1, n1, 1, p1);654M1 = Flm_mul_i(A11, B11, m1 + 1, n1 + 1, p1 + 1, l, li);655if (gc_needed(av, 1))656gerepileall(av, 5, &M2, &T3, &M3, &S3, &M1); /* destroy A11, B11 */657A12 = matslice(A, 1, m1, n1 + 1, n);658B21 = matslice(B, n1 + 1, n, 1, p1);659M4 = Flm_mul_i(A12, B21, m1 + 1, n2 + 1, p1 + 1, l, li);660if (gc_needed(av, 1))661gerepileall(av, 6, &M2, &T3, &M3, &S3, &M1, &M4); /* destroy A12, B21 */662add_slices_ip(m1, p1, M1, 0, m1, 0, p1, M4, 0, m1, 0, p1, C, 0, 0, l);663if (gc_needed(av, 1))664gerepileall(av, 5, &M2, &T3, &M3, &S3, &M1); /* destroy M4 */665M5 = Flm_mul_i(S3, T3, m1 + 1, n1 + 1, p1 + 1, l, li);666S4 = subtract_slices(m1, n2, A, 0, m1, n1, n2, S3, 0, m1, 0, n2, l);667if (gc_needed(av, 1))668gerepileall(av, 6, &M2, &T3, &M3, &M1, &M5, &S4); /* destroy S3 */669T4 = add_slices(n2, p1, B, n1, n2, 0, p1, T3, 0, n2, 0, p1, l);670if (gc_needed(av, 1))671gerepileall(av, 6, &M2, &M3, &M1, &M5, &S4, &T4); /* destroy T3 */672V1 = subtract_slices(m1, p1, M1, 0, m1, 0, p1, M5, 0, m1, 0, p1, l);673if (gc_needed(av, 1))674gerepileall(av, 5, &M2, &M3, &S4, &T4, &V1); /* destroy M1, M5 */675B22 = matslice(B, n1 + 1, n, p1 + 1, p);676M6 = Flm_mul_i(S4, B22, m1 + 1, n2 + 1, p2 + 1, l, li);677if (gc_needed(av, 1))678gerepileall(av, 5, &M2, &M3, &T4, &V1, &M6); /* destroy S4, B22 */679A22 = matslice(A, m1 + 1, m, n1 + 1, n);680M7 = Flm_mul_i(A22, T4, m2 + 1, n2 + 1, p1 + 1, l, li);681if (gc_needed(av, 1))682gerepileall(av, 5, &M2, &M3, &V1, &M6, &M7); /* destroy A22, T4 */683V3 = add_slices(m1, p2, V1, 0, m1, 0, p2, M3, 0, m2, 0, p2, l);684add_slices_ip(m1, p2, V3, 0, m1, 0, p2, M6, 0, m1, 0, p2, C, 0, p1, l);685if (gc_needed(av, 1))686gerepileall(av, 4, &M2, &M3, &V1, &M7); /* destroy V3, M6 */687V2 = add_slices(m2, p1, V1, 0, m2, 0, p1, M2, 0, m2, 0, p2, l);688if (gc_needed(av, 1))689gerepileall(av, 3, &M3, &M7, &V2); /* destroy V1, M2 */690add_slices_ip(m2, p1, V2, 0, m2, 0, p1, M7, 0, m2, 0, p1, C, m1, 0, l);691add_slices_ip(m2, p2, V2, 0, m2, 0, p2, M3, 0, m2, 0, p2, C, m1, p1, l);692set_avma(av); return C;693}694695/* Strassen-Winograd used for dim >= ZM_sw_bound */696static GEN697Flm_mul_i(GEN x, GEN y, long l, long lx, long ly, ulong p, ulong pi)698{699ulong e = expu(p);700#ifdef LONG_IS_64BIT701long ZM_sw_bound = e <= 29 ? 140: e <=62 ? 40: 70;702#else703long ZM_sw_bound = e <= 12 ? 230: e <=14 ? 170 : e <=17 ? 110: 120;704#endif705if (l <= ZM_sw_bound || lx <= ZM_sw_bound || ly <= ZM_sw_bound)706return Flm_mul_classical(x, y, l, lx, ly, p, pi);707else708return Flm_mul_sw(x, y, l - 1, lx - 1, ly - 1, p, pi);709}710711GEN712Flm_mul_pre(GEN x, GEN y, ulong p, ulong pi)713{714long lx=lg(x), ly=lg(y);715if (ly==1) return cgetg(1,t_MAT);716if (lx==1) return zero_Flm(0, ly-1);717return Flm_mul_i(x, y, lgcols(x), lx, ly, p, pi);718}719720GEN721Flm_mul(GEN x, GEN y, ulong p)722{723long lx=lg(x), ly=lg(y);724if (ly==1) return cgetg(1,t_MAT);725if (lx==1) return zero_Flm(0, ly-1);726return Flm_mul_i(x, y, lgcols(x), lx, ly, p, get_Fl_red(p));727}728729struct _Flm730{731ulong p;732long n;733};734735static GEN736_Flm_mul(void *E , GEN x, GEN y)737{ return Flm_mul(x,y,((struct _Flm*)E)->p); }738static GEN739_Flm_sqr(void *E, GEN x)740{ return Flm_mul(x,x,((struct _Flm*)E)->p); }741static GEN742_Flm_one(void *E)743{ return matid_Flm(((struct _Flm*)E)->n); }744GEN745Flm_powu(GEN x, ulong n, ulong p)746{747struct _Flm d;748if (!n) return matid(lg(x)-1);749d.p = p;750return gen_powu(x, n, (void*)&d, &_Flm_sqr, &_Flm_mul);751}752GEN753Flm_powers(GEN x, ulong n, ulong p)754{755struct _Flm d;756d.p = p;757d.n = lg(x)-1;758return gen_powers(x, n, 1, (void*)&d,759&_Flm_sqr, &_Flm_mul, &_Flm_one);760}761762static GEN763_FpM_mul(void *p , GEN x, GEN y)764{ return FpM_mul(x,y,(GEN)p); }765static GEN766_FpM_sqr(void *p, GEN x)767{ return FpM_mul(x,x,(GEN)p); }768GEN769FpM_powu(GEN x, ulong n, GEN p)770{771if (!n) return matid(lg(x)-1);772if (lgefint(p) == 3)773{774pari_sp av = avma;775ulong pp = uel(p,2);776GEN z;777if (pp == 2)778z = F2m_to_ZM(F2m_powu(ZM_to_F2m(x),n));779else780z = Flm_to_ZM(Flm_powu(ZM_to_Flm(x, pp), n, pp));781return gerepileupto(av, z);782}783return gen_powu(x, n, (void*)p, &_FpM_sqr, &_FpM_mul);784}785786/*Multiple a column vector by a line vector to make a matrix*/787GEN788Flc_Flv_mul(GEN x, GEN y, ulong p)789{790long i,j, lx=lg(x), ly=lg(y);791GEN z;792if (ly==1) return cgetg(1,t_MAT);793z = cgetg(ly, t_MAT);794for (j=1; j < ly; j++)795{796GEN zj = cgetg(lx,t_VECSMALL);797for (i=1; i<lx; i++)798uel(zj,i) = Fl_mul(uel(x,i), uel(y,j), p);799gel(z,j) = zj;800}801return z;802}803804/*Multiple a column vector by a line vector to make a matrix*/805GEN806FpC_FpV_mul(GEN x, GEN y, GEN p)807{808long i,j, lx=lg(x), ly=lg(y);809GEN z;810if (ly==1) return cgetg(1,t_MAT);811z = cgetg(ly,t_MAT);812for (j=1; j < ly; j++)813{814GEN zj = cgetg(lx,t_COL);815for (i=1; i<lx; i++) gel(zj,i) = Fp_mul(gel(x,i),gel(y,j), p);816gel(z, j) = zj;817}818return z;819}820821/* Multiply a line vector by a column and return a scalar (t_INT) */822GEN823FpV_dotproduct(GEN x, GEN y, GEN p)824{825long i, lx = lg(x);826pari_sp av;827GEN c;828if (lx == 1) return gen_0;829av = avma; c = mulii(gel(x,1),gel(y,1));830for (i=2; i<lx; i++) c = addii(c, mulii(gel(x,i),gel(y,i)));831return gerepileuptoint(av, modii(c,p));832}833GEN834FpV_dotsquare(GEN x, GEN p)835{836long i, lx = lg(x);837pari_sp av;838GEN c;839if (lx == 1) return gen_0;840av = avma; c = sqri(gel(x,1));841for (i=2; i<lx; i++) c = addii(c, sqri(gel(x,i)));842return gerepileuptoint(av, modii(c,p));843}844845INLINE ulong846Flv_dotproductspec_SMALL(GEN x, GEN y, ulong p, long lx)847{848ulong c = uel(x,0) * uel(y,0);849long k;850for (k = 1; k < lx; k++) {851c += uel(x,k) * uel(y,k);852if (c & HIGHBIT) c %= p;853}854return c % p;855}856857INLINE ulong858Flv_dotproductspec_i(GEN x, GEN y, ulong p, ulong pi, long lx)859{860ulong l0, l1, v1, h0, h1;861long i = 0;862LOCAL_OVERFLOW;863LOCAL_HIREMAINDER;864l1 = mulll(uel(x,i), uel(y,i)); h1 = hiremainder; v1 = 0;865while (++i < lx) {866l0 = mulll(uel(x,i), uel(y,i)); h0 = hiremainder;867l1 = addll(l0, l1); h1 = addllx(h0, h1); v1 += overflow;868}869if (v1 == 0) return remll_pre(h1, l1, p, pi);870else return remlll_pre(v1, h1, l1, p, pi);871}872873ulong874Flv_dotproduct(GEN x, GEN y, ulong p)875{876long lx = lg(x)-1;877if (lx == 0) return 0;878if (SMALL_ULONG(p))879return Flv_dotproductspec_SMALL(x+1, y+1, p, lx);880else881return Flv_dotproductspec_i(x+1, y+1, p, get_Fl_red(p), lx);882}883884ulong885Flv_dotproduct_pre(GEN x, GEN y, ulong p, ulong pi)886{887long lx = lg(x)-1;888if (lx == 0) return 0;889if (SMALL_ULONG(p))890return Flv_dotproductspec_SMALL(x+1, y+1, p, lx);891else892return Flv_dotproductspec_i(x+1, y+1, p, pi, lx);893}894895ulong896Flx_dotproduct(GEN x, GEN y, ulong p)897{898long lx = minss(lgpol(x), lgpol(y));899if (lx == 0) return 0;900if (SMALL_ULONG(p))901return Flv_dotproductspec_SMALL(x+2, y+2, p, lx);902else903return Flv_dotproductspec_i(x+2, y+2, p, get_Fl_red(p), lx);904}905906GEN907FpM_FpC_mul(GEN x, GEN y, GEN p)908{909long lx = lg(x);910return lx==1? cgetg(1,t_COL): FpM_FpC_mul_i(x, y, lx, lgcols(x), p);911}912GEN913Flm_Flc_mul(GEN x, GEN y, ulong p)914{915long l, lx = lg(x);916if (lx==1) return cgetg(1,t_VECSMALL);917l = lgcols(x);918if (p==2)919return Flm_Flc_mul_i_2(x, y, lx, l);920else if (SMALL_ULONG(p))921return Flm_Flc_mul_i_SMALL(x, y, lx, l, p);922else923return Flm_Flc_mul_i(x, y, lx, l, p, get_Fl_red(p));924}925926GEN927Flm_Flc_mul_pre(GEN x, GEN y, ulong p, ulong pi)928{929long l, lx = lg(x);930GEN z;931if (lx==1) return cgetg(1,t_VECSMALL);932l = lgcols(x);933z = cgetg(l, t_VECSMALL);934if (SMALL_ULONG(p))935__Flm_Flc_mul_i_SMALL(z, x, y, lx, l, p);936else937__Flm_Flc_mul_i(z, x, y, lx, l, p, pi);938return z;939}940941GEN942Flm_Flc_mul_pre_Flx(GEN x, GEN y, ulong p, ulong pi, long sv)943{944long l, lx = lg(x);945GEN z;946if (lx==1) return pol0_Flx(sv);947l = lgcols(x);948z = cgetg(l + 1, t_VECSMALL); z[1] = sv;949if (SMALL_ULONG(p))950__Flm_Flc_mul_i_SMALL(z + 1, x, y, lx, l, p);951else952__Flm_Flc_mul_i(z + 1, x, y, lx, l, p, pi);953return Flx_renormalize(z, l + 1);954}955956/* RgV_to_RgX(FpM_FpC_mul(x,y,p), v), p != NULL, memory clean */957GEN958FpM_FpC_mul_FpX(GEN x, GEN y, GEN p, long v)959{960long i, l, lx = lg(x);961GEN z;962if (lx==1) return pol_0(v);963l = lgcols(x);964z = new_chunk(l+1);965for (i=l-1; i; i--)966{967pari_sp av = avma;968GEN p1 = ZMrow_ZC_mul_i(x,y,lx,i);969p1 = modii(p1, p);970if (signe(p1))971{972if (i != l-1) stackdummy((pari_sp)(z + l+1), (pari_sp)(z + i+2));973gel(z,i+1) = gerepileuptoint(av, p1);974break;975}976set_avma(av);977}978if (!i) { set_avma((pari_sp)(z + l+1)); return pol_0(v); }979z[0] = evaltyp(t_POL) | evallg(i+2);980z[1] = evalsigne(1) | evalvarn(v);981for (; i; i--)982{983pari_sp av = avma;984GEN p1 = ZMrow_ZC_mul_i(x,y,lx,i);985gel(z,i+1) = gerepileuptoint(av, modii(p1,p));986}987return z;988}989990/********************************************************************/991/** **/992/** TRANSPOSITION **/993/** **/994/********************************************************************/995996/* == zm_transpose */997GEN998Flm_transpose(GEN x)999{1000long i, dx, lx = lg(x);1001GEN y;1002if (lx == 1) return cgetg(1,t_MAT);1003dx = lgcols(x); y = cgetg(dx,t_MAT);1004for (i=1; i<dx; i++) gel(y,i) = Flm_row(x,i);1005return y;1006}10071008/********************************************************************/1009/** **/1010/** SCALAR MATRICES **/1011/** **/1012/********************************************************************/10131014GEN1015gen_matid(long n, void *E, const struct bb_field *S)1016{1017GEN y = cgetg(n+1,t_MAT), _0, _1;1018long i;1019if (n < 0) pari_err_DOMAIN("gen_matid", "dimension","<",gen_0,stoi(n));1020_0 = S->s(E,0);1021_1 = S->s(E,1);1022for (i=1; i<=n; i++)1023{1024GEN z = const_col(n, _0); gel(z,i) = _1;1025gel(y, i) = z;1026}1027return y;1028}10291030GEN1031matid_F2xqM(long n, GEN T)1032{1033void *E;1034const struct bb_field *S = get_F2xq_field(&E, T);1035return gen_matid(n, E, S);1036}1037GEN1038matid_FlxqM(long n, GEN T, ulong p)1039{1040void *E;1041const struct bb_field *S = get_Flxq_field(&E, T, p);1042return gen_matid(n, E, S);1043}10441045GEN1046matid_Flm(long n)1047{1048GEN y = cgetg(n+1,t_MAT);1049long i;1050if (n < 0) pari_err_DOMAIN("matid_Flm", "dimension","<",gen_0,stoi(n));1051for (i=1; i<=n; i++) { gel(y,i) = zero_zv(n); ucoeff(y, i,i) = 1; }1052return y;1053}10541055GEN1056scalar_Flm(long s, long n)1057{1058long i;1059GEN y = cgetg(n+1,t_MAT);1060for (i=1; i<=n; i++) { gel(y,i) = zero_Flv(n); coeff(y, i,i) = s; }1061return y;1062}10631064/********************************************************************/1065/** **/1066/** CONVERSIONS **/1067/** **/1068/********************************************************************/1069GEN1070ZV_to_Flv(GEN x, ulong p)1071{ pari_APPLY_ulong(umodiu(gel(x,i), p)) }10721073GEN1074ZM_to_Flm(GEN x, ulong p)1075{ pari_APPLY_same(ZV_to_Flv(gel(x,i), p)) }10761077GEN1078ZMV_to_FlmV(GEN x, ulong m)1079{ pari_APPLY_type(t_VEC,ZM_to_Flm(gel(x,i), m)) }10801081/* TO INTMOD */1082static GEN1083to_intmod(GEN x, GEN p) { retmkintmod(modii(x, p), p); }1084static GEN1085Fl_to_intmod(ulong x, GEN p) { retmkintmod(utoi(x), p); }10861087GEN1088Fp_to_mod(GEN z, GEN p)1089{1090retmkintmod(modii(z, p), icopy(p));1091}10921093GEN1094FpX_to_mod_raw(GEN z, GEN p)1095{1096long i,l = lg(z);1097GEN x;1098if (l == 2)1099{1100x = cgetg(3,t_POL); x[1] = z[1];1101gel(x,2) = mkintmod(gen_0,p); return x;1102}1103x = cgetg(l,t_POL);1104for (i=2; i<l; i++) gel(x,i) = to_intmod(gel(z,i), p);1105x[1] = z[1]; return normalizepol_lg(x,l);1106}11071108/* z in Z[X], return z * Mod(1,p), normalized*/1109GEN1110FpX_to_mod(GEN z, GEN p)1111{1112long i,l = lg(z);1113GEN x;1114if (l == 2)1115{1116x = cgetg(3,t_POL); x[1] = z[1];1117gel(x,2) = mkintmod(gen_0,icopy(p)); return x;1118}1119x = cgetg(l,t_POL);1120if (l >2) p = icopy(p);1121for (i=2; i<l; i++) gel(x,i) = to_intmod(gel(z,i), p);1122x[1] = z[1]; return normalizepol_lg(x,l);1123}11241125GEN1126FpXC_to_mod(GEN z, GEN p)1127{1128long i,l = lg(z);1129GEN x = cgetg(l, t_COL);1130p = icopy(p);1131for (i=1; i<l; i++)1132gel(x,i) = FpX_to_mod_raw(gel(z,i), p);1133return x;1134}11351136static GEN1137FpXC_to_mod_raw(GEN x, GEN p)1138{ pari_APPLY_type(t_COL, FpX_to_mod_raw(gel(x,i), p)) }11391140GEN1141FpXM_to_mod(GEN z, GEN p)1142{1143long i,l = lg(z);1144GEN x = cgetg(l, t_MAT);1145p = icopy(p);1146for (i=1; i<l; i++)1147gel(x,i) = FpXC_to_mod_raw(gel(z,i), p);1148return x;1149}11501151/* z in Z^n, return z * Mod(1,p), normalized*/1152GEN1153FpV_to_mod(GEN z, GEN p)1154{1155long i,l = lg(z);1156GEN x = cgetg(l, t_VEC);1157if (l == 1) return x;1158p = icopy(p);1159for (i=1; i<l; i++) gel(x,i) = to_intmod(gel(z,i), p);1160return x;1161}1162/* z in Z^n, return z * Mod(1,p), normalized*/1163GEN1164FpC_to_mod(GEN z, GEN p)1165{1166long i, l = lg(z);1167GEN x = cgetg(l, t_COL);1168if (l == 1) return x;1169p = icopy(p);1170for (i=1; i<l; i++) gel(x,i) = to_intmod(gel(z,i), p);1171return x;1172}1173/* z in Mat m,n(Z), return z * Mod(1,p), normalized*/1174GEN1175FpM_to_mod(GEN z, GEN p)1176{1177long i, j, m, l = lg(z);1178GEN x = cgetg(l,t_MAT), y, zi;1179if (l == 1) return x;1180m = lgcols(z);1181p = icopy(p);1182for (i=1; i<l; i++)1183{1184gel(x,i) = cgetg(m,t_COL);1185y = gel(x,i); zi= gel(z,i);1186for (j=1; j<m; j++) gel(y,j) = to_intmod(gel(zi,j), p);1187}1188return x;1189}1190GEN1191Flc_to_mod(GEN z, ulong pp)1192{1193long i, l = lg(z);1194GEN p, x = cgetg(l, t_COL);1195if (l == 1) return x;1196p = utoipos(pp);1197for (i=1; i<l; i++) gel(x,i) = Fl_to_intmod(z[i], p);1198return x;1199}1200GEN1201Flm_to_mod(GEN z, ulong pp)1202{1203long i, j, m, l = lg(z);1204GEN p, x = cgetg(l,t_MAT), y, zi;1205if (l == 1) return x;1206m = lgcols(z);1207p = utoipos(pp);1208for (i=1; i<l; i++)1209{1210gel(x,i) = cgetg(m,t_COL);1211y = gel(x,i); zi= gel(z,i);1212for (j=1; j<m; j++) gel(y,j) = Fl_to_intmod(zi[j], p);1213}1214return x;1215}12161217GEN1218FpVV_to_mod(GEN z, GEN p)1219{1220long i, j, m, l = lg(z);1221GEN x = cgetg(l,t_VEC), y, zi;1222if (l == 1) return x;1223m = lgcols(z);1224p = icopy(p);1225for (i=1; i<l; i++)1226{1227gel(x,i) = cgetg(m,t_VEC);1228y = gel(x,i); zi= gel(z,i);1229for (j=1; j<m; j++) gel(y,j) = to_intmod(gel(zi,j), p);1230}1231return x;1232}12331234/* z in Z^n, return z * Mod(1,p), normalized*/1235GEN1236FpXQC_to_mod(GEN z, GEN T, GEN p)1237{1238long i,l = lg(z);1239GEN x = cgetg(l, t_COL);1240if (l == 1) return x;1241p = icopy(p);1242T = FpX_to_mod_raw(T, p);1243for (i=1; i<l; i++)1244gel(x,i) = mkpolmod(FpX_to_mod_raw(gel(z,i), p), T);1245return x;1246}12471248static GEN1249Fq_to_mod_raw(GEN x, GEN T, GEN p)1250{1251GEN z = typ(x)==t_INT ? mkintmod(modii(x, p), p): FpX_to_mod_raw(x, p);1252return mkpolmod(z, T);1253}12541255/* z in Z^n, return z * Mod(1,p), normalized*/1256GEN1257FqC_to_mod(GEN z, GEN T, GEN p)1258{1259GEN x;1260long i,l = lg(z);1261if (!T) return FpC_to_mod(z, p);1262x = cgetg(l, t_COL);1263if (l == 1) return x;1264p = icopy(p);1265T = FpX_to_mod_raw(T, p);1266for (i=1; i<l; i++)1267gel(x,i) = Fq_to_mod_raw(gel(z, i), T, p);1268return x;1269}12701271/* z in Z^n, return z * Mod(1,p), normalized*/1272static GEN1273FqC_to_mod_raw(GEN x, GEN T, GEN p)1274{ pari_APPLY_type(t_COL, Fq_to_mod_raw(gel(x,i), T, p)) }12751276/* z in Z^n, return z * Mod(1,p), normalized*/1277GEN1278FqM_to_mod(GEN z, GEN T, GEN p)1279{1280GEN x;1281long i,l = lg(z);1282if (!T) return FpM_to_mod(z, p);1283x = cgetg(l, t_MAT);1284if (l == 1) return x;1285p = icopy(p);1286T = FpX_to_mod_raw(T, p);1287for (i=1; i<l; i++)1288gel(x,i) = FqC_to_mod_raw(gel(z, i), T, p);1289return x;1290}12911292/********************************************************************/1293/* */1294/* Blackbox linear algebra */1295/* */1296/********************************************************************/12971298/* A sparse column (zCs) v is a t_COL with two components C and E which are1299* t_VECSMALL of the same length, representing sum_i E[i]*e_{C[i]}, where1300* (e_j) is the canonical basis.1301* A sparse matrix (zMs) is a t_VEC of zCs */13021303/* FpCs and FpMs are identical but E[i] is interpreted as a _signed_ C long1304* integer representing an element of Fp. This is important since p can be1305* large and p+E[i] would not fit in a C long. */13061307/* RgCs and RgMs are similar, except that the type of the component is1308* unspecified. Functions handling RgCs/RgMs must be independent of the type1309* of E. */13101311/* Most functions take an argument nbrow which is the number of lines of the1312* column/matrix, which cannot be derived from the data. */13131314GEN1315zCs_to_ZC(GEN R, long nbrow)1316{1317GEN C = gel(R,1), E = gel(R,2), c = zerocol(nbrow);1318long j, l = lg(C);1319for (j = 1; j < l; ++j) gel(c, C[j]) = stoi(E[j]);1320return c;1321}13221323GEN1324zMs_to_ZM(GEN x, long nbrow)1325{ pari_APPLY_same(zCs_to_ZC(gel(x, i), nbrow)) }13261327/* Solve equation f(X) = B (mod p) where B is a FpV, and f is an endomorphism.1328* Return either a solution as a t_COL, or a kernel vector as a t_VEC. */1329GEN1330gen_FpM_Wiedemann(void *E, GEN (*f)(void*, GEN), GEN B, GEN p)1331{1332pari_sp ltop = avma;1333long col = 0, n = lg(B)-1, m = 2*n+1;1334if (ZV_equal0(B)) return zerocol(n);1335while (++col <= n)1336{1337pari_sp btop = avma, av;1338long i, lQ;1339GEN V, Q, M, W = B;1340GEN b = cgetg(m+2, t_POL);1341b[1] = evalsigne(1)|evalvarn(0);1342gel(b, 2) = gel(W, col);1343for (i = 3; i<m+2; i++)1344gel(b, i) = cgeti(lgefint(p));1345av = avma;1346for (i = 3; i<m+2; i++)1347{1348W = f(E, W);1349affii(gel(W, col),gel(b, i));1350if (gc_needed(av,1))1351{1352if (DEBUGMEM>1) pari_warn(warnmem,"Wiedemann: first loop, %ld",i);1353W = gerepileupto(av, W);1354}1355}1356b = FpX_renormalize(b, m+2);1357if (lgpol(b)==0) {set_avma(btop); continue; }1358M = FpX_halfgcd(b, pol_xn(m, 0), p);1359Q = FpX_neg(FpX_normalize(gcoeff(M, 2, 1),p),p);1360W = B; lQ =lg(Q);1361if (DEBUGLEVEL) err_printf("Wiedemann: deg. minpoly: %ld\n",lQ-3);1362V = FpC_Fp_mul(W, gel(Q, lQ-2), p);1363av = avma;1364for (i = lQ-3; i > 1; i--)1365{1366W = f(E, W);1367V = ZC_lincomb(gen_1, gel(Q,i), V, W);1368if (gc_needed(av,1))1369{1370if (DEBUGMEM>1) pari_warn(warnmem,"Wiedemann: second loop, %ld",i);1371gerepileall(av, 2, &V, &W);1372}1373}1374V = FpC_red(V, p);1375W = FpC_sub(f(E,V), B, p);1376if (ZV_equal0(W)) return gerepilecopy(ltop, V);1377av = avma;1378for (i = 1; i <= n; ++i)1379{1380V = W;1381W = f(E, V);1382if (ZV_equal0(W))1383return gerepilecopy(ltop, shallowtrans(V));1384gerepileall(av, 2, &V, &W);1385}1386set_avma(btop);1387}1388return NULL;1389}13901391GEN1392zMs_ZC_mul(GEN M, GEN B)1393{1394long i, j;1395long n = lg(B)-1;1396GEN V = zerocol(n);1397for (i = 1; i <= n; ++i)1398if (signe(gel(B, i)))1399{1400GEN R = gel(M, i), C = gel(R, 1), E = gel(R, 2);1401long l = lg(C);1402for (j = 1; j < l; ++j)1403{1404long k = C[j];1405switch(E[j])1406{1407case 1:1408gel(V, k) = gel(V,k)==gen_0 ? gel(B,i) : addii(gel(V, k), gel(B,i));1409break;1410case -1:1411gel(V, k) = gel(V,k)==gen_0 ? negi(gel(B,i)) : subii(gel(V, k), gel(B,i));1412break;1413default:1414gel(V, k) = gel(V,k)==gen_0 ? mulis(gel(B, i), E[j]) : addii(gel(V, k), mulis(gel(B, i), E[j]));1415break;1416}1417}1418}1419return V;1420}14211422GEN1423FpMs_FpC_mul(GEN M, GEN B, GEN p) { return FpC_red(zMs_ZC_mul(M, B), p); }14241425GEN1426ZV_zMs_mul(GEN B, GEN M)1427{1428long i, j;1429long m = lg(M)-1;1430GEN V = cgetg(m+1,t_VEC);1431for (i = 1; i <= m; ++i)1432{1433GEN R = gel(M, i), C = gel(R, 1), E = gel(R, 2);1434long l = lg(C);1435GEN z;1436if (l == 1)1437{1438gel(V,i) = gen_0;1439continue;1440}1441z = mulis(gel(B, C[1]), E[1]);1442for (j = 2; j < l; ++j)1443{1444long k = C[j];1445switch(E[j])1446{1447case 1:1448z = addii(z, gel(B,k));1449break;1450case -1:1451z = subii(z, gel(B,k));1452break;1453default:1454z = addii(z, mulis(gel(B,k), E[j]));1455break;1456}1457}1458gel(V,i) = z;1459}1460return V;1461}14621463GEN1464FpV_FpMs_mul(GEN B, GEN M, GEN p) { return FpV_red(ZV_zMs_mul(B, M), p); }14651466GEN1467ZlM_gauss(GEN a, GEN b, ulong p, long e, GEN C)1468{1469pari_sp av = avma, av2;1470GEN xi, xb, pi = gen_1, P;1471long i;1472if (!C) {1473C = Flm_inv(ZM_to_Flm(a, p), p);1474if (!C) pari_err_INV("ZlM_gauss", a);1475}1476P = utoipos(p);1477av2 = avma;1478xi = Flm_mul(C, ZM_to_Flm(b, p), p);1479xb = Flm_to_ZM(xi);1480for (i = 2; i <= e; i++)1481{1482pi = muliu(pi, p); /* = p^(i-1) */1483b = ZM_Z_divexact(ZM_sub(b, ZM_nm_mul(a, xi)), P);1484if (gc_needed(av,2))1485{1486if(DEBUGMEM>1) pari_warn(warnmem,"ZlM_gauss. i=%ld/%ld",i,e);1487gerepileall(av2,3, &pi,&b,&xb);1488}1489xi = Flm_mul(C, ZM_to_Flm(b, p), p);1490xb = ZM_add(xb, nm_Z_mul(xi, pi));1491}1492return gerepileupto(av, xb);1493}14941495struct wrapper_modp_s {1496GEN (*f)(void*, GEN);1497void *E;1498GEN p;1499};15001501/* compute f(x) mod p */1502static GEN1503wrap_relcomb_modp(void *E, GEN x)1504{1505struct wrapper_modp_s *W = (struct wrapper_modp_s*)E;1506return FpC_red(W->f(W->E, x), W->p);1507}1508static GEN1509wrap_relcomb(void*E, GEN x) { return zMs_ZC_mul((GEN)E, x); }15101511static GEN1512wrap_relker(void*E, GEN x) { return ZV_zMs_mul(x, (GEN)E); }15131514/* Solve f(X) = B (mod p^e); blackbox version of ZlM_gauss */1515GEN1516gen_ZpM_Dixon_Wiedemann(void *E, GEN (*f)(void*, GEN), GEN B, GEN p, long e)1517{1518struct wrapper_modp_s W;1519pari_sp av = avma;1520GEN xb, xi, pi = gen_1;1521long i;1522W.E = E;1523W.f = f;1524W.p = p;1525xi = gen_FpM_Wiedemann((void*)&W, wrap_relcomb_modp, FpC_red(B, p), p); /* f^(-1) B */1526if (!xi || e == 1 || typ(xi) == t_VEC) return xi;1527xb = xi;1528for (i = 2; i <= e; i++)1529{1530pi = mulii(pi, p); /* = p^(i-1) */1531B = ZC_Z_divexact(ZC_sub(B, f(E, xi)), p);1532if (gc_needed(av,2))1533{1534if(DEBUGMEM>1) pari_warn(warnmem,"gen_ZpM_Dixon_Wiedemann. i=%ld",i);1535gerepileall(av,3, &pi,&B,&xb);1536}1537xi = gen_FpM_Wiedemann((void*)&W, wrap_relcomb_modp, FpC_red(B, p), p);1538if (!xi) return NULL;1539if (typ(xi) == t_VEC) return gerepileupto(av, xi);1540xb = ZC_add(xb, ZC_Z_mul(xi, pi));1541}1542return gerepileupto(av, xb);1543}15441545static GEN1546vecprow(GEN A, GEN prow)1547{1548return mkvec2(vecsmallpermute(prow,gel(A,1)), gel(A,2));1549}15501551/* Solve the equation MX = A. Return either a solution as a t_COL,1552* or the index of a column which is linearly dependent from the others as a1553* t_VECSMALL with a single component. */1554GEN1555ZpMs_ZpCs_solve(GEN M, GEN A, long nbrow, GEN p, long e)1556{1557pari_sp av = avma;1558GEN pcol, prow;1559long nbi=lg(M)-1, lR;1560long i, n;1561GEN Mp, Ap, Rp;1562pari_timer ti;1563if (DEBUGLEVEL) timer_start(&ti);1564RgMs_structelim(M, nbrow, gel(A, 1), &pcol, &prow);1565if (!pcol) return gc_NULL(av);1566if (DEBUGLEVEL)1567timer_printf(&ti,"structured elimination (%ld -> %ld)",nbi,lg(pcol)-1);1568n = lg(pcol)-1;1569Mp = cgetg(n+1, t_MAT);1570for(i=1; i<=n; i++)1571gel(Mp, i) = vecprow(gel(M,pcol[i]), prow);1572Ap = zCs_to_ZC(vecprow(A, prow), n);1573if (DEBUGLEVEL) timer_start(&ti);1574Rp = gen_ZpM_Dixon_Wiedemann((void*)Mp,wrap_relcomb, Ap, p, e);1575if (DEBUGLEVEL) timer_printf(&ti,"linear algebra");1576if (!Rp) return gc_NULL(av);1577lR = lg(Rp)-1;1578if (typ(Rp) == t_COL)1579{1580GEN R = zerocol(nbi+1);1581for (i=1; i<=lR; i++)1582gel(R,pcol[i]) = gel(Rp,i);1583return gerepilecopy(av, R);1584}1585for (i = 1; i <= lR; ++i)1586if (signe(gel(Rp, i)))1587return gerepileuptoleaf(av, mkvecsmall(pcol[i]));1588return NULL; /* LCOV_EXCL_LINE */1589}15901591GEN1592FpMs_FpCs_solve(GEN M, GEN A, long nbrow, GEN p)1593{1594return ZpMs_ZpCs_solve(M, A, nbrow, p, 1);1595}15961597GEN1598FpMs_FpCs_solve_safe(GEN M, GEN A, long nbrow, GEN p)1599{1600GEN res;1601pari_CATCH(e_INV)1602{1603GEN E = pari_err_last();1604GEN x = err_get_compo(E,2);1605if (typ(x) != t_INTMOD) pari_err(0,E);1606if (DEBUGLEVEL)1607pari_warn(warner,"FpMs_FpCs_solve_safe, impossible inverse %Ps", x);1608res = NULL;1609} pari_TRY {1610res = ZpMs_ZpCs_solve(M, A, nbrow, p, 1);1611} pari_ENDCATCH1612return res;1613}16141615static GEN1616FpMs_structelim_back(GEN M, GEN V, GEN prow, GEN p)1617{1618long i, j, oldf = 0, f = 0;1619long lrow = lg(prow), lM = lg(M);1620GEN W = const_vecsmall(lM-1,1);1621GEN R = cgetg(lrow, t_VEC);1622for (i=1; i<lrow; i++)1623gel(R,i) = prow[i]==0 ? NULL: gel(V,prow[i]);1624do1625{1626oldf = f;1627for(i=1; i<lM; i++)1628if (W[i])1629{1630GEN C = gel(M,i), F = gel(C,1), E = gel(C,2);1631long c=0, cj=0, lF = lg(F);1632for(j=1; j<lF; j++)1633if (!gel(R,F[j])) { c++; cj=j; }1634if (c>=2) continue;1635if (c==1)1636{1637pari_sp av = avma;1638GEN s = gen_0;1639for(j=1; j<lF; j++)1640if (j!=cj) s = Fp_add(s, mulis(gel(R,F[j]), E[j]), p);1641/* s /= -E[cj] mod p */1642s = E[cj] < 0? Fp_divu(s, -E[cj], p): Fp_divu(Fp_neg(s,p), E[cj], p);1643gel(R,F[cj]) = gerepileuptoint(av, s);1644f++;1645}1646W[i]=0;1647}1648} while(oldf!=f);1649for (i=1; i<lrow; i++)1650if (!gel(R,i)) gel(R,i) = gen_0;1651return R;1652}16531654/* Return a linear form Y such that YM is essentially 0 */1655GEN1656FpMs_leftkernel_elt_col(GEN M, long nbcol, long nbrow, GEN p)1657{1658pari_sp av = avma, av2;1659GEN pcol, prow;1660long i, n;1661GEN Mp, B, MB, R, Rp;1662pari_timer ti;1663struct wrapper_modp_s W;1664if (DEBUGLEVEL) timer_start(&ti);1665RgMs_structelim_col(M, nbcol, nbrow, cgetg(1,t_VECSMALL), &pcol, &prow);1666if (!pcol) return gc_NULL(av);1667if (DEBUGLEVEL)1668timer_printf(&ti,"structured elimination (%ld -> %ld)",nbcol,lg(pcol)-1);1669n = lg(pcol)-1;1670Mp = cgetg(n+1, t_MAT);1671for(i=1; i<=n; i++)1672gel(Mp, i) = vecprow(gel(M,pcol[i]), prow);1673W.E = (void*) Mp;1674W.f = wrap_relker;1675W.p = p;1676av2 = avma;1677for(;;)1678{1679set_avma(av2);1680B = random_FpV(n, p);1681MB = FpV_FpMs_mul(B, Mp, p);1682if (DEBUGLEVEL) timer_start(&ti);1683pari_CATCH(e_INV)1684{1685GEN E = pari_err_last();1686GEN x = err_get_compo(E,2);1687if (typ(x) != t_INTMOD) pari_err(0,E);1688if (DEBUGLEVEL)1689pari_warn(warner,"FpMs_leftkernel_elt, impossible inverse %Ps", x);1690Rp = NULL;1691} pari_TRY {1692Rp = gen_FpM_Wiedemann((void*)&W, wrap_relcomb_modp, MB, p);1693} pari_ENDCATCH1694if (!Rp) continue;1695if (typ(Rp)==t_COL)1696{1697Rp = FpC_sub(Rp,B,p);1698if (ZV_equal0(Rp)) continue;1699}1700R = FpMs_structelim_back(M, Rp, prow, p);1701if (DEBUGLEVEL) timer_printf(&ti,"Wiedemann left kernel");1702return gerepilecopy(av, R);1703}1704}17051706GEN1707FpMs_leftkernel_elt(GEN M, long nbrow, GEN p)1708{1709return FpMs_leftkernel_elt_col(M, lg(M)-1, nbrow, p);1710}171117121713