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"1617static int18check_ZV(GEN x, long l)19{20long i;21for (i=1; i<l; i++)22if (typ(gel(x,i)) != t_INT) return 0;23return 1;24}25void26RgV_check_ZV(GEN A, const char *s)27{28if (!RgV_is_ZV(A)) pari_err_TYPE(stack_strcat(s," [integer vector]"), A);29}30void31RgM_check_ZM(GEN A, const char *s)32{33long n = lg(A);34if (n != 1)35{36long j, m = lgcols(A);37for (j=1; j<n; j++)38if (!check_ZV(gel(A,j), m))39pari_err_TYPE(stack_strcat(s," [integer matrix]"), A);40}41}4243static long44ZV_max_lg_i(GEN x, long m)45{46long i, prec = 2;47for (i = 1; i < m; i++) prec = maxss(prec, lgefint(gel(x,i)));48return prec;49}50long51ZV_max_lg(GEN x) { return ZV_max_lg_i(x, lg(x)); }5253static long54ZM_max_lg_i(GEN x, long n, long m)55{56long j, prec = 2;57for (j = 1; j < n; j++) prec = maxss(prec, ZV_max_lg_i(gel(x,j), m));58return prec;59}60long61ZM_max_lg(GEN x)62{63long n = lg(x);64if (n == 1) return 2;65return ZM_max_lg_i(x, n, lgcols(x));66}6768GEN69ZM_supnorm(GEN x)70{71long i, j, h, lx = lg(x);72GEN s = gen_0;73if (lx == 1) return gen_1;74h = lgcols(x);75for (j=1; j<lx; j++)76{77GEN xj = gel(x,j);78for (i=1; i<h; i++)79{80GEN c = gel(xj,i);81if (abscmpii(c, s) > 0) s = c;82}83}84return absi(s);85}8687/********************************************************************/88/** **/89/** MULTIPLICATION **/90/** **/91/********************************************************************/92/* x nonempty ZM, y a compatible nc (dimension > 0). */93static GEN94ZM_nc_mul_i(GEN x, GEN y, long c, long l)95{96long i, j;97pari_sp av;98GEN z = cgetg(l,t_COL), s;99100for (i=1; i<l; i++)101{102av = avma; s = muliu(gcoeff(x,i,1),y[1]);103for (j=2; j<c; j++)104if (y[j]) s = addii(s, muliu(gcoeff(x,i,j),y[j]));105gel(z,i) = gerepileuptoint(av,s);106}107return z;108}109110/* x ZV, y a compatible zc. */111GEN112ZV_zc_mul(GEN x, GEN y)113{114long j, l = lg(x);115pari_sp av = avma;116GEN s = mulis(gel(x,1),y[1]);117for (j=2; j<l; j++)118if (y[j]) s = addii(s, mulis(gel(x,j),y[j]));119return gerepileuptoint(av,s);120}121122/* x nonempty ZM, y a compatible zc (dimension > 0). */123static GEN124ZM_zc_mul_i(GEN x, GEN y, long c, long l)125{126long i, j;127GEN z = cgetg(l,t_COL);128129for (i=1; i<l; i++)130{131pari_sp av = avma;132GEN s = mulis(gcoeff(x,i,1),y[1]);133for (j=2; j<c; j++)134if (y[j]) s = addii(s, mulis(gcoeff(x,i,j),y[j]));135gel(z,i) = gerepileuptoint(av,s);136}137return z;138}139GEN140ZM_zc_mul(GEN x, GEN y) {141long l = lg(x);142if (l == 1) return cgetg(1, t_COL);143return ZM_zc_mul_i(x,y, l, lgcols(x));144}145146/* y nonempty ZM, x a compatible zv (dimension > 0). */147GEN148zv_ZM_mul(GEN x, GEN y) {149long i,j, lx = lg(x), ly = lg(y);150GEN z;151if (lx == 1) return zerovec(ly-1);152z = cgetg(ly,t_VEC);153for (j=1; j<ly; j++)154{155pari_sp av = avma;156GEN s = mulsi(x[1], gcoeff(y,1,j));157for (i=2; i<lx; i++)158if (x[i]) s = addii(s, mulsi(x[i], gcoeff(y,i,j)));159gel(z,j) = gerepileuptoint(av,s);160}161return z;162}163164/* x ZM, y a compatible zm (dimension > 0). */165GEN166ZM_zm_mul(GEN x, GEN y)167{168long j, c, l = lg(x), ly = lg(y);169GEN z = cgetg(ly, t_MAT);170if (l == 1) return z;171c = lgcols(x);172for (j = 1; j < ly; j++) gel(z,j) = ZM_zc_mul_i(x, gel(y,j), l,c);173return z;174}175/* x ZM, y a compatible zn (dimension > 0). */176GEN177ZM_nm_mul(GEN x, GEN y)178{179long j, c, l = lg(x), ly = lg(y);180GEN z = cgetg(ly, t_MAT);181if (l == 1) return z;182c = lgcols(x);183for (j = 1; j < ly; j++) gel(z,j) = ZM_nc_mul_i(x, gel(y,j), l,c);184return z;185}186187/* Strassen-Winograd algorithm */188189/* Return A[ma+1..ma+da, na+1..na+ea] - B[mb+1..mb+db, nb+1..nb+eb]190* as an (m x n)-matrix, padding the input with zeroes as necessary. */191static GEN192add_slices(long m, long n,193GEN A, long ma, long da, long na, long ea,194GEN B, long mb, long db, long nb, long eb)195{196long min_d = minss(da, db), min_e = minss(ea, eb), i, j;197GEN M = cgetg(n + 1, t_MAT), C;198199for (j = 1; j <= min_e; j++) {200gel(M, j) = C = cgetg(m + 1, t_COL);201for (i = 1; i <= min_d; i++)202gel(C, i) = addii(gcoeff(A, ma + i, na + j),203gcoeff(B, mb + i, nb + j));204for (; i <= da; i++) gel(C, i) = gcoeff(A, ma + i, na + j);205for (; i <= db; i++) gel(C, i) = gcoeff(B, mb + i, nb + j);206for (; i <= m; i++) gel(C, i) = gen_0;207}208for (; j <= ea; j++) {209gel(M, j) = C = cgetg(m + 1, t_COL);210for (i = 1; i <= da; i++) gel(C, i) = gcoeff(A, ma + i, na + j);211for (; i <= m; i++) gel(C, i) = gen_0;212}213for (; j <= eb; j++) {214gel(M, j) = C = cgetg(m + 1, t_COL);215for (i = 1; i <= db; i++) gel(C, i) = gcoeff(B, mb + i, nb + j);216for (; i <= m; i++) gel(C, i) = gen_0;217}218for (; j <= n; j++) gel(M, j) = zerocol(m);219return M;220}221222/* Return A[ma+1..ma+da, na+1..na+ea] - B[mb+1..mb+db, nb+1..nb+eb]223* as an (m x n)-matrix, padding the input with zeroes as necessary. */224static GEN225subtract_slices(long m, long n,226GEN A, long ma, long da, long na, long ea,227GEN B, long mb, long db, long nb, long eb)228{229long min_d = minss(da, db), min_e = minss(ea, eb), i, j;230GEN M = cgetg(n + 1, t_MAT), C;231232for (j = 1; j <= min_e; j++) {233gel(M, j) = C = cgetg(m + 1, t_COL);234for (i = 1; i <= min_d; i++)235gel(C, i) = subii(gcoeff(A, ma + i, na + j),236gcoeff(B, mb + i, nb + j));237for (; i <= da; i++) gel(C, i) = gcoeff(A, ma + i, na + j);238for (; i <= db; i++) gel(C, i) = negi(gcoeff(B, mb + i, nb + j));239for (; i <= m; i++) gel(C, i) = gen_0;240}241for (; j <= ea; j++) {242gel(M, j) = C = cgetg(m + 1, t_COL);243for (i = 1; i <= da; i++) gel(C, i) = gcoeff(A, ma + i, na + j);244for (; i <= m; i++) gel(C, i) = gen_0;245}246for (; j <= eb; j++) {247gel(M, j) = C = cgetg(m + 1, t_COL);248for (i = 1; i <= db; i++) gel(C, i) = negi(gcoeff(B, mb + i, nb + j));249for (; i <= m; i++) gel(C, i) = gen_0;250}251for (; j <= n; j++) gel(M, j) = zerocol(m);252return M;253}254255static GEN ZM_mul_i(GEN x, GEN y, long l, long lx, long ly);256257/* Strassen-Winograd matrix product A (m x n) * B (n x p) */258static GEN259ZM_mul_sw(GEN A, GEN B, long m, long n, long p)260{261pari_sp av = avma;262long m1 = (m + 1)/2, m2 = m/2,263n1 = (n + 1)/2, n2 = n/2,264p1 = (p + 1)/2, p2 = p/2;265GEN A11, A12, A22, B11, B21, B22,266S1, S2, S3, S4, T1, T2, T3, T4,267M1, M2, M3, M4, M5, M6, M7,268V1, V2, V3, C11, C12, C21, C22, C;269270T2 = subtract_slices(n1, p2, B, 0, n1, p1, p2, B, n1, n2, p1, p2);271S1 = subtract_slices(m2, n1, A, m1, m2, 0, n1, A, 0, m2, 0, n1);272M2 = ZM_mul_i(S1, T2, m2 + 1, n1 + 1, p2 + 1);273if (gc_needed(av, 1))274gerepileall(av, 2, &T2, &M2); /* destroy S1 */275T3 = subtract_slices(n1, p1, T2, 0, n1, 0, p2, B, 0, n1, 0, p1);276if (gc_needed(av, 1))277gerepileall(av, 2, &M2, &T3); /* destroy T2 */278S2 = add_slices(m2, n1, A, m1, m2, 0, n1, A, m1, m2, n1, n2);279T1 = subtract_slices(n1, p1, B, 0, n1, p1, p2, B, 0, n1, 0, p2);280M3 = ZM_mul_i(S2, T1, m2 + 1, n1 + 1, p2 + 1);281if (gc_needed(av, 1))282gerepileall(av, 4, &M2, &T3, &S2, &M3); /* destroy T1 */283S3 = subtract_slices(m1, n1, S2, 0, m2, 0, n1, A, 0, m1, 0, n1);284if (gc_needed(av, 1))285gerepileall(av, 4, &M2, &T3, &M3, &S3); /* destroy S2 */286A11 = matslice(A, 1, m1, 1, n1);287B11 = matslice(B, 1, n1, 1, p1);288M1 = ZM_mul_i(A11, B11, m1 + 1, n1 + 1, p1 + 1);289if (gc_needed(av, 1))290gerepileall(av, 5, &M2, &T3, &M3, &S3, &M1); /* destroy A11, B11 */291A12 = matslice(A, 1, m1, n1 + 1, n);292B21 = matslice(B, n1 + 1, n, 1, p1);293M4 = ZM_mul_i(A12, B21, m1 + 1, n2 + 1, p1 + 1);294if (gc_needed(av, 1))295gerepileall(av, 6, &M2, &T3, &M3, &S3, &M1, &M4); /* destroy A12, B21 */296C11 = add_slices(m1, p1, M1, 0, m1, 0, p1, M4, 0, m1, 0, p1);297if (gc_needed(av, 1))298gerepileall(av, 6, &M2, &T3, &M3, &S3, &M1, &C11); /* destroy M4 */299M5 = ZM_mul_i(S3, T3, m1 + 1, n1 + 1, p1 + 1);300S4 = subtract_slices(m1, n2, A, 0, m1, n1, n2, S3, 0, m1, 0, n2);301if (gc_needed(av, 1))302gerepileall(av, 7, &M2, &T3, &M3, &M1, &C11, &M5, &S4); /* destroy S3 */303T4 = add_slices(n2, p1, B, n1, n2, 0, p1, T3, 0, n2, 0, p1);304if (gc_needed(av, 1))305gerepileall(av, 7, &M2, &M3, &M1, &C11, &M5, &S4, &T4); /* destroy T3 */306V1 = subtract_slices(m1, p1, M1, 0, m1, 0, p1, M5, 0, m1, 0, p1);307if (gc_needed(av, 1))308gerepileall(av, 6, &M2, &M3, &S4, &T4, &C11, &V1); /* destroy M1, M5 */309B22 = matslice(B, n1 + 1, n, p1 + 1, p);310M6 = ZM_mul_i(S4, B22, m1 + 1, n2 + 1, p2 + 1);311if (gc_needed(av, 1))312gerepileall(av, 6, &M2, &M3, &T4, &C11, &V1, &M6); /* destroy S4, B22 */313A22 = matslice(A, m1 + 1, m, n1 + 1, n);314M7 = ZM_mul_i(A22, T4, m2 + 1, n2 + 1, p1 + 1);315if (gc_needed(av, 1))316gerepileall(av, 6, &M2, &M3, &C11, &V1, &M6, &M7); /* destroy A22, T4 */317V3 = add_slices(m1, p2, V1, 0, m1, 0, p2, M3, 0, m2, 0, p2);318C12 = add_slices(m1, p2, V3, 0, m1, 0, p2, M6, 0, m1, 0, p2);319if (gc_needed(av, 1))320gerepileall(av, 6, &M2, &M3, &C11, &V1, &M7, &C12); /* destroy V3, M6 */321V2 = add_slices(m2, p1, V1, 0, m2, 0, p1, M2, 0, m2, 0, p2);322if (gc_needed(av, 1))323gerepileall(av, 5, &M3, &C11, &M7, &C12, &V2); /* destroy V1, M2 */324C21 = add_slices(m2, p1, V2, 0, m2, 0, p1, M7, 0, m2, 0, p1);325if (gc_needed(av, 1))326gerepileall(av, 5, &M3, &C11, &C12, &V2, &C21); /* destroy M7 */327C22 = add_slices(m2, p2, V2, 0, m2, 0, p2, M3, 0, m2, 0, p2);328if (gc_needed(av, 1))329gerepileall(av, 4, &C11, &C12, &C21, &C22); /* destroy V2, M3 */330C = shallowconcat(vconcat(C11, C21), vconcat(C12, C22));331return gerepilecopy(av, C);332}333334/* x[i,]*y. Assume lg(x) > 1 and 0 < i < lgcols(x) */335static GEN336ZMrow_ZC_mul_i(GEN x, GEN y, long i, long lx)337{338pari_sp av = avma;339GEN c = mulii(gcoeff(x,i,1), gel(y,1)), ZERO = gen_0;340long k;341for (k = 2; k < lx; k++)342{343GEN t = mulii(gcoeff(x,i,k), gel(y,k));344if (t != ZERO) c = addii(c, t);345}346return gerepileuptoint(av, c);347}348GEN349ZMrow_ZC_mul(GEN x, GEN y, long i)350{ return ZMrow_ZC_mul_i(x, y, i, lg(x)); }351352/* return x * y, 1 < lx = lg(x), l = lgcols(x) */353static GEN354ZM_ZC_mul_i(GEN x, GEN y, long lx, long l)355{356GEN z = cgetg(l,t_COL);357long i;358for (i=1; i<l; i++) gel(z,i) = ZMrow_ZC_mul_i(x,y,i,lx);359return z;360}361362static GEN363ZM_mul_classical(GEN x, GEN y, long l, long lx, long ly)364{365long j;366GEN z = cgetg(ly, t_MAT);367for (j = 1; j < ly; j++)368gel(z, j) = ZM_ZC_mul_i(x, gel(y, j), lx, l);369return z;370}371372static GEN373ZM_mul_slice(GEN A, GEN B, GEN P, GEN *mod)374{375pari_sp av = avma;376long i, n = lg(P)-1;377GEN H, T;378if (n == 1)379{380ulong p = uel(P,1);381GEN a = ZM_to_Flm(A, p);382GEN b = ZM_to_Flm(B, p);383GEN Hp = gerepileupto(av, Flm_to_ZM(Flm_mul(a, b, p)));384*mod = utoi(p); return Hp;385}386T = ZV_producttree(P);387A = ZM_nv_mod_tree(A, P, T);388B = ZM_nv_mod_tree(B, P, T);389H = cgetg(n+1, t_VEC);390for(i=1; i <= n; i++)391gel(H,i) = Flm_mul(gel(A,i),gel(B,i),P[i]);392H = nmV_chinese_center_tree_seq(H, P, T, ZV_chinesetree(P, T));393*mod = gmael(T, lg(T)-1, 1);394gerepileall(av, 2, &H, mod);395return H;396}397398GEN399ZM_mul_worker(GEN P, GEN A, GEN B)400{401GEN V = cgetg(3, t_VEC);402gel(V,1) = ZM_mul_slice(A, B, P, &gel(V,2));403return V;404}405406static GEN407ZM_mul_fast(GEN A, GEN B, long lA, long lB, long sA, long sB)408{409pari_sp av = avma;410forprime_t S;411GEN H, worker;412long h;413if (sA == 2 || sB == 2) return zeromat(nbrows(A),lB-1);414h = 1 + (sA + sB - 4) * BITS_IN_LONG + expu(lA-1);415init_modular_big(&S);416worker = snm_closure(is_entry("_ZM_mul_worker"), mkvec2(A,B));417H = gen_crt("ZM_mul", worker, &S, NULL, h, 0, NULL,418nmV_chinese_center, FpM_center);419return gerepileupto(av, H);420}421422/* s = min(log_BIL |x|, log_BIL |y|), use Strassen-Winograd when423* min(dims) > B */424static long425sw_bound(long s)426{ return s > 60 ? 2: s > 25 ? 4: s > 15 ? 8 : s > 8 ? 16 : 32; }427428static GEN429ZM_mul_i(GEN x, GEN y, long l, long lx, long ly)430{431long sx, sy, B;432if (lx==3 && l==3 && ly==3) return ZM2_mul(x,y);433sx = ZM_max_lg_i(x, lx, l);434sy = ZM_max_lg_i(y, ly, lx);435if ((lx > 70 && ly > 70 && l > 70) && (sx <= 10 * sy && sy <= 10 * sx))436return ZM_mul_fast(x, y, lx, ly, sx, sy);437438B = sw_bound(minss(sx, sy));439if (l <= B || lx <= B || ly <= B)440return ZM_mul_classical(x, y, l, lx, ly);441else442return ZM_mul_sw(x, y, l - 1, lx - 1, ly - 1);443}444445GEN446ZM_mul(GEN x, GEN y)447{448long lx = lg(x), ly = lg(y);449if (ly == 1) return cgetg(1,t_MAT);450if (lx == 1) return zeromat(0, ly-1);451return ZM_mul_i(x, y, lgcols(x), lx, ly);452}453454GEN455QM_mul(GEN x, GEN y)456{457GEN dx, nx = Q_primitive_part(x, &dx);458GEN dy, ny = Q_primitive_part(y, &dy);459GEN z = ZM_mul(nx, ny);460if (dx || dy)461{462GEN d = dx ? dy ? gmul(dx, dy): dx : dy;463if (!gequal1(d)) z = ZM_Q_mul(z, d);464}465return z;466}467468GEN469QM_sqr(GEN x)470{471GEN dx, nx = Q_primitive_part(x, &dx);472GEN z = ZM_sqr(nx);473if (dx)474z = ZM_Q_mul(z, gsqr(dx));475return z;476}477478GEN479QM_QC_mul(GEN x, GEN y)480{481GEN dx, nx = Q_primitive_part(x, &dx);482GEN dy, ny = Q_primitive_part(y, &dy);483GEN z = ZM_ZC_mul(nx, ny);484if (dx || dy)485{486GEN d = dx ? dy ? gmul(dx, dy): dx : dy;487if (!gequal1(d)) z = ZC_Q_mul(z, d);488}489return z;490}491492/* assume result is symmetric */493GEN494ZM_multosym(GEN x, GEN y)495{496long j, lx, ly = lg(y);497GEN M;498if (ly == 1) return cgetg(1,t_MAT);499lx = lg(x); /* = lgcols(y) */500if (lx == 1) return cgetg(1,t_MAT);501/* ly = lgcols(x) */502M = cgetg(ly, t_MAT);503for (j=1; j<ly; j++)504{505GEN z = cgetg(ly,t_COL), yj = gel(y,j);506long i;507for (i=1; i<j; i++) gel(z,i) = gcoeff(M,j,i);508for (i=j; i<ly; i++)gel(z,i) = ZMrow_ZC_mul_i(x,yj,i,lx);509gel(M,j) = z;510}511return M;512}513514/* compute m*diagonal(d), assume lg(d) = lg(m). Shallow */515GEN516ZM_mul_diag(GEN m, GEN d)517{518long j, l;519GEN y = cgetg_copy(m, &l);520for (j=1; j<l; j++)521{522GEN c = gel(d,j);523gel(y,j) = equali1(c)? gel(m,j): ZC_Z_mul(gel(m,j), c);524}525return y;526}527/* compute diagonal(d)*m, assume lg(d) = lg(m~). Shallow */528GEN529ZM_diag_mul(GEN d, GEN m)530{531long i, j, l = lg(d), lm = lg(m);532GEN y = cgetg(lm, t_MAT);533for (j=1; j<lm; j++) gel(y,j) = cgetg(l, t_COL);534for (i=1; i<l; i++)535{536GEN c = gel(d,i);537if (equali1(c))538for (j=1; j<lm; j++) gcoeff(y,i,j) = gcoeff(m,i,j);539else540for (j=1; j<lm; j++) gcoeff(y,i,j) = mulii(gcoeff(m,i,j), c);541}542return y;543}544545/* assume lx > 1 is lg(x) = lg(y) */546static GEN547ZV_dotproduct_i(GEN x, GEN y, long lx)548{549pari_sp av = avma;550GEN c = mulii(gel(x,1), gel(y,1));551long i;552for (i = 2; i < lx; i++)553{554GEN t = mulii(gel(x,i), gel(y,i));555if (t != gen_0) c = addii(c, t);556}557return gerepileuptoint(av, c);558}559560/* x~ * y, assuming result is symmetric */561GEN562ZM_transmultosym(GEN x, GEN y)563{564long i, j, l, ly = lg(y);565GEN M;566if (ly == 1) return cgetg(1,t_MAT);567/* lg(x) = ly */568l = lgcols(y); /* = lgcols(x) */569M = cgetg(ly, t_MAT);570for (i=1; i<ly; i++)571{572GEN xi = gel(x,i), c = cgetg(ly,t_COL);573gel(M,i) = c;574for (j=1; j<i; j++)575gcoeff(M,i,j) = gel(c,j) = ZV_dotproduct_i(xi,gel(y,j),l);576gel(c,i) = ZV_dotproduct_i(xi,gel(y,i),l);577}578return M;579}580/* x~ * y */581GEN582ZM_transmul(GEN x, GEN y)583{584long i, j, l, lx, ly = lg(y);585GEN M;586if (ly == 1) return cgetg(1,t_MAT);587lx = lg(x);588l = lgcols(y);589if (lgcols(x) != l) pari_err_OP("operation 'ZM_transmul'", x,y);590M = cgetg(ly, t_MAT);591for (i=1; i<ly; i++)592{593GEN yi = gel(y,i), c = cgetg(lx,t_COL);594gel(M,i) = c;595for (j=1; j<lx; j++) gel(c,j) = ZV_dotproduct_i(yi,gel(x,j),l);596}597return M;598}599600static GEN601ZM_sqr_i(GEN x, long l)602{603long s = ZM_max_lg_i(x, l, l);604if (l > 70) return ZM_mul_fast(x, x, l, l, s, s);605if (l <= sw_bound(s))606return ZM_mul_classical(x, x, l, l, l);607else608return ZM_mul_sw(x, x, l - 1, l - 1, l - 1);609}610611GEN612ZM_sqr(GEN x)613{614long lx=lg(x);615if (lx==1) return cgetg(1,t_MAT);616return ZM_sqr_i(x, lx);617}618GEN619ZM_ZC_mul(GEN x, GEN y)620{621long lx = lg(x);622return lx==1? cgetg(1,t_COL): ZM_ZC_mul_i(x, y, lx, lgcols(x));623}624625GEN626ZC_Z_div(GEN x, GEN c)627{ pari_APPLY_type(t_COL, Qdivii(gel(x,i), c)) }628629GEN630ZM_Z_div(GEN x, GEN c)631{ pari_APPLY_same(ZC_Z_div(gel(x, i), c)) }632633GEN634ZC_Q_mul(GEN A, GEN z)635{636pari_sp av = avma;637long i, l = lg(A);638GEN d, n, Ad, B, u;639if (typ(z)==t_INT) return ZC_Z_mul(A,z);640n = gel(z, 1); d = gel(z, 2);641Ad = FpC_red(A, d);642u = gcdii(d, FpV_factorback(Ad, NULL, d));643B = cgetg(l, t_COL);644if (equali1(u))645{646for(i=1; i<l; i++)647gel(B, i) = mkfrac(mulii(n, gel(A,i)), d);648} else649{650for(i=1; i<l; i++)651{652GEN di = gcdii(gel(Ad, i), u), ni = mulii(n, diviiexact(gel(A,i), di));653if (equalii(d, di))654gel(B, i) = ni;655else656gel(B, i) = mkfrac(ni, diviiexact(d, di));657}658}659return gerepilecopy(av, B);660}661662GEN663ZM_Q_mul(GEN x, GEN z)664{665if (typ(z)==t_INT) return ZM_Z_mul(x,z);666pari_APPLY_same(ZC_Q_mul(gel(x, i), z));667}668669long670zv_dotproduct(GEN x, GEN y)671{672long i, lx = lg(x);673ulong c;674if (lx == 1) return 0;675c = uel(x,1)*uel(y,1);676for (i = 2; i < lx; i++)677c += uel(x,i)*uel(y,i);678return c;679}680681GEN682ZV_ZM_mul(GEN x, GEN y)683{684long i, lx = lg(x), ly = lg(y);685GEN z;686if (lx == 1) return zerovec(ly-1);687z = cgetg(ly, t_VEC);688for (i = 1; i < ly; i++) gel(z,i) = ZV_dotproduct_i(x, gel(y,i), lx);689return z;690}691692GEN693ZC_ZV_mul(GEN x, GEN y)694{695long i,j, lx=lg(x), ly=lg(y);696GEN z;697if (ly==1) return cgetg(1,t_MAT);698z = cgetg(ly,t_MAT);699for (j=1; j < ly; j++)700{701gel(z,j) = cgetg(lx,t_COL);702for (i=1; i<lx; i++) gcoeff(z,i,j) = mulii(gel(x,i),gel(y,j));703}704return z;705}706707GEN708ZV_dotsquare(GEN x)709{710long i, lx;711pari_sp av;712GEN z;713lx = lg(x);714if (lx == 1) return gen_0;715av = avma; z = sqri(gel(x,1));716for (i=2; i<lx; i++) z = addii(z, sqri(gel(x,i)));717return gerepileuptoint(av,z);718}719720GEN721ZV_dotproduct(GEN x,GEN y)722{723long lx;724if (x == y) return ZV_dotsquare(x);725lx = lg(x);726if (lx == 1) return gen_0;727return ZV_dotproduct_i(x, y, lx);728}729730static GEN731_ZM_mul(void *data /*ignored*/, GEN x, GEN y)732{ (void)data; return ZM_mul(x,y); }733static GEN734_ZM_sqr(void *data /*ignored*/, GEN x)735{ (void)data; return ZM_sqr(x); }736GEN737ZM_pow(GEN x, GEN n)738{739pari_sp av = avma;740if (!signe(n)) return matid(lg(x)-1);741return gerepilecopy(av, gen_pow_i(x, n, NULL, &_ZM_sqr, &_ZM_mul));742}743GEN744ZM_powu(GEN x, ulong n)745{746pari_sp av = avma;747if (!n) return matid(lg(x)-1);748return gerepilecopy(av, gen_powu_i(x, n, NULL, &_ZM_sqr, &_ZM_mul));749}750/********************************************************************/751/** **/752/** ADD, SUB **/753/** **/754/********************************************************************/755static GEN756ZC_add_i(GEN x, GEN y, long lx)757{758GEN A = cgetg(lx, t_COL);759long i;760for (i=1; i<lx; i++) gel(A,i) = addii(gel(x,i), gel(y,i));761return A;762}763GEN764ZC_add(GEN x, GEN y) { return ZC_add_i(x, y, lg(x)); }765GEN766ZC_Z_add(GEN x, GEN y)767{768long k, lx = lg(x);769GEN z = cgetg(lx, t_COL);770if (lx == 1) pari_err_TYPE2("+",x,y);771gel(z,1) = addii(y,gel(x,1));772for (k = 2; k < lx; k++) gel(z,k) = icopy(gel(x,k));773return z;774}775776static GEN777ZC_sub_i(GEN x, GEN y, long lx)778{779long i;780GEN A = cgetg(lx, t_COL);781for (i=1; i<lx; i++) gel(A,i) = subii(gel(x,i), gel(y,i));782return A;783}784GEN785ZC_sub(GEN x, GEN y) { return ZC_sub_i(x, y, lg(x)); }786GEN787ZC_Z_sub(GEN x, GEN y)788{789long k, lx = lg(x);790GEN z = cgetg(lx, t_COL);791if (lx == 1) pari_err_TYPE2("+",x,y);792gel(z,1) = subii(gel(x,1), y);793for (k = 2; k < lx; k++) gel(z,k) = icopy(gel(x,k));794return z;795}796GEN797Z_ZC_sub(GEN a, GEN x)798{799long k, lx = lg(x);800GEN z = cgetg(lx, t_COL);801if (lx == 1) pari_err_TYPE2("-",a,x);802gel(z,1) = subii(a, gel(x,1));803for (k = 2; k < lx; k++) gel(z,k) = negi(gel(x,k));804return z;805}806807GEN808ZM_add(GEN x, GEN y)809{810long lx = lg(x), l, j;811GEN z;812if (lx == 1) return cgetg(1, t_MAT);813z = cgetg(lx, t_MAT); l = lgcols(x);814for (j = 1; j < lx; j++) gel(z,j) = ZC_add_i(gel(x,j), gel(y,j), l);815return z;816}817GEN818ZM_sub(GEN x, GEN y)819{820long lx = lg(x), l, j;821GEN z;822if (lx == 1) return cgetg(1, t_MAT);823z = cgetg(lx, t_MAT); l = lgcols(x);824for (j = 1; j < lx; j++) gel(z,j) = ZC_sub_i(gel(x,j), gel(y,j), l);825return z;826}827/********************************************************************/828/** **/829/** LINEAR COMBINATION **/830/** **/831/********************************************************************/832/* return X/c assuming division is exact */833GEN834ZC_Z_divexact(GEN x, GEN c)835{ pari_APPLY_type(t_COL, diviiexact(gel(x,i), c)) }836GEN837ZC_divexactu(GEN x, ulong c)838{ pari_APPLY_type(t_COL, diviuexact(gel(x,i), c)) }839840GEN841ZM_Z_divexact(GEN x, GEN c)842{ pari_APPLY_same(ZC_Z_divexact(gel(x,i), c)) }843844GEN845ZM_divexactu(GEN x, ulong c)846{ pari_APPLY_same(ZC_divexactu(gel(x,i), c)) }847848GEN849ZC_Z_mul(GEN x, GEN c)850{851if (!signe(c)) return zerocol(lg(x)-1);852if (is_pm1(c)) return (signe(c) > 0)? ZC_copy(x): ZC_neg(x);853pari_APPLY_type(t_COL, mulii(gel(x,i), c))854}855856GEN857ZC_z_mul(GEN x, long c)858{859if (!c) return zerocol(lg(x)-1);860if (c == 1) return ZC_copy(x);861if (c ==-1) return ZC_neg(x);862pari_APPLY_type(t_COL, mulsi(c, gel(x,i)))863}864865GEN866zv_z_mul(GEN x, long n)867{ pari_APPLY_long(x[i]*n) }868869/* return a ZM */870GEN871nm_Z_mul(GEN X, GEN c)872{873long i, j, h, l = lg(X), s = signe(c);874GEN A;875if (l == 1) return cgetg(1, t_MAT);876h = lgcols(X);877if (!s) return zeromat(h-1, l-1);878if (is_pm1(c)) {879if (s > 0) return Flm_to_ZM(X);880X = Flm_to_ZM(X); ZM_togglesign(X); return X;881}882A = cgetg(l, t_MAT);883for (j = 1; j < l; j++)884{885GEN a = cgetg(h, t_COL), x = gel(X, j);886for (i = 1; i < h; i++) gel(a,i) = muliu(c, x[i]);887gel(A,j) = a;888}889return A;890}891GEN892ZM_Z_mul(GEN X, GEN c)893{894long i, j, h, l = lg(X);895GEN A;896if (l == 1) return cgetg(1, t_MAT);897h = lgcols(X);898if (!signe(c)) return zeromat(h-1, l-1);899if (is_pm1(c)) return (signe(c) > 0)? ZM_copy(X): ZM_neg(X);900A = cgetg(l, t_MAT);901for (j = 1; j < l; j++)902{903GEN a = cgetg(h, t_COL), x = gel(X, j);904for (i = 1; i < h; i++) gel(a,i) = mulii(c, gel(x,i));905gel(A,j) = a;906}907return A;908}909void910ZC_lincomb1_inplace_i(GEN X, GEN Y, GEN v, long n)911{912long i;913for (i = n; i; i--) gel(X,i) = addmulii_inplace(gel(X,i), gel(Y,i), v);914}915/* X <- X + v Y (elementary col operation) */916void917ZC_lincomb1_inplace(GEN X, GEN Y, GEN v)918{919if (lgefint(v) != 2) return ZC_lincomb1_inplace_i(X, Y, v, lg(X)-1);920}921void922Flc_lincomb1_inplace(GEN X, GEN Y, ulong v, ulong q)923{924long i;925if (!v) return; /* v = 0 */926for (i = lg(X)-1; i; i--) X[i] = Fl_add(X[i], Fl_mul(Y[i], v, q), q);927}928929/* X + v Y, wasteful if (v = 0) */930static GEN931ZC_lincomb1(GEN v, GEN x, GEN y)932{ pari_APPLY_type(t_COL, addmulii(gel(x,i), gel(y,i), v)) }933934/* -X + vY */935static GEN936ZC_lincomb_1(GEN v, GEN x, GEN y)937{ pari_APPLY_type(t_COL, mulsubii(gel(y,i), v, gel(x,i))) }938939/* X,Y compatible ZV; u,v in Z. Returns A = u*X + v*Y */940GEN941ZC_lincomb(GEN u, GEN v, GEN X, GEN Y)942{943long su, sv;944GEN A;945946su = signe(u); if (!su) return ZC_Z_mul(Y, v);947sv = signe(v); if (!sv) return ZC_Z_mul(X, u);948if (is_pm1(v))949{950if (is_pm1(u))951{952if (su != sv) A = ZC_sub(X, Y);953else A = ZC_add(X, Y);954if (su < 0) ZV_togglesign(A); /* in place but was created above */955}956else957{958if (sv > 0) A = ZC_lincomb1 (u, Y, X);959else A = ZC_lincomb_1(u, Y, X);960}961}962else if (is_pm1(u))963{964if (su > 0) A = ZC_lincomb1 (v, X, Y);965else A = ZC_lincomb_1(v, X, Y);966}967else968{ /* not cgetg_copy: x may be a t_VEC */969long i, lx = lg(X);970A = cgetg(lx,t_COL);971for (i=1; i<lx; i++) gel(A,i) = lincombii(u,v,gel(X,i),gel(Y,i));972}973return A;974}975976/********************************************************************/977/** **/978/** CONVERSIONS **/979/** **/980/********************************************************************/981GEN982ZV_to_nv(GEN x)983{ pari_APPLY_ulong(itou(gel(x,i))) }984985GEN986zm_to_ZM(GEN x)987{ pari_APPLY_type(t_MAT, zc_to_ZC(gel(x,i))) }988989GEN990zmV_to_ZMV(GEN x)991{ pari_APPLY_type(t_VEC, zm_to_ZM(gel(x,i))) }992993/* same as Flm_to_ZM but do not assume positivity */994GEN995ZM_to_zm(GEN x)996{ pari_APPLY_same(ZV_to_zv(gel(x,i))) }997998GEN999zv_to_Flv(GEN x, ulong p)1000{ pari_APPLY_ulong(umodsu(x[i], p)) }10011002GEN1003zm_to_Flm(GEN x, ulong p)1004{ pari_APPLY_same(zv_to_Flv(gel(x,i),p)) }10051006GEN1007ZMV_to_zmV(GEN x)1008{ pari_APPLY_type(t_VEC, ZM_to_zm(gel(x,i))) }10091010/********************************************************************/1011/** **/1012/** COPY, NEGATION **/1013/** **/1014/********************************************************************/1015GEN1016ZC_copy(GEN x)1017{1018long i, lx = lg(x);1019GEN y = cgetg(lx, t_COL);1020for (i=1; i<lx; i++)1021{1022GEN c = gel(x,i);1023gel(y,i) = lgefint(c) == 2? gen_0: icopy(c);1024}1025return y;1026}10271028GEN1029ZM_copy(GEN x)1030{ pari_APPLY_same(ZC_copy(gel(x,i))) }10311032void1033ZV_neg_inplace(GEN M)1034{1035long l = lg(M);1036while (--l > 0) gel(M,l) = negi(gel(M,l));1037}1038GEN1039ZC_neg(GEN x)1040{ pari_APPLY_type(t_COL, negi(gel(x,i))) }10411042GEN1043zv_neg(GEN x)1044{ pari_APPLY_long(-x[i]) }1045GEN1046zv_neg_inplace(GEN M)1047{1048long l = lg(M);1049while (--l > 0) M[l] = -M[l];1050return M;1051}1052GEN1053zv_abs(GEN x)1054{ pari_APPLY_ulong(labs(x[i])) }1055GEN1056ZM_neg(GEN x)1057{ pari_APPLY_same(ZC_neg(gel(x,i))) }10581059void1060ZV_togglesign(GEN M)1061{1062long l = lg(M);1063while (--l > 0) togglesign_safe(&gel(M,l));1064}1065void1066ZM_togglesign(GEN M)1067{1068long l = lg(M);1069while (--l > 0) ZV_togglesign(gel(M,l));1070}10711072/********************************************************************/1073/** **/1074/** "DIVISION" mod HNF **/1075/** **/1076/********************************************************************/1077/* Reduce ZC x modulo ZM y in HNF, may return x itself (not a copy) */1078GEN1079ZC_hnfremdiv(GEN x, GEN y, GEN *Q)1080{1081long i, l = lg(x);1082GEN q;10831084if (Q) *Q = cgetg(l,t_COL);1085if (l == 1) return cgetg(1,t_COL);1086for (i = l-1; i>0; i--)1087{1088q = diviiround(gel(x,i), gcoeff(y,i,i));1089if (signe(q)) {1090togglesign(q);1091x = ZC_lincomb(gen_1, q, x, gel(y,i));1092}1093if (Q) gel(*Q, i) = q;1094}1095return x;1096}10971098/* x = y Q + R, may return some columns of x (not copies) */1099GEN1100ZM_hnfdivrem(GEN x, GEN y, GEN *Q)1101{1102long lx = lg(x), i;1103GEN R = cgetg(lx, t_MAT);1104if (Q)1105{1106GEN q = cgetg(lx, t_MAT); *Q = q;1107for (i=1; i<lx; i++) gel(R,i) = ZC_hnfremdiv(gel(x,i),y,(GEN*)(q+i));1108}1109else1110for (i=1; i<lx; i++)1111{1112pari_sp av = avma;1113GEN z = ZC_hnfrem(gel(x,i),y);1114gel(R,i) = (avma == av)? ZC_copy(z): gerepileupto(av, z);1115}1116return R;1117}11181119/********************************************************************/1120/** **/1121/** TESTS **/1122/** **/1123/********************************************************************/1124int1125zv_equal0(GEN V)1126{1127long l = lg(V);1128while (--l > 0)1129if (V[l]) return 0;1130return 1;1131}11321133int1134ZV_equal0(GEN V)1135{1136long l = lg(V);1137while (--l > 0)1138if (signe(gel(V,l))) return 0;1139return 1;1140}1141int1142ZMrow_equal0(GEN V, long i)1143{1144long l = lg(V);1145while (--l > 0)1146if (signe(gcoeff(V,i,l))) return 0;1147return 1;1148}11491150static int1151ZV_equal_lg(GEN V, GEN W, long l)1152{1153while (--l > 0)1154if (!equalii(gel(V,l), gel(W,l))) return 0;1155return 1;1156}1157int1158ZV_equal(GEN V, GEN W)1159{1160long l = lg(V);1161if (lg(W) != l) return 0;1162return ZV_equal_lg(V, W, l);1163}1164int1165ZM_equal(GEN A, GEN B)1166{1167long i, m, l = lg(A);1168if (lg(B) != l) return 0;1169if (l == 1) return 1;1170m = lgcols(A);1171if (lgcols(B) != m) return 0;1172for (i = 1; i < l; i++)1173if (!ZV_equal_lg(gel(A,i), gel(B,i), m)) return 0;1174return 1;1175}1176int1177ZM_equal0(GEN A)1178{1179long i, j, m, l = lg(A);1180if (l == 1) return 1;1181m = lgcols(A);1182for (j = 1; j < l; j++)1183for (i = 1; i < m; i++)1184if (signe(gcoeff(A,i,j))) return 0;1185return 1;1186}1187int1188zv_equal(GEN V, GEN W)1189{1190long l = lg(V);1191if (lg(W) != l) return 0;1192while (--l > 0)1193if (V[l] != W[l]) return 0;1194return 1;1195}11961197int1198zvV_equal(GEN V, GEN W)1199{1200long l = lg(V);1201if (lg(W) != l) return 0;1202while (--l > 0)1203if (!zv_equal(gel(V,l),gel(W,l))) return 0;1204return 1;1205}12061207int1208ZM_ishnf(GEN x)1209{1210long i,j, lx = lg(x);1211for (i=1; i<lx; i++)1212{1213GEN xii = gcoeff(x,i,i);1214if (signe(xii) <= 0) return 0;1215for (j=1; j<i; j++)1216if (signe(gcoeff(x,i,j))) return 0;1217for (j=i+1; j<lx; j++)1218{1219GEN xij = gcoeff(x,i,j);1220if (signe(xij)<0 || cmpii(xij,xii)>=0) return 0;1221}1222}1223return 1;1224}1225int1226ZM_isidentity(GEN x)1227{1228long i,j, lx = lg(x);12291230if (lx == 1) return 1;1231if (lx != lgcols(x)) return 0;1232for (j=1; j<lx; j++)1233{1234GEN c = gel(x,j);1235for (i=1; i<j; )1236if (signe(gel(c,i++))) return 0;1237/* i = j */1238if (!equali1(gel(c,i++))) return 0;1239for ( ; i<lx; )1240if (signe(gel(c,i++))) return 0;1241}1242return 1;1243}1244int1245ZM_isdiagonal(GEN x)1246{1247long i,j, lx = lg(x);1248if (lx == 1) return 1;1249if (lx != lgcols(x)) return 0;12501251for (j=1; j<lx; j++)1252{1253GEN c = gel(x,j);1254for (i=1; i<j; i++)1255if (signe(gel(c,i))) return 0;1256for (i++; i<lx; i++)1257if (signe(gel(c,i))) return 0;1258}1259return 1;1260}1261int1262ZM_isscalar(GEN x, GEN s)1263{1264long i, j, lx = lg(x);12651266if (lx == 1) return 1;1267if (!s) s = gcoeff(x,1,1);1268if (equali1(s)) return ZM_isidentity(x);1269if (lx != lgcols(x)) return 0;1270for (j=1; j<lx; j++)1271{1272GEN c = gel(x,j);1273for (i=1; i<j; )1274if (signe(gel(c,i++))) return 0;1275/* i = j */1276if (!equalii(gel(c,i++), s)) return 0;1277for ( ; i<lx; )1278if (signe(gel(c,i++))) return 0;1279}1280return 1;1281}12821283long1284ZC_is_ei(GEN x)1285{1286long i, j = 0, l = lg(x);1287for (i = 1; i < l; i++)1288{1289GEN c = gel(x,i);1290long s = signe(c);1291if (!s) continue;1292if (s < 0 || !is_pm1(c) || j) return 0;1293j = i;1294}1295return j;1296}12971298/********************************************************************/1299/** **/1300/** MISCELLANEOUS **/1301/** **/1302/********************************************************************/1303/* assume lg(x) = lg(y), x,y in Z^n */1304int1305ZV_cmp(GEN x, GEN y)1306{1307long fl,i, lx = lg(x);1308for (i=1; i<lx; i++)1309if (( fl = cmpii(gel(x,i), gel(y,i)) )) return fl;1310return 0;1311}1312/* assume lg(x) = lg(y), x,y in Z^n */1313int1314ZV_abscmp(GEN x, GEN y)1315{1316long fl,i, lx = lg(x);1317for (i=1; i<lx; i++)1318if (( fl = abscmpii(gel(x,i), gel(y,i)) )) return fl;1319return 0;1320}13211322long1323zv_content(GEN x)1324{1325long i, s, l = lg(x);1326if (l == 1) return 0;1327s = labs(x[1]);1328for (i = 2; i < l && s != 1; i++) s = ugcd(s, labs(x[i]));1329return s;1330}1331GEN1332ZV_content(GEN x)1333{1334long i, l = lg(x);1335pari_sp av = avma;1336GEN c;1337if (l == 1) return gen_0;1338if (l == 2) return absi(gel(x,1));1339c = gel(x,1);1340for (i = 2; i < l; i++)1341{1342c = gcdii(c, gel(x,i));1343if (is_pm1(c)) { set_avma(av); return gen_1; }1344}1345return gerepileuptoint(av, c);1346}13471348GEN1349ZM_det_triangular(GEN mat)1350{1351pari_sp av;1352long i,l = lg(mat);1353GEN s;13541355if (l<3) return l<2? gen_1: icopy(gcoeff(mat,1,1));1356av = avma; s = gcoeff(mat,1,1);1357for (i=2; i<l; i++) s = mulii(s,gcoeff(mat,i,i));1358return gerepileuptoint(av,s);1359}13601361/* assumes no overflow */1362long1363zv_prod(GEN v)1364{1365long n, i, l = lg(v);1366if (l == 1) return 1;1367n = v[1]; for (i = 2; i < l; i++) n *= v[i];1368return n;1369}13701371static GEN1372_mulii(void *E, GEN a, GEN b)1373{ (void) E; return mulii(a, b); }13741375/* product of ulongs */1376GEN1377zv_prod_Z(GEN v)1378{1379pari_sp av = avma;1380long k, m, n = lg(v)-1;1381GEN V;1382switch(n) {1383case 0: return gen_1;1384case 1: return utoi(v[1]);1385case 2: return muluu(v[1], v[2]);1386}1387m = n >> 1;1388V = cgetg(m + (odd(n)? 2: 1), t_VEC);1389for (k = 1; k <= m; k++) gel(V,k) = muluu(v[k<<1], v[(k<<1)-1]);1390if (odd(n)) gel(V,k) = utoipos(v[n]);1391return gerepileuptoint(av, gen_product(V, NULL, &_mulii));1392}1393GEN1394vecsmall_prod(GEN v)1395{1396pari_sp av = avma;1397long k, m, n = lg(v)-1;1398GEN V;1399switch (n) {1400case 0: return gen_1;1401case 1: return stoi(v[1]);1402case 2: return mulss(v[1], v[2]);1403}1404m = n >> 1;1405V = cgetg(m + (odd(n)? 2: 1), t_VEC);1406for (k = 1; k <= m; k++) gel(V,k) = mulss(v[k<<1], v[(k<<1)-1]);1407if (odd(n)) gel(V,k) = stoi(v[n]);1408return gerepileuptoint(av, gen_product(V, NULL, &_mulii));1409}14101411GEN1412ZV_prod(GEN v)1413{1414pari_sp av = avma;1415long i, l = lg(v);1416GEN n;1417if (l == 1) return gen_1;1418if (l > 7) return gerepileuptoint(av, gen_product(v, NULL, _mulii));1419n = gel(v,1);1420if (l == 2) return icopy(n);1421for (i = 2; i < l; i++) n = mulii(n, gel(v,i));1422return gerepileuptoint(av, n);1423}1424/* assumes no overflow */1425long1426zv_sum(GEN v)1427{1428long n, i, l = lg(v);1429if (l == 1) return 0;1430n = v[1]; for (i = 2; i < l; i++) n += v[i];1431return n;1432}1433/* assumes no overflow and 0 <= n <= #v */1434long1435zv_sumpart(GEN v, long n)1436{1437long i, p;1438if (!n) return 0;1439p = v[1]; for (i = 2; i <= n; i++) p += v[i];1440return p;1441}1442GEN1443ZV_sum(GEN v)1444{1445pari_sp av = avma;1446long i, l = lg(v);1447GEN n;1448if (l == 1) return gen_0;1449n = gel(v,1);1450if (l == 2) return icopy(n);1451for (i = 2; i < l; i++) n = addii(n, gel(v,i));1452return gerepileuptoint(av, n);1453}14541455/********************************************************************/1456/** **/1457/** GRAM SCHMIDT REDUCTION (integer matrices) **/1458/** **/1459/********************************************************************/14601461/* L[k,] += q * L[l,], l < k. Inefficient if q = 0 */1462static void1463Zupdate_row(long k, long l, GEN q, GEN L, GEN B)1464{1465long i, qq = itos_or_0(q);1466if (!qq)1467{1468for(i=1;i<l;i++) gcoeff(L,k,i) = addii(gcoeff(L,k,i),mulii(q,gcoeff(L,l,i)));1469gcoeff(L,k,l) = addii(gcoeff(L,k,l), mulii(q,B));1470return;1471}1472if (qq == 1) {1473for (i=1;i<l; i++) gcoeff(L,k,i) = addii(gcoeff(L,k,i),gcoeff(L,l,i));1474gcoeff(L,k,l) = addii(gcoeff(L,k,l), B);1475} else if (qq == -1) {1476for (i=1;i<l; i++) gcoeff(L,k,i) = subii(gcoeff(L,k,i),gcoeff(L,l,i));1477gcoeff(L,k,l) = addii(gcoeff(L,k,l), negi(B));1478} else {1479for(i=1;i<l;i++) gcoeff(L,k,i) = addii(gcoeff(L,k,i),mulsi(qq,gcoeff(L,l,i)));1480gcoeff(L,k,l) = addii(gcoeff(L,k,l), mulsi(qq,B));1481}1482}14831484/* update L[k,] */1485static void1486ZRED(long k, long l, GEN x, GEN L, GEN B)1487{1488GEN q = truedivii(addii(B,shifti(gcoeff(L,k,l),1)), shifti(B,1));1489if (!signe(q)) return;1490q = negi(q);1491Zupdate_row(k,l,q,L,B);1492gel(x,k) = ZC_lincomb(gen_1, q, gel(x,k), gel(x,l));1493}14941495/* Gram-Schmidt reduction, x a ZM */1496static void1497ZincrementalGS(GEN x, GEN L, GEN B, long k)1498{1499long i, j;1500for (j=1; j<=k; j++)1501{1502pari_sp av = avma;1503GEN u = ZV_dotproduct(gel(x,k), gel(x,j));1504for (i=1; i<j; i++)1505{1506u = subii(mulii(gel(B,i+1), u), mulii(gcoeff(L,k,i), gcoeff(L,j,i)));1507u = diviiexact(u, gel(B,i));1508}1509gcoeff(L,k,j) = gerepileuptoint(av, u);1510}1511gel(B,k+1) = gcoeff(L,k,k); gcoeff(L,k,k) = gen_1;1512}15131514/* Variant reducemodinvertible(ZC v, ZM y), when y singular.1515* Very inefficient if y is not LLL-reduced of maximal rank */1516static GEN1517ZC_reducemodmatrix_i(GEN v, GEN y)1518{1519GEN B, L, x = shallowconcat(y, v);1520long k, lx = lg(x), nx = lx-1;15211522B = scalarcol_shallow(gen_1, lx);1523L = zeromatcopy(nx, nx);1524for (k=1; k <= nx; k++) ZincrementalGS(x, L, B, k);1525for (k = nx-1; k >= 1; k--) ZRED(nx,k, x,L,gel(B,k+1));1526return gel(x,nx);1527}1528GEN1529ZC_reducemodmatrix(GEN v, GEN y) {1530pari_sp av = avma;1531return gerepilecopy(av, ZC_reducemodmatrix_i(v,y));1532}15331534/* Variant reducemodinvertible(ZM v, ZM y), when y singular.1535* Very inefficient if y is not LLL-reduced of maximal rank */1536static GEN1537ZM_reducemodmatrix_i(GEN v, GEN y)1538{1539GEN B, L, V;1540long j, k, lv = lg(v), nx = lg(y), lx = nx+1;15411542V = cgetg(lv, t_MAT);1543B = scalarcol_shallow(gen_1, lx);1544L = zeromatcopy(nx, nx);1545for (k=1; k < nx; k++) ZincrementalGS(y, L, B, k);1546for (j = 1; j < lg(v); j++)1547{1548GEN x = shallowconcat(y, gel(v,j));1549ZincrementalGS(x, L, B, nx); /* overwrite last */1550for (k = nx-1; k >= 1; k--) ZRED(nx,k, x,L,gel(B,k+1));1551gel(V,j) = gel(x,nx);1552}1553return V;1554}1555GEN1556ZM_reducemodmatrix(GEN v, GEN y) {1557pari_sp av = avma;1558return gerepilecopy(av, ZM_reducemodmatrix_i(v,y));1559}15601561GEN1562ZC_reducemodlll(GEN x,GEN y)1563{1564pari_sp av = avma;1565GEN z = ZC_reducemodmatrix(x, ZM_lll(y, 0.75, LLL_INPLACE));1566return gerepilecopy(av, z);1567}1568GEN1569ZM_reducemodlll(GEN x,GEN y)1570{1571pari_sp av = avma;1572GEN z = ZM_reducemodmatrix(x, ZM_lll(y, 0.75, LLL_INPLACE));1573return gerepilecopy(av, z);1574}157515761577