Line data Source code
1 : /* Copyright (C) 2000 The PARI group.
2 :
3 : This file is part of the PARI/GP package.
4 :
5 : PARI/GP is free software; you can redistribute it and/or modify it under the
6 : terms of the GNU General Public License as published by the Free Software
7 : Foundation; either version 2 of the License, or (at your option) any later
8 : version. It is distributed in the hope that it will be useful, but WITHOUT
9 : ANY WARRANTY WHATSOEVER.
10 :
11 : Check the License for details. You should have received a copy of it, along
12 : with the package; see the file 'COPYING'. If not, write to the Free Software
13 : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
14 :
15 : #include "pari.h"
16 : #include "paripriv.h"
17 :
18 : static int
19 1931230 : check_ZV(GEN x, long l)
20 : {
21 : long i;
22 13037423 : for (i=1; i<l; i++)
23 11106249 : if (typ(gel(x,i)) != t_INT) return 0;
24 1931174 : return 1;
25 : }
26 : void
27 1422011 : RgV_check_ZV(GEN A, const char *s)
28 : {
29 1422011 : if (!RgV_is_ZV(A)) pari_err_TYPE(stack_strcat(s," [integer vector]"), A);
30 1422003 : }
31 : void
32 650549 : RgM_check_ZM(GEN A, const char *s)
33 : {
34 650549 : long n = lg(A);
35 650549 : if (n != 1)
36 : {
37 650416 : long j, m = lgcols(A);
38 2581590 : for (j=1; j<n; j++)
39 1931230 : if (!check_ZV(gel(A,j), m))
40 56 : pari_err_TYPE(stack_strcat(s," [integer matrix]"), A);
41 : }
42 650493 : }
43 :
44 : /* assume m > 1 */
45 : static long
46 108656164 : ZV_max_lg_i(GEN x, long m)
47 : {
48 108656164 : long i, l = lgefint(gel(x,1));
49 875436146 : for (i = 2; i < m; i++) l = maxss(l, lgefint(gel(x,i)));
50 108656698 : return l;
51 : }
52 : long
53 10640 : ZV_max_lg(GEN x)
54 : {
55 10640 : long m = lg(x);
56 10640 : return m == 1? 2: ZV_max_lg_i(x, m);
57 : }
58 :
59 : /* assume n > 1 and m > 1 */
60 : static long
61 26762354 : ZM_max_lg_i(GEN x, long n, long m)
62 : {
63 26762354 : long j, l = ZV_max_lg_i(gel(x,1), m);
64 108646405 : for (j = 2; j < n; j++) l = maxss(l, ZV_max_lg_i(gel(x,j), m));
65 26762690 : return l;
66 : }
67 : long
68 22536 : ZM_max_lg(GEN x)
69 : {
70 22536 : long n = lg(x), m;
71 22536 : if (n == 1) return 2;
72 22536 : m = lgcols(x); return m == 1? 2: ZM_max_lg_i(x, n, m);
73 : }
74 :
75 : /* assume m > 1 */
76 : static long
77 0 : ZV_max_expi_i(GEN x, long m)
78 : {
79 0 : long i, prec = expi(gel(x,1));
80 0 : for (i = 2; i < m; i++) prec = maxss(prec, expi(gel(x,i)));
81 0 : return prec;
82 : }
83 : long
84 0 : ZV_max_expi(GEN x)
85 : {
86 0 : long m = lg(x);
87 0 : return m == 1? 2: ZV_max_expi_i(x, m);
88 : }
89 :
90 : /* assume n > 1 and m > 1 */
91 : static long
92 0 : ZM_max_expi_i(GEN x, long n, long m)
93 : {
94 0 : long j, prec = ZV_max_expi_i(gel(x,1), m);
95 0 : for (j = 2; j < n; j++) prec = maxss(prec, ZV_max_expi_i(gel(x,j), m));
96 0 : return prec;
97 : }
98 : long
99 0 : ZM_max_expi(GEN x)
100 : {
101 0 : long n = lg(x), m;
102 0 : if (n == 1) return 2;
103 0 : m = lgcols(x); return m == 1? 2: ZM_max_expi_i(x, n, m);
104 : }
105 :
106 : GEN
107 4246 : ZM_supnorm(GEN x)
108 : {
109 4246 : long i, j, h, lx = lg(x);
110 4246 : GEN s = gen_0;
111 4246 : if (lx == 1) return gen_1;
112 4246 : h = lgcols(x);
113 26144 : for (j=1; j<lx; j++)
114 : {
115 21898 : GEN xj = gel(x,j);
116 294710 : for (i=1; i<h; i++)
117 : {
118 272812 : GEN c = gel(xj,i);
119 272812 : if (abscmpii(c, s) > 0) s = c;
120 : }
121 : }
122 4246 : return absi(s);
123 : }
124 :
125 : /********************************************************************/
126 : /** **/
127 : /** MULTIPLICATION **/
128 : /** **/
129 : /********************************************************************/
130 : /* x nonempty ZM, y a compatible nc (dimension > 0). */
131 : static GEN
132 0 : ZM_nc_mul_i(GEN x, GEN y, long c, long l)
133 : {
134 : long i, j;
135 : pari_sp av;
136 0 : GEN z = cgetg(l,t_COL), s;
137 :
138 0 : for (i=1; i<l; i++)
139 : {
140 0 : av = avma; s = muliu(gcoeff(x,i,1),y[1]);
141 0 : for (j=2; j<c; j++)
142 0 : if (y[j]) s = addii(s, muliu(gcoeff(x,i,j),y[j]));
143 0 : gel(z,i) = gerepileuptoint(av,s);
144 : }
145 0 : return z;
146 : }
147 :
148 : /* x ZV, y a compatible zc. */
149 : GEN
150 2229528 : ZV_zc_mul(GEN x, GEN y)
151 : {
152 2229528 : long j, l = lg(x);
153 2229528 : pari_sp av = avma;
154 2229528 : GEN s = mulis(gel(x,1),y[1]);
155 50292998 : for (j=2; j<l; j++)
156 48063470 : if (y[j]) s = addii(s, mulis(gel(x,j),y[j]));
157 2229528 : return gerepileuptoint(av,s);
158 : }
159 :
160 : /* x nonempty ZM, y a compatible zc (dimension > 0). */
161 : static GEN
162 19278801 : ZM_zc_mul_i(GEN x, GEN y, long c, long l)
163 : {
164 : long i, j;
165 19278801 : GEN z = cgetg(l,t_COL);
166 :
167 121745616 : for (i=1; i<l; i++)
168 : {
169 102467595 : pari_sp av = avma;
170 102467595 : GEN s = mulis(gcoeff(x,i,1),y[1]);
171 1170839459 : for (j=2; j<c; j++)
172 1068364300 : if (y[j]) s = addii(s, mulis(gcoeff(x,i,j),y[j]));
173 102475159 : gel(z,i) = gerepileuptoint(av,s);
174 : }
175 19278021 : return z;
176 : }
177 : GEN
178 17214554 : ZM_zc_mul(GEN x, GEN y) {
179 17214554 : long l = lg(x);
180 17214554 : if (l == 1) return cgetg(1, t_COL);
181 17214554 : return ZM_zc_mul_i(x,y, l, lgcols(x));
182 : }
183 :
184 : /* y nonempty ZM, x a compatible zv (dimension > 0). */
185 : GEN
186 1736 : zv_ZM_mul(GEN x, GEN y) {
187 1736 : long i,j, lx = lg(x), ly = lg(y);
188 : GEN z;
189 1736 : if (lx == 1) return zerovec(ly-1);
190 1736 : z = cgetg(ly,t_VEC);
191 4046 : for (j=1; j<ly; j++)
192 : {
193 2310 : pari_sp av = avma;
194 2310 : GEN s = mulsi(x[1], gcoeff(y,1,j));
195 3990 : for (i=2; i<lx; i++)
196 1680 : if (x[i]) s = addii(s, mulsi(x[i], gcoeff(y,i,j)));
197 2310 : gel(z,j) = gerepileuptoint(av,s);
198 : }
199 1736 : return z;
200 : }
201 :
202 : /* x ZM, y a compatible zm (dimension > 0). */
203 : GEN
204 975604 : ZM_zm_mul(GEN x, GEN y)
205 : {
206 975604 : long j, c, l = lg(x), ly = lg(y);
207 975604 : GEN z = cgetg(ly, t_MAT);
208 975603 : if (l == 1) return z;
209 975596 : c = lgcols(x);
210 3039853 : for (j = 1; j < ly; j++) gel(z,j) = ZM_zc_mul_i(x, gel(y,j), l,c);
211 975596 : return z;
212 : }
213 : /* x ZM, y a compatible zn (dimension > 0). */
214 : GEN
215 0 : ZM_nm_mul(GEN x, GEN y)
216 : {
217 0 : long j, c, l = lg(x), ly = lg(y);
218 0 : GEN z = cgetg(ly, t_MAT);
219 0 : if (l == 1) return z;
220 0 : c = lgcols(x);
221 0 : for (j = 1; j < ly; j++) gel(z,j) = ZM_nc_mul_i(x, gel(y,j), l,c);
222 0 : return z;
223 : }
224 :
225 : /* Strassen-Winograd algorithm */
226 :
227 : /* Return A[ma+1..ma+da, na+1..na+ea] - B[mb+1..mb+db, nb+1..nb+eb]
228 : * as an (m x n)-matrix, padding the input with zeroes as necessary. */
229 : static GEN
230 593352 : add_slices(long m, long n,
231 : GEN A, long ma, long da, long na, long ea,
232 : GEN B, long mb, long db, long nb, long eb)
233 : {
234 593352 : long min_d = minss(da, db), min_e = minss(ea, eb), i, j;
235 593352 : GEN M = cgetg(n + 1, t_MAT), C;
236 :
237 5019733 : for (j = 1; j <= min_e; j++) {
238 4426381 : gel(M, j) = C = cgetg(m + 1, t_COL);
239 86173582 : for (i = 1; i <= min_d; i++)
240 81747201 : gel(C, i) = addii(gcoeff(A, ma + i, na + j),
241 81747201 : gcoeff(B, mb + i, nb + j));
242 4492759 : for (; i <= da; i++) gel(C, i) = gcoeff(A, ma + i, na + j);
243 4426381 : for (; i <= db; i++) gel(C, i) = gcoeff(B, mb + i, nb + j);
244 4426381 : for (; i <= m; i++) gel(C, i) = gen_0;
245 : }
246 652910 : for (; j <= ea; j++) {
247 59558 : gel(M, j) = C = cgetg(m + 1, t_COL);
248 208567 : for (i = 1; i <= da; i++) gel(C, i) = gcoeff(A, ma + i, na + j);
249 59558 : for (; i <= m; i++) gel(C, i) = gen_0;
250 : }
251 593352 : for (; j <= eb; j++) {
252 0 : gel(M, j) = C = cgetg(m + 1, t_COL);
253 0 : for (i = 1; i <= db; i++) gel(C, i) = gcoeff(B, mb + i, nb + j);
254 0 : for (; i <= m; i++) gel(C, i) = gen_0;
255 : }
256 593352 : for (; j <= n; j++) gel(M, j) = zerocol(m);
257 593352 : return M;
258 : }
259 :
260 : /* Return A[ma+1..ma+da, na+1..na+ea] - B[mb+1..mb+db, nb+1..nb+eb]
261 : * as an (m x n)-matrix, padding the input with zeroes as necessary. */
262 : static GEN
263 519183 : subtract_slices(long m, long n,
264 : GEN A, long ma, long da, long na, long ea,
265 : GEN B, long mb, long db, long nb, long eb)
266 : {
267 519183 : long min_d = minss(da, db), min_e = minss(ea, eb), i, j;
268 519183 : GEN M = cgetg(n + 1, t_MAT), C;
269 :
270 4408587 : for (j = 1; j <= min_e; j++) {
271 3889404 : gel(M, j) = C = cgetg(m + 1, t_COL);
272 77110189 : for (i = 1; i <= min_d; i++)
273 73220785 : gel(C, i) = subii(gcoeff(A, ma + i, na + j),
274 73220785 : gcoeff(B, mb + i, nb + j));
275 3946514 : for (; i <= da; i++) gel(C, i) = gcoeff(A, ma + i, na + j);
276 3977547 : for (; i <= db; i++) gel(C, i) = negi(gcoeff(B, mb + i, nb + j));
277 3889404 : for (; i <= m; i++) gel(C, i) = gen_0;
278 : }
279 519183 : for (; j <= ea; j++) {
280 0 : gel(M, j) = C = cgetg(m + 1, t_COL);
281 0 : for (i = 1; i <= da; i++) gel(C, i) = gcoeff(A, ma + i, na + j);
282 0 : for (; i <= m; i++) gel(C, i) = gen_0;
283 : }
284 554253 : for (; j <= eb; j++) {
285 35070 : gel(M, j) = C = cgetg(m + 1, t_COL);
286 127205 : for (i = 1; i <= db; i++) gel(C, i) = negi(gcoeff(B, mb + i, nb + j));
287 35070 : for (; i <= m; i++) gel(C, i) = gen_0;
288 : }
289 554253 : for (; j <= n; j++) gel(M, j) = zerocol(m);
290 519183 : return M;
291 : }
292 :
293 : static GEN ZM_mul_i(GEN x, GEN y, long l, long lx, long ly);
294 :
295 : /* Strassen-Winograd matrix product A (m x n) * B (n x p) */
296 : static GEN
297 74169 : ZM_mul_sw(GEN A, GEN B, long m, long n, long p)
298 : {
299 74169 : pari_sp av = avma;
300 74169 : long m1 = (m + 1)/2, m2 = m/2,
301 74169 : n1 = (n + 1)/2, n2 = n/2,
302 74169 : p1 = (p + 1)/2, p2 = p/2;
303 : GEN A11, A12, A22, B11, B21, B22,
304 : S1, S2, S3, S4, T1, T2, T3, T4,
305 : M1, M2, M3, M4, M5, M6, M7,
306 : V1, V2, V3, C11, C12, C21, C22, C;
307 :
308 74169 : T2 = subtract_slices(n1, p2, B, 0, n1, p1, p2, B, n1, n2, p1, p2);
309 74169 : S1 = subtract_slices(m2, n1, A, m1, m2, 0, n1, A, 0, m2, 0, n1);
310 74169 : M2 = ZM_mul_i(S1, T2, m2 + 1, n1 + 1, p2 + 1);
311 74169 : if (gc_needed(av, 1))
312 0 : gerepileall(av, 2, &T2, &M2); /* destroy S1 */
313 74169 : T3 = subtract_slices(n1, p1, T2, 0, n1, 0, p2, B, 0, n1, 0, p1);
314 74169 : if (gc_needed(av, 1))
315 0 : gerepileall(av, 2, &M2, &T3); /* destroy T2 */
316 74169 : S2 = add_slices(m2, n1, A, m1, m2, 0, n1, A, m1, m2, n1, n2);
317 74169 : T1 = subtract_slices(n1, p1, B, 0, n1, p1, p2, B, 0, n1, 0, p2);
318 74169 : M3 = ZM_mul_i(S2, T1, m2 + 1, n1 + 1, p2 + 1);
319 74169 : if (gc_needed(av, 1))
320 0 : gerepileall(av, 4, &M2, &T3, &S2, &M3); /* destroy T1 */
321 74169 : S3 = subtract_slices(m1, n1, S2, 0, m2, 0, n1, A, 0, m1, 0, n1);
322 74169 : if (gc_needed(av, 1))
323 0 : gerepileall(av, 4, &M2, &T3, &M3, &S3); /* destroy S2 */
324 74169 : A11 = matslice(A, 1, m1, 1, n1);
325 74169 : B11 = matslice(B, 1, n1, 1, p1);
326 74169 : M1 = ZM_mul_i(A11, B11, m1 + 1, n1 + 1, p1 + 1);
327 74169 : if (gc_needed(av, 1))
328 0 : gerepileall(av, 5, &M2, &T3, &M3, &S3, &M1); /* destroy A11, B11 */
329 74169 : A12 = matslice(A, 1, m1, n1 + 1, n);
330 74169 : B21 = matslice(B, n1 + 1, n, 1, p1);
331 74169 : M4 = ZM_mul_i(A12, B21, m1 + 1, n2 + 1, p1 + 1);
332 74169 : if (gc_needed(av, 1))
333 0 : gerepileall(av, 6, &M2, &T3, &M3, &S3, &M1, &M4); /* destroy A12, B21 */
334 74169 : C11 = add_slices(m1, p1, M1, 0, m1, 0, p1, M4, 0, m1, 0, p1);
335 74169 : if (gc_needed(av, 1))
336 0 : gerepileall(av, 6, &M2, &T3, &M3, &S3, &M1, &C11); /* destroy M4 */
337 74169 : M5 = ZM_mul_i(S3, T3, m1 + 1, n1 + 1, p1 + 1);
338 74169 : S4 = subtract_slices(m1, n2, A, 0, m1, n1, n2, S3, 0, m1, 0, n2);
339 74169 : if (gc_needed(av, 1))
340 5 : gerepileall(av, 7, &M2, &T3, &M3, &M1, &C11, &M5, &S4); /* destroy S3 */
341 74169 : T4 = add_slices(n2, p1, B, n1, n2, 0, p1, T3, 0, n2, 0, p1);
342 74169 : if (gc_needed(av, 1))
343 0 : gerepileall(av, 7, &M2, &M3, &M1, &C11, &M5, &S4, &T4); /* destroy T3 */
344 74169 : V1 = subtract_slices(m1, p1, M1, 0, m1, 0, p1, M5, 0, m1, 0, p1);
345 74169 : if (gc_needed(av, 1))
346 1 : gerepileall(av, 6, &M2, &M3, &S4, &T4, &C11, &V1); /* destroy M1, M5 */
347 74169 : B22 = matslice(B, n1 + 1, n, p1 + 1, p);
348 74169 : M6 = ZM_mul_i(S4, B22, m1 + 1, n2 + 1, p2 + 1);
349 74169 : if (gc_needed(av, 1))
350 6 : gerepileall(av, 6, &M2, &M3, &T4, &C11, &V1, &M6); /* destroy S4, B22 */
351 74169 : A22 = matslice(A, m1 + 1, m, n1 + 1, n);
352 74169 : M7 = ZM_mul_i(A22, T4, m2 + 1, n2 + 1, p1 + 1);
353 74169 : if (gc_needed(av, 1))
354 0 : gerepileall(av, 6, &M2, &M3, &C11, &V1, &M6, &M7); /* destroy A22, T4 */
355 74169 : V3 = add_slices(m1, p2, V1, 0, m1, 0, p2, M3, 0, m2, 0, p2);
356 74169 : C12 = add_slices(m1, p2, V3, 0, m1, 0, p2, M6, 0, m1, 0, p2);
357 74169 : if (gc_needed(av, 1))
358 6 : gerepileall(av, 6, &M2, &M3, &C11, &V1, &M7, &C12); /* destroy V3, M6 */
359 74169 : V2 = add_slices(m2, p1, V1, 0, m2, 0, p1, M2, 0, m2, 0, p2);
360 74169 : if (gc_needed(av, 1))
361 0 : gerepileall(av, 5, &M3, &C11, &M7, &C12, &V2); /* destroy V1, M2 */
362 74169 : C21 = add_slices(m2, p1, V2, 0, m2, 0, p1, M7, 0, m2, 0, p1);
363 74169 : if (gc_needed(av, 1))
364 6 : gerepileall(av, 5, &M3, &C11, &C12, &V2, &C21); /* destroy M7 */
365 74169 : C22 = add_slices(m2, p2, V2, 0, m2, 0, p2, M3, 0, m2, 0, p2);
366 74169 : if (gc_needed(av, 1))
367 0 : gerepileall(av, 4, &C11, &C12, &C21, &C22); /* destroy V2, M3 */
368 74169 : C = shallowconcat(vconcat(C11, C21), vconcat(C12, C22));
369 74169 : return gerepilecopy(av, C);
370 : }
371 :
372 : /* x[i,]*y. Assume lg(x) > 1 and 0 < i < lgcols(x) */
373 : static GEN
374 577766077 : ZMrow_ZC_mul_i(GEN x, GEN y, long i, long lx)
375 : {
376 577766077 : pari_sp av = avma;
377 577766077 : GEN c = mulii(gcoeff(x,i,1), gel(y,1)), ZERO = gen_0;
378 : long k;
379 6692079327 : for (k = 2; k < lx; k++)
380 : {
381 6114697869 : GEN t = mulii(gcoeff(x,i,k), gel(y,k));
382 6113139826 : if (t != ZERO) c = addii(c, t);
383 : }
384 577381458 : return gerepileuptoint(av, c);
385 : }
386 : GEN
387 135180938 : ZMrow_ZC_mul(GEN x, GEN y, long i)
388 135180938 : { return ZMrow_ZC_mul_i(x, y, i, lg(x)); }
389 :
390 : /* return x * y, 1 < lx = lg(x), l = lgcols(x) */
391 : static GEN
392 74765357 : ZM_ZC_mul_i(GEN x, GEN y, long lx, long l)
393 : {
394 74765357 : GEN z = cgetg(l,t_COL);
395 : long i;
396 517370872 : for (i=1; i<l; i++) gel(z,i) = ZMrow_ZC_mul_i(x,y,i,lx);
397 74724053 : return z;
398 : }
399 :
400 : static GEN
401 13306967 : ZM_mul_classical(GEN x, GEN y, long l, long lx, long ly)
402 : {
403 : long j;
404 13306967 : GEN z = cgetg(ly, t_MAT);
405 64019860 : for (j = 1; j < ly; j++)
406 50715562 : gel(z, j) = ZM_ZC_mul_i(x, gel(y, j), lx, l);
407 13304298 : return z;
408 : }
409 :
410 : static GEN
411 1105 : ZM_mul_slice(GEN A, GEN B, GEN P, GEN *mod)
412 : {
413 1105 : pari_sp av = avma;
414 1105 : long i, n = lg(P)-1;
415 : GEN H, T;
416 1105 : if (n == 1)
417 : {
418 0 : ulong p = uel(P,1);
419 0 : GEN a = ZM_to_Flm(A, p);
420 0 : GEN b = ZM_to_Flm(B, p);
421 0 : GEN Hp = gerepileupto(av, Flm_to_ZM(Flm_mul(a, b, p)));
422 0 : *mod = utoi(p); return Hp;
423 : }
424 1105 : T = ZV_producttree(P);
425 1105 : A = ZM_nv_mod_tree(A, P, T);
426 1105 : B = ZM_nv_mod_tree(B, P, T);
427 1105 : H = cgetg(n+1, t_VEC);
428 7878 : for(i=1; i <= n; i++)
429 6773 : gel(H,i) = Flm_mul(gel(A,i),gel(B,i),P[i]);
430 1105 : H = nmV_chinese_center_tree_seq(H, P, T, ZV_chinesetree(P, T));
431 1105 : *mod = gmael(T, lg(T)-1, 1); return gc_all(av, 2, &H, mod);
432 : }
433 :
434 : GEN
435 1105 : ZM_mul_worker(GEN P, GEN A, GEN B)
436 : {
437 1105 : GEN V = cgetg(3, t_VEC);
438 1105 : gel(V,1) = ZM_mul_slice(A, B, P, &gel(V,2));
439 1105 : return V;
440 : }
441 :
442 : static GEN
443 839 : ZM_mul_fast(GEN A, GEN B, long lA, long lB, long sA, long sB)
444 : {
445 839 : pari_sp av = avma;
446 : forprime_t S;
447 : GEN H, worker;
448 : long h;
449 839 : if (sA == 2 || sB == 2) return zeromat(nbrows(A),lB-1);
450 827 : h = 1 + (sA + sB - 4) * BITS_IN_LONG + expu(lA-1);
451 827 : init_modular_big(&S);
452 827 : worker = snm_closure(is_entry("_ZM_mul_worker"), mkvec2(A,B));
453 827 : H = gen_crt("ZM_mul", worker, &S, NULL, h, 0, NULL,
454 : nmV_chinese_center, FpM_center);
455 827 : return gerepileupto(av, H);
456 : }
457 :
458 : /* s = min(log_BIL |x|, log_BIL |y|), use Strassen-Winograd when
459 : * min(dims) > B */
460 : static long
461 13381153 : sw_bound(long s)
462 13381153 : { return s > 60 ? 2: s > 25 ? 4: s > 15 ? 8 : s > 8 ? 16 : 32; }
463 :
464 : /* assume lx > 1 and ly > 1; x is (l-1) x (lx-1), y is (lx-1) x (ly-1).
465 : * Return x * y */
466 : static GEN
467 20688128 : ZM_mul_i(GEN x, GEN y, long l, long lx, long ly)
468 : {
469 : long sx, sy, B;
470 : #ifdef LONG_IS_64BIT /* From Flm_mul_i */
471 17558232 : long Flm_sw_bound = 70;
472 : #else
473 3129896 : long Flm_sw_bound = 120;
474 : #endif
475 20688128 : if (l == 1) return zeromat(0, ly-1);
476 20686242 : if (lx==2 && l==2 && ly==2)
477 359530 : { retmkmat(mkcol(mulii(gcoeff(x,1,1), gcoeff(y,1,1)))); }
478 20326712 : if (lx==3 && l==3 && ly==3) return ZM2_mul(x, y);
479 13357501 : sx = ZM_max_lg_i(x, lx, l);
480 13358331 : sy = ZM_max_lg_i(y, ly, lx);
481 : /* Use modular reconstruction if Flm_mul would use Strassen and the input
482 : * sizes look balanced */
483 13358379 : if (lx > Flm_sw_bound && ly > Flm_sw_bound && l > Flm_sw_bound
484 851 : && sx <= 10 * sy && sy <= 10 * sx) return ZM_mul_fast(x,y, lx,ly, sx,sy);
485 :
486 13357540 : B = sw_bound(minss(sx, sy));
487 13357518 : if (l <= B || lx <= B || ly <= B)
488 13283379 : return ZM_mul_classical(x, y, l, lx, ly);
489 : else
490 74139 : return ZM_mul_sw(x, y, l - 1, lx - 1, ly - 1);
491 : }
492 :
493 : GEN
494 20305699 : ZM_mul(GEN x, GEN y)
495 : {
496 20305699 : long lx = lg(x), ly = lg(y);
497 20305699 : if (ly == 1) return cgetg(1,t_MAT);
498 20170512 : if (lx == 1) return zeromat(0, ly-1);
499 20168930 : return ZM_mul_i(x, y, lgcols(x), lx, ly);
500 : }
501 :
502 : static GEN
503 0 : ZM_sqr_slice(GEN A, GEN P, GEN *mod)
504 : {
505 0 : pari_sp av = avma;
506 0 : long i, n = lg(P)-1;
507 : GEN H, T;
508 0 : if (n == 1)
509 : {
510 0 : ulong p = uel(P,1);
511 0 : GEN a = ZM_to_Flm(A, p);
512 0 : GEN Hp = gerepileupto(av, Flm_to_ZM(Flm_sqr(a, p)));
513 0 : *mod = utoi(p); return Hp;
514 : }
515 0 : T = ZV_producttree(P);
516 0 : A = ZM_nv_mod_tree(A, P, T);
517 0 : H = cgetg(n+1, t_VEC);
518 0 : for(i=1; i <= n; i++)
519 0 : gel(H,i) = Flm_sqr(gel(A,i), P[i]);
520 0 : H = nmV_chinese_center_tree_seq(H, P, T, ZV_chinesetree(P, T));
521 0 : *mod = gmael(T, lg(T)-1, 1); return gc_all(av, 2, &H, mod);
522 : }
523 :
524 : GEN
525 0 : ZM_sqr_worker(GEN P, GEN A)
526 : {
527 0 : GEN V = cgetg(3, t_VEC);
528 0 : gel(V,1) = ZM_sqr_slice(A, P, &gel(V,2));
529 0 : return V;
530 : }
531 :
532 : static GEN
533 0 : ZM_sqr_fast(GEN A, long l, long s)
534 : {
535 0 : pari_sp av = avma;
536 : forprime_t S;
537 : GEN H, worker;
538 : long h;
539 0 : if (s == 2) return zeromat(l-1,l-1);
540 0 : h = 1 + (2*s - 4) * BITS_IN_LONG + expu(l-1);
541 0 : init_modular_big(&S);
542 0 : worker = snm_closure(is_entry("_ZM_sqr_worker"), mkvec(A));
543 0 : H = gen_crt("ZM_sqr", worker, &S, NULL, h, 0, NULL,
544 : nmV_chinese_center, FpM_center);
545 0 : return gerepileupto(av, H);
546 : }
547 :
548 : GEN
549 558609 : QM_mul(GEN x, GEN y)
550 : {
551 558609 : GEN dx, nx = Q_primitive_part(x, &dx);
552 558609 : GEN dy, ny = Q_primitive_part(y, &dy);
553 558609 : GEN z = ZM_mul(nx, ny);
554 558609 : if (dx || dy)
555 : {
556 474494 : GEN d = dx ? dy ? gmul(dx, dy): dx : dy;
557 474494 : if (!gequal1(d)) z = ZM_Q_mul(z, d);
558 : }
559 558609 : return z;
560 : }
561 :
562 : GEN
563 700 : QM_sqr(GEN x)
564 : {
565 700 : GEN dx, nx = Q_primitive_part(x, &dx);
566 700 : GEN z = ZM_sqr(nx);
567 700 : if (dx)
568 700 : z = ZM_Q_mul(z, gsqr(dx));
569 700 : return z;
570 : }
571 :
572 : GEN
573 126853 : QM_QC_mul(GEN x, GEN y)
574 : {
575 126853 : GEN dx, nx = Q_primitive_part(x, &dx);
576 126853 : GEN dy, ny = Q_primitive_part(y, &dy);
577 126853 : GEN z = ZM_ZC_mul(nx, ny);
578 126853 : if (dx || dy)
579 : {
580 126832 : GEN d = dx ? dy ? gmul(dx, dy): dx : dy;
581 126832 : if (!gequal1(d)) z = ZC_Q_mul(z, d);
582 : }
583 126853 : return z;
584 : }
585 :
586 : /* assume result is symmetric */
587 : GEN
588 0 : ZM_multosym(GEN x, GEN y)
589 : {
590 0 : long j, lx, ly = lg(y);
591 : GEN M;
592 0 : if (ly == 1) return cgetg(1,t_MAT);
593 0 : lx = lg(x); /* = lgcols(y) */
594 0 : if (lx == 1) return cgetg(1,t_MAT);
595 : /* ly = lgcols(x) */
596 0 : M = cgetg(ly, t_MAT);
597 0 : for (j=1; j<ly; j++)
598 : {
599 0 : GEN z = cgetg(ly,t_COL), yj = gel(y,j);
600 : long i;
601 0 : for (i=1; i<j; i++) gel(z,i) = gcoeff(M,j,i);
602 0 : for (i=j; i<ly; i++)gel(z,i) = ZMrow_ZC_mul_i(x,yj,i,lx);
603 0 : gel(M,j) = z;
604 : }
605 0 : return M;
606 : }
607 :
608 : /* compute m*diagonal(d), assume lg(d) = lg(m). Shallow */
609 : GEN
610 21 : ZM_mul_diag(GEN m, GEN d)
611 : {
612 : long j, l;
613 21 : GEN y = cgetg_copy(m, &l);
614 77 : for (j=1; j<l; j++)
615 : {
616 56 : GEN c = gel(d,j);
617 56 : gel(y,j) = equali1(c)? gel(m,j): ZC_Z_mul(gel(m,j), c);
618 : }
619 21 : return y;
620 : }
621 : /* compute diagonal(d)*m, assume lg(d) = lg(m~). Shallow */
622 : GEN
623 534028 : ZM_diag_mul(GEN d, GEN m)
624 : {
625 534028 : long i, j, l = lg(d), lm = lg(m);
626 534028 : GEN y = cgetg(lm, t_MAT);
627 1521600 : for (j=1; j<lm; j++) gel(y,j) = cgetg(l, t_COL);
628 1648123 : for (i=1; i<l; i++)
629 : {
630 1114404 : GEN c = gel(d,i);
631 1114404 : if (equali1(c))
632 236838 : for (j=1; j<lm; j++) gcoeff(y,i,j) = gcoeff(m,i,j);
633 : else
634 3400628 : for (j=1; j<lm; j++) gcoeff(y,i,j) = mulii(gcoeff(m,i,j), c);
635 : }
636 533719 : return y;
637 : }
638 :
639 : /* assume lx > 1 is lg(x) = lg(y) */
640 : static GEN
641 19188582 : ZV_dotproduct_i(GEN x, GEN y, long lx)
642 : {
643 19188582 : pari_sp av = avma;
644 19188582 : GEN c = mulii(gel(x,1), gel(y,1));
645 : long i;
646 144271731 : for (i = 2; i < lx; i++)
647 : {
648 125084378 : GEN t = mulii(gel(x,i), gel(y,i));
649 125083852 : if (t != gen_0) c = addii(c, t);
650 : }
651 19187353 : return gerepileuptoint(av, c);
652 : }
653 :
654 : /* x~ * y, assuming result is symmetric */
655 : GEN
656 532221 : ZM_transmultosym(GEN x, GEN y)
657 : {
658 532221 : long i, j, l, ly = lg(y);
659 : GEN M;
660 532221 : if (ly == 1) return cgetg(1,t_MAT);
661 : /* lg(x) = ly */
662 532221 : l = lgcols(y); /* = lgcols(x) */
663 532221 : M = cgetg(ly, t_MAT);
664 2708730 : for (i=1; i<ly; i++)
665 : {
666 2176509 : GEN xi = gel(x,i), c = cgetg(ly,t_COL);
667 2176509 : gel(M,i) = c;
668 7109000 : for (j=1; j<i; j++)
669 4932491 : gcoeff(M,i,j) = gel(c,j) = ZV_dotproduct_i(xi,gel(y,j),l);
670 2176509 : gel(c,i) = ZV_dotproduct_i(xi,gel(y,i),l);
671 : }
672 532221 : return M;
673 : }
674 : /* x~ * y */
675 : GEN
676 2289 : ZM_transmul(GEN x, GEN y)
677 : {
678 2289 : long i, j, l, lx, ly = lg(y);
679 : GEN M;
680 2289 : if (ly == 1) return cgetg(1,t_MAT);
681 2289 : lx = lg(x);
682 2289 : l = lgcols(y);
683 2289 : if (lgcols(x) != l) pari_err_OP("operation 'ZM_transmul'", x,y);
684 2289 : M = cgetg(ly, t_MAT);
685 6993 : for (i=1; i<ly; i++)
686 : {
687 4704 : GEN yi = gel(y,i), c = cgetg(lx,t_COL);
688 4704 : gel(M,i) = c;
689 12229 : for (j=1; j<lx; j++) gel(c,j) = ZV_dotproduct_i(yi,gel(x,j),l);
690 : }
691 2289 : return M;
692 : }
693 :
694 : /* assume l > 1; x is (l-1) x (l-1), return x^2.
695 : * FIXME: we ultimately rely on Strassen-Winograd which uses 7M + 15A.
696 : * Should use Bodrato's variant of Winograd, using 3M + 4S + 11A */
697 : static GEN
698 24402 : ZM_sqr_i(GEN x, long l)
699 : {
700 : long s;
701 24402 : if (l == 2) { retmkmat(mkcol(sqri(gcoeff(x,1,1)))); }
702 24402 : if (l == 3) return ZM2_sqr(x);
703 23625 : s = ZM_max_lg_i(x, l, l);
704 23625 : if (l > 70) return ZM_sqr_fast(x, l, s);
705 23625 : if (l <= sw_bound(s))
706 23595 : return ZM_mul_classical(x, x, l, l, l);
707 : else
708 30 : return ZM_mul_sw(x, x, l - 1, l - 1, l - 1);
709 : }
710 :
711 : GEN
712 24402 : ZM_sqr(GEN x)
713 : {
714 24402 : long lx=lg(x);
715 24402 : if (lx==1) return cgetg(1,t_MAT);
716 24402 : return ZM_sqr_i(x, lx);
717 : }
718 : GEN
719 24138569 : ZM_ZC_mul(GEN x, GEN y)
720 : {
721 24138569 : long lx = lg(x);
722 24138569 : return lx==1? cgetg(1,t_COL): ZM_ZC_mul_i(x, y, lx, lgcols(x));
723 : }
724 :
725 : GEN
726 3675996 : ZC_Z_div(GEN x, GEN c)
727 17413454 : { pari_APPLY_type(t_COL, Qdivii(gel(x,i), c)) }
728 :
729 : GEN
730 19308 : ZM_Z_div(GEN x, GEN c)
731 221027 : { pari_APPLY_same(ZC_Z_div(gel(x, i), c)) }
732 :
733 : GEN
734 2726897 : ZC_Q_mul(GEN A, GEN z)
735 : {
736 2726897 : pari_sp av = avma;
737 2726897 : long i, l = lg(A);
738 : GEN d, n, Ad, B, u;
739 2726897 : if (typ(z)==t_INT) return ZC_Z_mul(A,z);
740 2722522 : n = gel(z, 1); d = gel(z, 2);
741 2722522 : Ad = FpC_red(A, d);
742 2722452 : u = gcdii(d, FpV_factorback(Ad, NULL, d));
743 2722443 : B = cgetg(l, t_COL);
744 2722437 : if (equali1(u))
745 : {
746 416942 : for(i=1; i<l; i++)
747 352920 : gel(B, i) = mkfrac(mulii(n, gel(A,i)), d);
748 : } else
749 : {
750 18817349 : for(i=1; i<l; i++)
751 : {
752 16159004 : GEN di = gcdii(gel(Ad, i), u), ni = mulii(n, diviiexact(gel(A,i), di));
753 16158917 : if (equalii(d, di))
754 11188954 : gel(B, i) = ni;
755 : else
756 4969945 : gel(B, i) = mkfrac(ni, diviiexact(d, di));
757 : }
758 : }
759 2722367 : return gerepilecopy(av, B);
760 : }
761 :
762 : GEN
763 1091922 : ZM_Q_mul(GEN x, GEN z)
764 : {
765 1091922 : if (typ(z)==t_INT) return ZM_Z_mul(x,z);
766 3259581 : pari_APPLY_same(ZC_Q_mul(gel(x, i), z));
767 : }
768 :
769 : long
770 196399119 : zv_dotproduct(GEN x, GEN y)
771 : {
772 196399119 : long i, lx = lg(x);
773 : ulong c;
774 196399119 : if (lx == 1) return 0;
775 196399119 : c = uel(x,1)*uel(y,1);
776 3047359154 : for (i = 2; i < lx; i++)
777 2850960035 : c += uel(x,i)*uel(y,i);
778 196399119 : return c;
779 : }
780 :
781 : GEN
782 230629 : ZV_ZM_mul(GEN x, GEN y)
783 : {
784 230629 : long i, lx = lg(x), ly = lg(y);
785 : GEN z;
786 230629 : if (lx == 1) return zerovec(ly-1);
787 230510 : z = cgetg(ly, t_VEC);
788 882063 : for (i = 1; i < ly; i++) gel(z,i) = ZV_dotproduct_i(x, gel(y,i), lx);
789 230510 : return z;
790 : }
791 :
792 : GEN
793 0 : ZC_ZV_mul(GEN x, GEN y)
794 : {
795 0 : long i,j, lx=lg(x), ly=lg(y);
796 : GEN z;
797 0 : if (ly==1) return cgetg(1,t_MAT);
798 0 : z = cgetg(ly,t_MAT);
799 0 : for (j=1; j < ly; j++)
800 : {
801 0 : gel(z,j) = cgetg(lx,t_COL);
802 0 : for (i=1; i<lx; i++) gcoeff(z,i,j) = mulii(gel(x,i),gel(y,j));
803 : }
804 0 : return z;
805 : }
806 :
807 : GEN
808 6662547 : ZV_dotsquare(GEN x)
809 : {
810 : long i, lx;
811 : pari_sp av;
812 : GEN z;
813 6662547 : lx = lg(x);
814 6662547 : if (lx == 1) return gen_0;
815 6662547 : av = avma; z = sqri(gel(x,1));
816 26379044 : for (i=2; i<lx; i++) z = addii(z, sqri(gel(x,i)));
817 6660815 : return gerepileuptoint(av,z);
818 : }
819 :
820 : GEN
821 16202468 : ZV_dotproduct(GEN x,GEN y)
822 : {
823 : long lx;
824 16202468 : if (x == y) return ZV_dotsquare(x);
825 11419720 : lx = lg(x);
826 11419720 : if (lx == 1) return gen_0;
827 11419720 : return ZV_dotproduct_i(x, y, lx);
828 : }
829 :
830 : static GEN
831 280 : _ZM_mul(void *data /*ignored*/, GEN x, GEN y)
832 280 : { (void)data; return ZM_mul(x,y); }
833 : static GEN
834 23331 : _ZM_sqr(void *data /*ignored*/, GEN x)
835 23331 : { (void)data; return ZM_sqr(x); }
836 : /* FIXME: Using Bodrato's squaring, precomputations attached to fixed
837 : * multiplicand should be reused. And some postcomputations can be fused */
838 : GEN
839 0 : ZM_pow(GEN x, GEN n)
840 : {
841 0 : pari_sp av = avma;
842 0 : if (!signe(n)) return matid(lg(x)-1);
843 0 : return gerepilecopy(av, gen_pow_i(x, n, NULL, &_ZM_sqr, &_ZM_mul));
844 : }
845 : GEN
846 22806 : ZM_powu(GEN x, ulong n)
847 : {
848 22806 : pari_sp av = avma;
849 22806 : if (!n) return matid(lg(x)-1);
850 22806 : return gerepilecopy(av, gen_powu_i(x, n, NULL, &_ZM_sqr, &_ZM_mul));
851 : }
852 : /********************************************************************/
853 : /** **/
854 : /** ADD, SUB **/
855 : /** **/
856 : /********************************************************************/
857 : static GEN
858 36064610 : ZC_add_i(GEN x, GEN y, long lx)
859 : {
860 36064610 : GEN A = cgetg(lx, t_COL);
861 : long i;
862 484298850 : for (i=1; i<lx; i++) gel(A,i) = addii(gel(x,i), gel(y,i));
863 36056293 : return A;
864 : }
865 : GEN
866 26037153 : ZC_add(GEN x, GEN y) { return ZC_add_i(x, y, lg(x)); }
867 : GEN
868 365423 : ZC_Z_add(GEN x, GEN y)
869 : {
870 365423 : long k, lx = lg(x);
871 365423 : GEN z = cgetg(lx, t_COL);
872 365423 : if (lx == 1) pari_err_TYPE2("+",x,y);
873 365423 : gel(z,1) = addii(y,gel(x,1));
874 2513426 : for (k = 2; k < lx; k++) gel(z,k) = icopy(gel(x,k));
875 365438 : return z;
876 : }
877 :
878 : static GEN
879 9085778 : ZC_sub_i(GEN x, GEN y, long lx)
880 : {
881 : long i;
882 9085778 : GEN A = cgetg(lx, t_COL);
883 64507266 : for (i=1; i<lx; i++) gel(A,i) = subii(gel(x,i), gel(y,i));
884 9084997 : return A;
885 : }
886 : GEN
887 8757262 : ZC_sub(GEN x, GEN y) { return ZC_sub_i(x, y, lg(x)); }
888 : GEN
889 0 : ZC_Z_sub(GEN x, GEN y)
890 : {
891 0 : long k, lx = lg(x);
892 0 : GEN z = cgetg(lx, t_COL);
893 0 : if (lx == 1) pari_err_TYPE2("+",x,y);
894 0 : gel(z,1) = subii(gel(x,1), y);
895 0 : for (k = 2; k < lx; k++) gel(z,k) = icopy(gel(x,k));
896 0 : return z;
897 : }
898 : GEN
899 545305 : Z_ZC_sub(GEN a, GEN x)
900 : {
901 545305 : long k, lx = lg(x);
902 545305 : GEN z = cgetg(lx, t_COL);
903 545303 : if (lx == 1) pari_err_TYPE2("-",a,x);
904 545303 : gel(z,1) = subii(a, gel(x,1));
905 1514845 : for (k = 2; k < lx; k++) gel(z,k) = negi(gel(x,k));
906 545305 : return z;
907 : }
908 :
909 : GEN
910 843752 : ZM_add(GEN x, GEN y)
911 : {
912 843752 : long lx = lg(x), l, j;
913 : GEN z;
914 843752 : if (lx == 1) return cgetg(1, t_MAT);
915 803964 : z = cgetg(lx, t_MAT); l = lgcols(x);
916 10831433 : for (j = 1; j < lx; j++) gel(z,j) = ZC_add_i(gel(x,j), gel(y,j), l);
917 803964 : return z;
918 : }
919 : GEN
920 155533 : ZM_sub(GEN x, GEN y)
921 : {
922 155533 : long lx = lg(x), l, j;
923 : GEN z;
924 155533 : if (lx == 1) return cgetg(1, t_MAT);
925 115745 : z = cgetg(lx, t_MAT); l = lgcols(x);
926 444228 : for (j = 1; j < lx; j++) gel(z,j) = ZC_sub_i(gel(x,j), gel(y,j), l);
927 115745 : return z;
928 : }
929 : /********************************************************************/
930 : /** **/
931 : /** LINEAR COMBINATION **/
932 : /** **/
933 : /********************************************************************/
934 : /* return X/c assuming division is exact */
935 : GEN
936 5497788 : ZC_Z_divexact(GEN x, GEN c)
937 78263031 : { pari_APPLY_type(t_COL, diviiexact(gel(x,i), c)) }
938 : GEN
939 2261 : ZC_divexactu(GEN x, ulong c)
940 11375 : { pari_APPLY_type(t_COL, diviuexact(gel(x,i), c)) }
941 :
942 : GEN
943 490825 : ZM_Z_divexact(GEN x, GEN c)
944 3917216 : { pari_APPLY_same(ZC_Z_divexact(gel(x,i), c)) }
945 :
946 : GEN
947 441 : ZM_divexactu(GEN x, ulong c)
948 2702 : { pari_APPLY_same(ZC_divexactu(gel(x,i), c)) }
949 :
950 : GEN
951 34675546 : ZC_Z_mul(GEN x, GEN c)
952 : {
953 34675546 : if (!signe(c)) return zerocol(lg(x)-1);
954 33338312 : if (is_pm1(c)) return (signe(c) > 0)? ZC_copy(x): ZC_neg(x);
955 253731514 : pari_APPLY_type(t_COL, mulii(gel(x,i), c))
956 : }
957 :
958 : GEN
959 62247 : ZC_z_mul(GEN x, long c)
960 : {
961 62247 : if (!c) return zerocol(lg(x)-1);
962 52368 : if (c == 1) return ZC_copy(x);
963 47453 : if (c ==-1) return ZC_neg(x);
964 477380 : pari_APPLY_type(t_COL, mulsi(c, gel(x,i)))
965 : }
966 :
967 : GEN
968 28500 : zv_z_mul(GEN x, long n)
969 432141 : { pari_APPLY_long(x[i]*n) }
970 :
971 : /* return a ZM */
972 : GEN
973 0 : nm_Z_mul(GEN X, GEN c)
974 : {
975 0 : long i, j, h, l = lg(X), s = signe(c);
976 : GEN A;
977 0 : if (l == 1) return cgetg(1, t_MAT);
978 0 : h = lgcols(X);
979 0 : if (!s) return zeromat(h-1, l-1);
980 0 : if (is_pm1(c)) {
981 0 : if (s > 0) return Flm_to_ZM(X);
982 0 : X = Flm_to_ZM(X); ZM_togglesign(X); return X;
983 : }
984 0 : A = cgetg(l, t_MAT);
985 0 : for (j = 1; j < l; j++)
986 : {
987 0 : GEN a = cgetg(h, t_COL), x = gel(X, j);
988 0 : for (i = 1; i < h; i++) gel(a,i) = muliu(c, x[i]);
989 0 : gel(A,j) = a;
990 : }
991 0 : return A;
992 : }
993 : GEN
994 3268268 : ZM_Z_mul(GEN X, GEN c)
995 : {
996 3268268 : long i, j, h, l = lg(X);
997 : GEN A;
998 3268268 : if (l == 1) return cgetg(1, t_MAT);
999 3268268 : h = lgcols(X);
1000 3268264 : if (!signe(c)) return zeromat(h-1, l-1);
1001 3222919 : if (is_pm1(c)) return (signe(c) > 0)? ZM_copy(X): ZM_neg(X);
1002 2822922 : A = cgetg(l, t_MAT);
1003 15760038 : for (j = 1; j < l; j++)
1004 : {
1005 12937441 : GEN a = cgetg(h, t_COL), x = gel(X, j);
1006 219994779 : for (i = 1; i < h; i++) gel(a,i) = mulii(c, gel(x,i));
1007 12937129 : gel(A,j) = a;
1008 : }
1009 2822597 : return A;
1010 : }
1011 : void
1012 73140811 : ZC_lincomb1_inplace_i(GEN X, GEN Y, GEN v, long n)
1013 : {
1014 : long i;
1015 1619063609 : for (i = n; i; i--) gel(X,i) = addmulii_inplace(gel(X,i), gel(Y,i), v);
1016 73089247 : }
1017 : /* X <- X + v Y (elementary col operation) */
1018 : void
1019 62189856 : ZC_lincomb1_inplace(GEN X, GEN Y, GEN v)
1020 : {
1021 62189856 : if (lgefint(v) != 2) return ZC_lincomb1_inplace_i(X, Y, v, lg(X)-1);
1022 : }
1023 : void
1024 29656476 : Flc_lincomb1_inplace(GEN X, GEN Y, ulong v, ulong q)
1025 : {
1026 : long i;
1027 29656476 : if (!v) return; /* v = 0 */
1028 743720409 : for (i = lg(X)-1; i; i--) X[i] = Fl_add(X[i], Fl_mul(Y[i], v, q), q);
1029 : }
1030 :
1031 : /* X + v Y, wasteful if (v = 0) */
1032 : static GEN
1033 15389890 : ZC_lincomb1(GEN v, GEN x, GEN y)
1034 121950220 : { pari_APPLY_type(t_COL, addmulii(gel(x,i), gel(y,i), v)) }
1035 :
1036 : /* -X + vY */
1037 : static GEN
1038 741083 : ZC_lincomb_1(GEN v, GEN x, GEN y)
1039 4588821 : { pari_APPLY_type(t_COL, mulsubii(gel(y,i), v, gel(x,i))) }
1040 :
1041 : /* X,Y compatible ZV; u,v in Z. Returns A = u*X + v*Y */
1042 : GEN
1043 33292143 : ZC_lincomb(GEN u, GEN v, GEN X, GEN Y)
1044 : {
1045 : long su, sv;
1046 : GEN A;
1047 :
1048 33292143 : su = signe(u); if (!su) return ZC_Z_mul(Y, v);
1049 33291191 : sv = signe(v); if (!sv) return ZC_Z_mul(X, u);
1050 32861371 : if (is_pm1(v))
1051 : {
1052 11087963 : if (is_pm1(u))
1053 : {
1054 9761250 : if (su != sv) A = ZC_sub(X, Y);
1055 2710845 : else A = ZC_add(X, Y);
1056 9760584 : if (su < 0) ZV_togglesign(A); /* in place but was created above */
1057 : }
1058 : else
1059 : {
1060 1326703 : if (sv > 0) A = ZC_lincomb1 (u, Y, X);
1061 606580 : else A = ZC_lincomb_1(u, Y, X);
1062 : }
1063 : }
1064 21774030 : else if (is_pm1(u))
1065 : {
1066 14803791 : if (su > 0) A = ZC_lincomb1 (v, X, Y);
1067 133911 : else A = ZC_lincomb_1(v, X, Y);
1068 : }
1069 : else
1070 : { /* not cgetg_copy: x may be a t_VEC */
1071 6974240 : long i, lx = lg(X);
1072 6974240 : A = cgetg(lx,t_COL);
1073 44823038 : for (i=1; i<lx; i++) gel(A,i) = lincombii(u,v,gel(X,i),gel(Y,i));
1074 : }
1075 32851940 : return A;
1076 : }
1077 :
1078 : /********************************************************************/
1079 : /** **/
1080 : /** CONVERSIONS **/
1081 : /** **/
1082 : /********************************************************************/
1083 : GEN
1084 466573 : ZV_to_nv(GEN x)
1085 883931 : { pari_APPLY_ulong(itou(gel(x,i))) }
1086 :
1087 : GEN
1088 156610 : zm_to_ZM(GEN x)
1089 912516 : { pari_APPLY_type(t_MAT, zc_to_ZC(gel(x,i))) }
1090 :
1091 : GEN
1092 126 : zmV_to_ZMV(GEN x)
1093 791 : { pari_APPLY_type(t_VEC, zm_to_ZM(gel(x,i))) }
1094 :
1095 : /* same as Flm_to_ZM but do not assume positivity */
1096 : GEN
1097 1022 : ZM_to_zm(GEN x)
1098 17472 : { pari_APPLY_same(ZV_to_zv(gel(x,i))) }
1099 :
1100 : GEN
1101 366646 : zv_to_Flv(GEN x, ulong p)
1102 5418812 : { pari_APPLY_ulong(umodsu(x[i], p)) }
1103 :
1104 : GEN
1105 22694 : zm_to_Flm(GEN x, ulong p)
1106 351750 : { pari_APPLY_same(zv_to_Flv(gel(x,i),p)) }
1107 :
1108 : GEN
1109 49 : ZMV_to_zmV(GEN x)
1110 399 : { pari_APPLY_type(t_VEC, ZM_to_zm(gel(x,i))) }
1111 :
1112 : /********************************************************************/
1113 : /** **/
1114 : /** COPY, NEGATION **/
1115 : /** **/
1116 : /********************************************************************/
1117 : GEN
1118 17406086 : ZC_copy(GEN x)
1119 : {
1120 17406086 : long i, lx = lg(x);
1121 17406086 : GEN y = cgetg(lx, t_COL);
1122 135134978 : for (i=1; i<lx; i++)
1123 : {
1124 117728595 : GEN c = gel(x,i);
1125 117728595 : gel(y,i) = lgefint(c) == 2? gen_0: icopy(c);
1126 : }
1127 17406383 : return y;
1128 : }
1129 :
1130 : GEN
1131 654098 : ZM_copy(GEN x)
1132 5982648 : { pari_APPLY_same(ZC_copy(gel(x,i))) }
1133 :
1134 : void
1135 345679 : ZV_neg_inplace(GEN M)
1136 : {
1137 345679 : long l = lg(M);
1138 1303173 : while (--l > 0) gel(M,l) = negi(gel(M,l));
1139 345744 : }
1140 : GEN
1141 6345770 : ZC_neg(GEN x)
1142 37180738 : { pari_APPLY_type(t_COL, negi(gel(x,i))) }
1143 :
1144 : GEN
1145 51884 : zv_neg(GEN x)
1146 662823 : { pari_APPLY_long(-x[i]) }
1147 : GEN
1148 126 : zv_neg_inplace(GEN M)
1149 : {
1150 126 : long l = lg(M);
1151 429 : while (--l > 0) M[l] = -M[l];
1152 126 : return M;
1153 : }
1154 : GEN
1155 77 : zv_abs(GEN x)
1156 5446 : { pari_APPLY_ulong(labs(x[i])) }
1157 : GEN
1158 1722584 : ZM_neg(GEN x)
1159 5157434 : { pari_APPLY_same(ZC_neg(gel(x,i))) }
1160 :
1161 : void
1162 5038460 : ZV_togglesign(GEN M)
1163 : {
1164 5038460 : long l = lg(M);
1165 75531461 : while (--l > 0) togglesign_safe(&gel(M,l));
1166 5038530 : }
1167 : void
1168 0 : ZM_togglesign(GEN M)
1169 : {
1170 0 : long l = lg(M);
1171 0 : while (--l > 0) ZV_togglesign(gel(M,l));
1172 0 : }
1173 :
1174 : /********************************************************************/
1175 : /** **/
1176 : /** "DIVISION" mod HNF **/
1177 : /** **/
1178 : /********************************************************************/
1179 : /* Reduce ZC x modulo ZM y in HNF */
1180 : static GEN
1181 10751199 : ZC_hnfdivrem_i(GEN x, GEN y, GEN *Q, GEN (*div)(GEN,GEN))
1182 : {
1183 10751199 : long i, l = lg(x);
1184 10751199 : pari_sp av = avma;
1185 :
1186 10751199 : if (Q) *Q = cgetg(l,t_COL);
1187 10751199 : if (l == 1) return cgetg(1,t_COL);
1188 61758286 : for (i = l-1; i>0; i--)
1189 : {
1190 51010241 : GEN q = div(gel(x,i), gcoeff(y,i,i));
1191 51009403 : if (signe(q)) x = ZC_lincomb(gen_1, negi(q), x, gel(y,i));
1192 51007087 : if (Q) gel(*Q, i) = q;
1193 : }
1194 10748045 : if (avma == av) return ZC_copy(x);
1195 8004020 : if (!Q) return gerepileupto(av, x);
1196 57964 : gerepileall(av, 2, &x, Q); return x;
1197 : }
1198 : GEN
1199 10193063 : ZC_hnfdivrem(GEN x, GEN y, GEN *Q)
1200 10193063 : { return ZC_hnfdivrem_i(x, y, Q, diviiround); }
1201 : GEN
1202 28 : ZC_modhnf(GEN x, GEN y, GEN *Q)
1203 28 : { return ZC_hnfdivrem_i(x, y, Q, truedivii); }
1204 :
1205 : /* Return R such that x = y Q + R, y integral HNF */
1206 : static GEN
1207 454460 : ZM_hnfdivrem_i(GEN x, GEN y, GEN *Q, GEN (*div)(GEN,GEN))
1208 : {
1209 : long l, i;
1210 454460 : GEN R = cgetg_copy(x, &l);
1211 454507 : if (Q)
1212 : {
1213 127640 : GEN q = cgetg(l, t_MAT); *Q = q;
1214 188596 : for (i = 1; i < l; i++)
1215 60951 : gel(R,i) = ZC_hnfdivrem_i(gel(x,i),y,&gel(q,i),div);
1216 : }
1217 : else
1218 824195 : for (i = 1; i < l; i++)
1219 497334 : gel(R,i) = ZC_hnfdivrem_i(gel(x,i),y,NULL,div);
1220 454506 : return R;
1221 : }
1222 : GEN
1223 454443 : ZM_hnfdivrem(GEN x, GEN y, GEN *Q)
1224 454443 : { return ZM_hnfdivrem_i(x, y, Q, diviiround); }
1225 : GEN
1226 14 : ZM_modhnf(GEN x, GEN y, GEN *Q)
1227 14 : { return ZM_hnfdivrem_i(x, y, Q, truedivii); }
1228 :
1229 : static GEN
1230 7 : ZV_ZV_divrem(GEN x, GEN y, GEN *pQ)
1231 : {
1232 7 : long i, l = lg(x), tx = typ(x);
1233 : GEN Q, R;
1234 :
1235 7 : if (!pQ) return ZV_ZV_mod(x, y);
1236 0 : Q = cgetg(l,tx);
1237 0 : R = cgetg(l,tx);
1238 0 : for (i = 1; i < l; i++) gel(Q,i) = truedvmdii(gel(x,i), gel(y,i), &gel(R,i));
1239 0 : *pQ = Q; return R;
1240 : }
1241 : static GEN
1242 0 : ZM_ZV_divrem(GEN x, GEN y, GEN *Q)
1243 0 : { if (!Q) return ZM_ZV_mod(x, y);
1244 0 : pari_APPLY_same(ZV_ZV_divrem(gel(x,i), y, Q)); }
1245 :
1246 : static int
1247 42 : RgM_issquare(GEN x) { long l = lg(x); return l == 1 || lg(gel(x,1)) == l; }
1248 : static void
1249 98 : matmodhnf_check(GEN x)
1250 : {
1251 98 : switch(typ(x))
1252 : {
1253 42 : case t_VEC: case t_COL:
1254 42 : if (!RgV_is_ZV(x)) pari_err_TYPE("matmodhnf", x);
1255 42 : break;
1256 56 : case t_MAT:
1257 56 : if (!RgM_is_ZM(x)) pari_err_TYPE("matmodhnf", x);
1258 56 : break;
1259 0 : default: pari_err_TYPE("matmodhnf", x);
1260 : }
1261 98 : }
1262 : GEN
1263 49 : matmodhnf(GEN x, GEN y, GEN *Q)
1264 : {
1265 49 : long tx = typ(x), ty = typ(y), ly, lx;
1266 49 : matmodhnf_check(x); lx = lg(x);
1267 49 : matmodhnf_check(y); ly = lg(y);
1268 49 : if (ty == t_MAT && !RgM_issquare(y)) pari_err_TYPE("matmodhnf", y);
1269 49 : if (tx == t_MAT && lx == 1)
1270 : {
1271 0 : if (ly != 1) pari_err_DIM("matmodhnf");
1272 0 : if (!Q) *Q = cgetg(1, t_MAT);
1273 0 : return cgetg(1, t_MAT);
1274 : }
1275 49 : if (is_vec_t(ty))
1276 7 : return tx == t_MAT? ZM_ZV_divrem(x, y, Q): ZV_ZV_divrem(x, y, Q);
1277 : /* ty = t_MAT */
1278 42 : if (tx == t_MAT) return ZM_modhnf(x, y, Q);
1279 28 : x = ZC_modhnf(x, y, Q);
1280 28 : if (tx == t_VEC) { settyp(x, tx); if (Q) settyp(*Q, tx); }
1281 28 : return x;
1282 : }
1283 :
1284 : /********************************************************************/
1285 : /** **/
1286 : /** TESTS **/
1287 : /** **/
1288 : /********************************************************************/
1289 : int
1290 23103520 : zv_equal0(GEN V)
1291 : {
1292 23103520 : long l = lg(V);
1293 37579987 : while (--l > 0)
1294 30805728 : if (V[l]) return 0;
1295 6774259 : return 1;
1296 : }
1297 :
1298 : int
1299 14397145 : ZV_equal0(GEN V)
1300 : {
1301 14397145 : long l = lg(V);
1302 25490196 : while (--l > 0)
1303 25062878 : if (signe(gel(V,l))) return 0;
1304 427318 : return 1;
1305 : }
1306 : int
1307 16231 : ZMrow_equal0(GEN V, long i)
1308 : {
1309 16231 : long l = lg(V);
1310 25183 : while (--l > 0)
1311 21679 : if (signe(gcoeff(V,i,l))) return 0;
1312 3504 : return 1;
1313 : }
1314 :
1315 : static int
1316 6227028 : ZV_equal_lg(GEN V, GEN W, long l)
1317 : {
1318 25809428 : while (--l > 0)
1319 20071281 : if (!equalii(gel(V,l), gel(W,l))) return 0;
1320 5738147 : return 1;
1321 : }
1322 : int
1323 293632 : ZV_equal(GEN V, GEN W)
1324 : {
1325 293632 : long l = lg(V);
1326 293632 : if (lg(W) != l) return 0;
1327 293611 : return ZV_equal_lg(V, W, l);
1328 : }
1329 : int
1330 3327314 : ZM_equal(GEN A, GEN B)
1331 : {
1332 3327314 : long i, m, l = lg(A);
1333 3327314 : if (lg(B) != l) return 0;
1334 3327260 : if (l == 1) return 1;
1335 3327260 : m = lgcols(A);
1336 3327254 : if (lgcols(B) != m) return 0;
1337 8990080 : for (i = 1; i < l; i++)
1338 5933399 : if (!ZV_equal_lg(gel(A,i), gel(B,i), m)) return 0;
1339 3056681 : return 1;
1340 : }
1341 : int
1342 73962 : ZM_equal0(GEN A)
1343 : {
1344 73962 : long i, j, m, l = lg(A);
1345 73962 : if (l == 1) return 1;
1346 73962 : m = lgcols(A);
1347 129561 : for (j = 1; j < l; j++)
1348 2763485 : for (i = 1; i < m; i++)
1349 2707886 : if (signe(gcoeff(A,i,j))) return 0;
1350 17146 : return 1;
1351 : }
1352 : int
1353 30987605 : zv_equal(GEN V, GEN W)
1354 : {
1355 30987605 : long l = lg(V);
1356 30987605 : if (lg(W) != l) return 0;
1357 263091908 : while (--l > 0)
1358 233220556 : if (V[l] != W[l]) return 0;
1359 29871352 : return 1;
1360 : }
1361 :
1362 : int
1363 1638 : zvV_equal(GEN V, GEN W)
1364 : {
1365 1638 : long l = lg(V);
1366 1638 : if (lg(W) != l) return 0;
1367 80388 : while (--l > 0)
1368 79912 : if (!zv_equal(gel(V,l),gel(W,l))) return 0;
1369 476 : return 1;
1370 : }
1371 :
1372 : int
1373 752217 : ZM_ishnf(GEN x)
1374 : {
1375 752217 : long i,j, lx = lg(x);
1376 2260016 : for (i=1; i<lx; i++)
1377 : {
1378 1620614 : GEN xii = gcoeff(x,i,i);
1379 1620614 : if (signe(xii) <= 0) return 0;
1380 3152914 : for (j=1; j<i; j++)
1381 1569768 : if (signe(gcoeff(x,i,j))) return 0;
1382 3211436 : for (j=i+1; j<lx; j++)
1383 : {
1384 1703637 : GEN xij = gcoeff(x,i,j);
1385 1703637 : if (signe(xij)<0 || cmpii(xij,xii)>=0) return 0;
1386 : }
1387 : }
1388 639402 : return 1;
1389 : }
1390 : int
1391 658005 : ZM_isidentity(GEN x)
1392 : {
1393 658005 : long i,j, lx = lg(x);
1394 :
1395 658005 : if (lx == 1) return 1;
1396 657998 : if (lx != lgcols(x)) return 0;
1397 3215429 : for (j=1; j<lx; j++)
1398 : {
1399 2557488 : GEN c = gel(x,j);
1400 8023691 : for (i=1; i<j; )
1401 5466225 : if (signe(gel(c,i++))) return 0;
1402 : /* i = j */
1403 2557466 : if (!equali1(gel(c,i++))) return 0;
1404 8023706 : for ( ; i<lx; )
1405 5466261 : if (signe(gel(c,i++))) return 0;
1406 : }
1407 657941 : return 1;
1408 : }
1409 : int
1410 556686 : ZM_isdiagonal(GEN x)
1411 : {
1412 556686 : long i,j, lx = lg(x);
1413 556686 : if (lx == 1) return 1;
1414 556686 : if (lx != lgcols(x)) return 0;
1415 :
1416 1439173 : for (j=1; j<lx; j++)
1417 : {
1418 1209362 : GEN c = gel(x,j);
1419 1698099 : for (i=1; i<j; i++)
1420 815566 : if (signe(gel(c,i))) return 0;
1421 2022093 : for (i++; i<lx; i++)
1422 1139595 : if (signe(gel(c,i))) return 0;
1423 : }
1424 229811 : return 1;
1425 : }
1426 : int
1427 159995 : ZM_isscalar(GEN x, GEN s)
1428 : {
1429 159995 : long i, j, lx = lg(x);
1430 :
1431 159995 : if (lx == 1) return 1;
1432 159995 : if (!s) s = gcoeff(x,1,1);
1433 159995 : if (equali1(s)) return ZM_isidentity(x);
1434 158538 : if (lx != lgcols(x)) return 0;
1435 712929 : for (j=1; j<lx; j++)
1436 : {
1437 615776 : GEN c = gel(x,j);
1438 2832284 : for (i=1; i<j; )
1439 2275794 : if (signe(gel(c,i++))) return 0;
1440 : /* i = j */
1441 556490 : if (!equalii(gel(c,i++), s)) return 0;
1442 2860386 : for ( ; i<lx; )
1443 2305995 : if (signe(gel(c,i++))) return 0;
1444 : }
1445 97153 : return 1;
1446 : }
1447 :
1448 : long
1449 174398 : ZC_is_ei(GEN x)
1450 : {
1451 174398 : long i, j = 0, l = lg(x);
1452 1789913 : for (i = 1; i < l; i++)
1453 : {
1454 1615516 : GEN c = gel(x,i);
1455 1615516 : long s = signe(c);
1456 1615516 : if (!s) continue;
1457 174392 : if (s < 0 || !is_pm1(c) || j) return 0;
1458 174391 : j = i;
1459 : }
1460 174397 : return j;
1461 : }
1462 :
1463 : /********************************************************************/
1464 : /** **/
1465 : /** MISCELLANEOUS **/
1466 : /** **/
1467 : /********************************************************************/
1468 : /* assume lg(x) = lg(y), x,y in Z^n */
1469 : int
1470 3147281 : ZV_cmp(GEN x, GEN y)
1471 : {
1472 3147281 : long fl,i, lx = lg(x);
1473 6376287 : for (i=1; i<lx; i++)
1474 5080679 : if (( fl = cmpii(gel(x,i), gel(y,i)) )) return fl;
1475 1295608 : return 0;
1476 : }
1477 : /* assume lg(x) = lg(y), x,y in Z^n */
1478 : int
1479 19740 : ZV_abscmp(GEN x, GEN y)
1480 : {
1481 19740 : long fl,i, lx = lg(x);
1482 54201 : for (i=1; i<lx; i++)
1483 54074 : if (( fl = abscmpii(gel(x,i), gel(y,i)) )) return fl;
1484 127 : return 0;
1485 : }
1486 :
1487 : long
1488 18452708 : zv_content(GEN x)
1489 : {
1490 18452708 : long i, s, l = lg(x);
1491 18452708 : if (l == 1) return 0;
1492 18452701 : s = labs(x[1]);
1493 42002344 : for (i = 2; i < l && s != 1; i++) s = ugcd(s, labs(x[i]));
1494 18452767 : return s;
1495 : }
1496 : GEN
1497 297390 : ZV_content(GEN x)
1498 : {
1499 297390 : long i, l = lg(x);
1500 297390 : pari_sp av = avma;
1501 : GEN c;
1502 297390 : if (l == 1) return gen_0;
1503 297390 : if (l == 2) return absi(gel(x,1));
1504 204864 : c = gel(x,1);
1505 558021 : for (i = 2; i < l; i++)
1506 : {
1507 403684 : c = gcdii(c, gel(x,i));
1508 403684 : if (is_pm1(c)) { set_avma(av); return gen_1; }
1509 : }
1510 154337 : return gerepileuptoint(av, c);
1511 : }
1512 :
1513 : GEN
1514 3953088 : ZM_det_triangular(GEN mat)
1515 : {
1516 : pari_sp av;
1517 3953088 : long i,l = lg(mat);
1518 : GEN s;
1519 :
1520 3953088 : if (l<3) return l<2? gen_1: icopy(gcoeff(mat,1,1));
1521 3522308 : av = avma; s = gcoeff(mat,1,1);
1522 9539977 : for (i=2; i<l; i++) s = mulii(s,gcoeff(mat,i,i));
1523 3521874 : return gerepileuptoint(av,s);
1524 : }
1525 :
1526 : /* assumes no overflow */
1527 : long
1528 950721 : zv_prod(GEN v)
1529 : {
1530 950721 : long n, i, l = lg(v);
1531 950721 : if (l == 1) return 1;
1532 960749 : n = v[1]; for (i = 2; i < l; i++) n *= v[i];
1533 772242 : return n;
1534 : }
1535 :
1536 : static GEN
1537 319043558 : _mulii(void *E, GEN a, GEN b)
1538 319043558 : { (void) E; return mulii(a, b); }
1539 :
1540 : /* product of ulongs */
1541 : GEN
1542 1864644 : zv_prod_Z(GEN v)
1543 : {
1544 : pari_sp av;
1545 1864644 : long k, m, n = lg(v)-1;
1546 1864644 : int stop = 0;
1547 : GEN V;
1548 1864644 : switch(n) {
1549 21035 : case 0: return gen_1;
1550 130095 : case 1: return utoi(v[1]);
1551 976951 : case 2: return muluu(v[1], v[2]);
1552 : }
1553 736563 : av = avma; m = n >> 1;
1554 736563 : V = cgetg(m + (odd(n)? 2: 1), t_VEC);
1555 154564630 : for (k = n; k; k--) /* start from the end: v is usually sorted */
1556 153829133 : if (v[k] & HIGHMASK) { stop = 1; break; }
1557 2420919 : while (!stop)
1558 : { /* HACK: handle V as a t_VECSMALL; gain a few iterations */
1559 83105044 : for (k = 1; k <= m; k++)
1560 : {
1561 80812326 : V[k] = uel(v,k<<1) * uel(v,(k<<1)-1);
1562 80812326 : if (V[k] & HIGHMASK) stop = 1; /* last "free" iteration */
1563 : }
1564 2292718 : if (odd(n))
1565 : {
1566 1351357 : if (n == 1) { set_avma(av); return utoi(v[1]); }
1567 742995 : V[++m] = v[n];
1568 : }
1569 1684356 : v = V; n = m; m = n >> 1;
1570 : }
1571 : /* n > 1; m > 0 */
1572 128201 : if (n == 2) { set_avma(av); return muluu(v[1], v[2]); }
1573 47265809 : for (k = 1; k <= m; k++) gel(V,k) = muluu(v[k<<1], v[(k<<1)-1]);
1574 88065 : if (odd(n)) gel(V, ++m) = utoipos(v[n]);
1575 88065 : setlg(V, m+1); /* HACK: now V is a bona fide t_VEC */
1576 88065 : return gerepileuptoint(av, gen_product(V, NULL, &_mulii));
1577 : }
1578 : GEN
1579 14694393 : vecsmall_prod(GEN v)
1580 : {
1581 14694393 : pari_sp av = avma;
1582 14694393 : long k, m, n = lg(v)-1;
1583 : GEN V;
1584 14694393 : switch (n) {
1585 0 : case 0: return gen_1;
1586 0 : case 1: return stoi(v[1]);
1587 21 : case 2: return mulss(v[1], v[2]);
1588 : }
1589 14694372 : m = n >> 1;
1590 14694372 : V = cgetg(m + (odd(n)? 2: 1), t_VEC);
1591 161556906 : for (k = 1; k <= m; k++) gel(V,k) = mulss(v[k<<1], v[(k<<1)-1]);
1592 14694372 : if (odd(n)) gel(V,k) = stoi(v[n]);
1593 14694372 : return gerepileuptoint(av, gen_product(V, NULL, &_mulii));
1594 : }
1595 :
1596 : GEN
1597 8829748 : ZV_prod(GEN v)
1598 : {
1599 8829748 : pari_sp av = avma;
1600 8829748 : long i, l = lg(v);
1601 : GEN n;
1602 8829748 : if (l == 1) return gen_1;
1603 8644792 : if (l > 7) return gerepileuptoint(av, gen_product(v, NULL, _mulii));
1604 1301798 : n = gel(v,1);
1605 1301798 : if (l == 2) return icopy(n);
1606 2104061 : for (i = 2; i < l; i++) n = mulii(n, gel(v,i));
1607 846316 : return gerepileuptoint(av, n);
1608 : }
1609 : /* assumes no overflow */
1610 : long
1611 16540 : zv_sum(GEN v)
1612 : {
1613 16540 : long n, i, l = lg(v);
1614 16540 : if (l == 1) return 0;
1615 95939 : n = v[1]; for (i = 2; i < l; i++) n += v[i];
1616 16519 : return n;
1617 : }
1618 : /* assumes no overflow and 0 <= n <= #v */
1619 : long
1620 0 : zv_sumpart(GEN v, long n)
1621 : {
1622 : long i, p;
1623 0 : if (!n) return 0;
1624 0 : p = v[1]; for (i = 2; i <= n; i++) p += v[i];
1625 0 : return p;
1626 : }
1627 : GEN
1628 77 : ZV_sum(GEN v)
1629 : {
1630 77 : pari_sp av = avma;
1631 77 : long i, l = lg(v);
1632 : GEN n;
1633 77 : if (l == 1) return gen_0;
1634 77 : n = gel(v,1);
1635 77 : if (l == 2) return icopy(n);
1636 581 : for (i = 2; i < l; i++) n = addii(n, gel(v,i));
1637 77 : return gerepileuptoint(av, n);
1638 : }
1639 :
1640 : /********************************************************************/
1641 : /** **/
1642 : /** GRAM SCHMIDT REDUCTION (integer matrices) **/
1643 : /** **/
1644 : /********************************************************************/
1645 :
1646 : /* L[k,] += q * L[l,], l < k. Inefficient if q = 0 */
1647 : static void
1648 360340 : Zupdate_row(long k, long l, GEN q, GEN L, GEN B)
1649 : {
1650 360340 : long i, qq = itos_or_0(q);
1651 360344 : if (!qq)
1652 : {
1653 36730 : for(i=1;i<l;i++) gcoeff(L,k,i) = addii(gcoeff(L,k,i),mulii(q,gcoeff(L,l,i)));
1654 7587 : gcoeff(L,k,l) = addii(gcoeff(L,k,l), mulii(q,B));
1655 7587 : return;
1656 : }
1657 352757 : if (qq == 1) {
1658 148961 : for (i=1;i<l; i++) gcoeff(L,k,i) = addii(gcoeff(L,k,i),gcoeff(L,l,i));
1659 102354 : gcoeff(L,k,l) = addii(gcoeff(L,k,l), B);
1660 250401 : } else if (qq == -1) {
1661 150651 : for (i=1;i<l; i++) gcoeff(L,k,i) = subii(gcoeff(L,k,i),gcoeff(L,l,i));
1662 88194 : gcoeff(L,k,l) = addii(gcoeff(L,k,l), negi(B));
1663 : } else {
1664 287197 : for(i=1;i<l;i++) gcoeff(L,k,i) = addii(gcoeff(L,k,i),mulsi(qq,gcoeff(L,l,i)));
1665 162226 : gcoeff(L,k,l) = addii(gcoeff(L,k,l), mulsi(qq,B));
1666 : }
1667 : }
1668 :
1669 : /* update L[k,] */
1670 : static void
1671 1037680 : ZRED(long k, long l, GEN x, GEN L, GEN B)
1672 : {
1673 1037680 : GEN q = truedivii(addii(B,shifti(gcoeff(L,k,l),1)), shifti(B,1));
1674 1037786 : if (!signe(q)) return;
1675 360316 : q = negi(q);
1676 360333 : Zupdate_row(k,l,q,L,B);
1677 360322 : gel(x,k) = ZC_lincomb(gen_1, q, gel(x,k), gel(x,l));
1678 : }
1679 :
1680 : /* Gram-Schmidt reduction, x a ZM */
1681 : static void
1682 1189977 : ZincrementalGS(GEN x, GEN L, GEN B, long k)
1683 : {
1684 : long i, j;
1685 3802070 : for (j=1; j<=k; j++)
1686 : {
1687 2612987 : pari_sp av = avma;
1688 2612987 : GEN u = ZV_dotproduct(gel(x,k), gel(x,j));
1689 5653681 : for (i=1; i<j; i++)
1690 : {
1691 3042446 : u = subii(mulii(gel(B,i+1), u), mulii(gcoeff(L,k,i), gcoeff(L,j,i)));
1692 3041826 : u = diviiexact(u, gel(B,i));
1693 : }
1694 2611235 : gcoeff(L,k,j) = gerepileuptoint(av, u);
1695 : }
1696 1189083 : gel(B,k+1) = gcoeff(L,k,k); gcoeff(L,k,k) = gen_1;
1697 1189083 : }
1698 :
1699 : /* Variant reducemodinvertible(ZC v, ZM y), when y singular.
1700 : * Very inefficient if y is not LLL-reduced of maximal rank */
1701 : static GEN
1702 112233 : ZC_reducemodmatrix_i(GEN v, GEN y)
1703 : {
1704 112233 : GEN B, L, x = shallowconcat(y, v);
1705 112234 : long k, lx = lg(x), nx = lx-1;
1706 :
1707 112234 : B = scalarcol_shallow(gen_1, lx);
1708 112234 : L = zeromatcopy(nx, nx);
1709 461660 : for (k=1; k <= nx; k++) ZincrementalGS(x, L, B, k);
1710 349428 : for (k = nx-1; k >= 1; k--) ZRED(nx,k, x,L,gel(B,k+1));
1711 112221 : return gel(x,nx);
1712 : }
1713 : GEN
1714 112233 : ZC_reducemodmatrix(GEN v, GEN y) {
1715 112233 : pari_sp av = avma;
1716 112233 : return gerepilecopy(av, ZC_reducemodmatrix_i(v,y));
1717 : }
1718 :
1719 : /* Variant reducemodinvertible(ZM v, ZM y), when y singular.
1720 : * Very inefficient if y is not LLL-reduced of maximal rank */
1721 : static GEN
1722 227040 : ZM_reducemodmatrix_i(GEN v, GEN y)
1723 : {
1724 : GEN B, L, V;
1725 227040 : long j, k, lv = lg(v), nx = lg(y), lx = nx+1;
1726 :
1727 227040 : V = cgetg(lv, t_MAT);
1728 227058 : B = scalarcol_shallow(gen_1, lx);
1729 227065 : L = zeromatcopy(nx, nx);
1730 601587 : for (k=1; k < nx; k++) ZincrementalGS(y, L, B, k);
1731 693076 : for (j = 1; j < lg(v); j++)
1732 : {
1733 466056 : GEN x = shallowconcat(y, gel(v,j));
1734 466084 : ZincrementalGS(x, L, B, nx); /* overwrite last */
1735 1266596 : for (k = nx-1; k >= 1; k--) ZRED(nx,k, x,L,gel(B,k+1));
1736 466029 : gel(V,j) = gel(x,nx);
1737 : }
1738 227020 : return V;
1739 : }
1740 : GEN
1741 227039 : ZM_reducemodmatrix(GEN v, GEN y) {
1742 227039 : pari_sp av = avma;
1743 227039 : return gerepilecopy(av, ZM_reducemodmatrix_i(v,y));
1744 : }
1745 :
1746 : GEN
1747 97928 : ZC_reducemodlll(GEN x,GEN y)
1748 : {
1749 97928 : pari_sp av = avma;
1750 97928 : GEN z = ZC_reducemodmatrix(x, ZM_lll(y, 0.75, LLL_INPLACE));
1751 97932 : return gerepilecopy(av, z);
1752 : }
1753 : GEN
1754 0 : ZM_reducemodlll(GEN x,GEN y)
1755 : {
1756 0 : pari_sp av = avma;
1757 0 : GEN z = ZM_reducemodmatrix(x, ZM_lll(y, 0.75, LLL_INPLACE));
1758 0 : return gerepilecopy(av, z);
1759 : }
|