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 1655377 : check_ZV(GEN x, long l)
20 : {
21 : long i;
22 11173952 : for (i=1; i<l; i++)
23 9518623 : if (typ(gel(x,i)) != t_INT) return 0;
24 1655329 : return 1;
25 : }
26 : void
27 1218857 : RgV_check_ZV(GEN A, const char *s)
28 : {
29 1218857 : if (!RgV_is_ZV(A)) pari_err_TYPE(stack_strcat(s," [integer vector]"), A);
30 1218851 : }
31 : void
32 557677 : RgM_check_ZM(GEN A, const char *s)
33 : {
34 557677 : long n = lg(A);
35 557677 : if (n != 1)
36 : {
37 557563 : long j, m = lgcols(A);
38 2212892 : for (j=1; j<n; j++)
39 1655377 : if (!check_ZV(gel(A,j), m))
40 48 : pari_err_TYPE(stack_strcat(s," [integer matrix]"), A);
41 : }
42 557629 : }
43 :
44 : /* assume m > 1 */
45 : static long
46 93286699 : ZV_max_lg_i(GEN x, long m)
47 : {
48 93286699 : long i, l = lgefint(gel(x,1));
49 752602088 : for (i = 2; i < m; i++) l = maxss(l, lgefint(gel(x,i)));
50 93286699 : return l;
51 : }
52 : long
53 9120 : ZV_max_lg(GEN x)
54 : {
55 9120 : long m = lg(x);
56 9120 : return m == 1? 2: ZV_max_lg_i(x, m);
57 : }
58 :
59 : /* assume n > 1 and m > 1 */
60 : static long
61 22981835 : ZM_max_lg_i(GEN x, long n, long m)
62 : {
63 22981835 : long j, l = ZV_max_lg_i(gel(x,1), m);
64 93277579 : for (j = 2; j < n; j++) l = maxss(l, ZV_max_lg_i(gel(x,j), m));
65 22981835 : return l;
66 : }
67 : long
68 19317 : ZM_max_lg(GEN x)
69 : {
70 19317 : long n = lg(x), m;
71 19317 : if (n == 1) return 2;
72 19317 : 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 3639 : ZM_supnorm(GEN x)
108 : {
109 3639 : long i, j, h, lx = lg(x);
110 3639 : GEN s = gen_0;
111 3639 : if (lx == 1) return gen_1;
112 3639 : h = lgcols(x);
113 22407 : for (j=1; j<lx; j++)
114 : {
115 18768 : GEN xj = gel(x,j);
116 252600 : for (i=1; i<h; i++)
117 : {
118 233832 : GEN c = gel(xj,i);
119 233832 : if (abscmpii(c, s) > 0) s = c;
120 : }
121 : }
122 3639 : 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 1911024 : ZV_zc_mul(GEN x, GEN y)
151 : {
152 1911024 : long j, l = lg(x);
153 1911024 : pari_sp av = avma;
154 1911024 : GEN s = mulis(gel(x,1),y[1]);
155 43108284 : for (j=2; j<l; j++)
156 41197260 : if (y[j]) s = addii(s, mulis(gel(x,j),y[j]));
157 1911024 : return gerepileuptoint(av,s);
158 : }
159 :
160 : /* x nonempty ZM, y a compatible zc (dimension > 0). */
161 : static GEN
162 16524656 : ZM_zc_mul_i(GEN x, GEN y, long c, long l)
163 : {
164 : long i, j;
165 16524656 : GEN z = cgetg(l,t_COL);
166 :
167 104354230 : for (i=1; i<l; i++)
168 : {
169 87829574 : pari_sp av = avma;
170 87829574 : GEN s = mulis(gcoeff(x,i,1),y[1]);
171 1003455880 : for (j=2; j<c; j++)
172 915626306 : if (y[j]) s = addii(s, mulis(gcoeff(x,i,j),y[j]));
173 87829574 : gel(z,i) = gerepileuptoint(av,s);
174 : }
175 16524656 : return z;
176 : }
177 : GEN
178 14755614 : ZM_zc_mul(GEN x, GEN y) {
179 14755614 : long l = lg(x);
180 14755614 : if (l == 1) return cgetg(1, t_COL);
181 14755614 : 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 1488 : zv_ZM_mul(GEN x, GEN y) {
187 1488 : long i,j, lx = lg(x), ly = lg(y);
188 : GEN z;
189 1488 : if (lx == 1) return zerovec(ly-1);
190 1488 : z = cgetg(ly,t_VEC);
191 3468 : for (j=1; j<ly; j++)
192 : {
193 1980 : pari_sp av = avma;
194 1980 : GEN s = mulsi(x[1], gcoeff(y,1,j));
195 3420 : for (i=2; i<lx; i++)
196 1440 : if (x[i]) s = addii(s, mulsi(x[i], gcoeff(y,i,j)));
197 1980 : gel(z,j) = gerepileuptoint(av,s);
198 : }
199 1488 : return z;
200 : }
201 :
202 : /* x ZM, y a compatible zm (dimension > 0). */
203 : GEN
204 836248 : ZM_zm_mul(GEN x, GEN y)
205 : {
206 836248 : long j, c, l = lg(x), ly = lg(y);
207 836248 : GEN z = cgetg(ly, t_MAT);
208 836248 : if (l == 1) return z;
209 836242 : c = lgcols(x);
210 2605284 : for (j = 1; j < ly; j++) gel(z,j) = ZM_zc_mul_i(x, gel(y,j), l,c);
211 836242 : 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 538328 : 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 538328 : long min_d = minss(da, db), min_e = minss(ea, eb), i, j;
235 538328 : GEN M = cgetg(n + 1, t_MAT), C;
236 :
237 4415961 : for (j = 1; j <= min_e; j++) {
238 3877633 : gel(M, j) = C = cgetg(m + 1, t_COL);
239 74934615 : for (i = 1; i <= min_d; i++)
240 71056982 : gel(C, i) = addii(gcoeff(A, ma + i, na + j),
241 71056982 : gcoeff(B, mb + i, nb + j));
242 3939252 : for (; i <= da; i++) gel(C, i) = gcoeff(A, ma + i, na + j);
243 3877633 : for (; i <= db; i++) gel(C, i) = gcoeff(B, mb + i, nb + j);
244 3877633 : for (; i <= m; i++) gel(C, i) = gen_0;
245 : }
246 593727 : for (; j <= ea; j++) {
247 55399 : gel(M, j) = C = cgetg(m + 1, t_COL);
248 192775 : for (i = 1; i <= da; i++) gel(C, i) = gcoeff(A, ma + i, na + j);
249 55399 : for (; i <= m; i++) gel(C, i) = gen_0;
250 : }
251 538328 : 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 538328 : for (; j <= n; j++) gel(M, j) = zerocol(m);
257 538328 : 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 471037 : 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 471037 : long min_d = minss(da, db), min_e = minss(ea, eb), i, j;
268 471037 : GEN M = cgetg(n + 1, t_MAT), C;
269 :
270 3881312 : for (j = 1; j <= min_e; j++) {
271 3410275 : gel(M, j) = C = cgetg(m + 1, t_COL);
272 67155889 : for (i = 1; i <= min_d; i++)
273 63745614 : gel(C, i) = subii(gcoeff(A, ma + i, na + j),
274 63745614 : gcoeff(B, mb + i, nb + j));
275 3463892 : for (; i <= da; i++) gel(C, i) = gcoeff(A, ma + i, na + j);
276 3492719 : for (; i <= db; i++) gel(C, i) = negi(gcoeff(B, mb + i, nb + j));
277 3410275 : for (; i <= m; i++) gel(C, i) = gen_0;
278 : }
279 471037 : 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 503261 : for (; j <= eb; j++) {
285 32224 : gel(M, j) = C = cgetg(m + 1, t_COL);
286 118173 : for (i = 1; i <= db; i++) gel(C, i) = negi(gcoeff(B, mb + i, nb + j));
287 32224 : for (; i <= m; i++) gel(C, i) = gen_0;
288 : }
289 503261 : for (; j <= n; j++) gel(M, j) = zerocol(m);
290 471037 : 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 67291 : ZM_mul_sw(GEN A, GEN B, long m, long n, long p)
298 : {
299 67291 : pari_sp av = avma;
300 67291 : long m1 = (m + 1)/2, m2 = m/2,
301 67291 : n1 = (n + 1)/2, n2 = n/2,
302 67291 : 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 67291 : T2 = subtract_slices(n1, p2, B, 0, n1, p1, p2, B, n1, n2, p1, p2);
309 67291 : S1 = subtract_slices(m2, n1, A, m1, m2, 0, n1, A, 0, m2, 0, n1);
310 67291 : M2 = ZM_mul_i(S1, T2, m2 + 1, n1 + 1, p2 + 1);
311 67291 : if (gc_needed(av, 1))
312 0 : gerepileall(av, 2, &T2, &M2); /* destroy S1 */
313 67291 : T3 = subtract_slices(n1, p1, T2, 0, n1, 0, p2, B, 0, n1, 0, p1);
314 67291 : if (gc_needed(av, 1))
315 0 : gerepileall(av, 2, &M2, &T3); /* destroy T2 */
316 67291 : S2 = add_slices(m2, n1, A, m1, m2, 0, n1, A, m1, m2, n1, n2);
317 67291 : T1 = subtract_slices(n1, p1, B, 0, n1, p1, p2, B, 0, n1, 0, p2);
318 67291 : M3 = ZM_mul_i(S2, T1, m2 + 1, n1 + 1, p2 + 1);
319 67291 : if (gc_needed(av, 1))
320 0 : gerepileall(av, 4, &M2, &T3, &S2, &M3); /* destroy T1 */
321 67291 : S3 = subtract_slices(m1, n1, S2, 0, m2, 0, n1, A, 0, m1, 0, n1);
322 67291 : if (gc_needed(av, 1))
323 0 : gerepileall(av, 4, &M2, &T3, &M3, &S3); /* destroy S2 */
324 67291 : A11 = matslice(A, 1, m1, 1, n1);
325 67291 : B11 = matslice(B, 1, n1, 1, p1);
326 67291 : M1 = ZM_mul_i(A11, B11, m1 + 1, n1 + 1, p1 + 1);
327 67291 : if (gc_needed(av, 1))
328 0 : gerepileall(av, 5, &M2, &T3, &M3, &S3, &M1); /* destroy A11, B11 */
329 67291 : A12 = matslice(A, 1, m1, n1 + 1, n);
330 67291 : B21 = matslice(B, n1 + 1, n, 1, p1);
331 67291 : M4 = ZM_mul_i(A12, B21, m1 + 1, n2 + 1, p1 + 1);
332 67291 : if (gc_needed(av, 1))
333 0 : gerepileall(av, 6, &M2, &T3, &M3, &S3, &M1, &M4); /* destroy A12, B21 */
334 67291 : C11 = add_slices(m1, p1, M1, 0, m1, 0, p1, M4, 0, m1, 0, p1);
335 67291 : if (gc_needed(av, 1))
336 0 : gerepileall(av, 6, &M2, &T3, &M3, &S3, &M1, &C11); /* destroy M4 */
337 67291 : M5 = ZM_mul_i(S3, T3, m1 + 1, n1 + 1, p1 + 1);
338 67291 : S4 = subtract_slices(m1, n2, A, 0, m1, n1, n2, S3, 0, m1, 0, n2);
339 67291 : if (gc_needed(av, 1))
340 5 : gerepileall(av, 7, &M2, &T3, &M3, &M1, &C11, &M5, &S4); /* destroy S3 */
341 67291 : T4 = add_slices(n2, p1, B, n1, n2, 0, p1, T3, 0, n2, 0, p1);
342 67291 : if (gc_needed(av, 1))
343 0 : gerepileall(av, 7, &M2, &M3, &M1, &C11, &M5, &S4, &T4); /* destroy T3 */
344 67291 : V1 = subtract_slices(m1, p1, M1, 0, m1, 0, p1, M5, 0, m1, 0, p1);
345 67291 : if (gc_needed(av, 1))
346 0 : gerepileall(av, 6, &M2, &M3, &S4, &T4, &C11, &V1); /* destroy M1, M5 */
347 67291 : B22 = matslice(B, n1 + 1, n, p1 + 1, p);
348 67291 : M6 = ZM_mul_i(S4, B22, m1 + 1, n2 + 1, p2 + 1);
349 67291 : if (gc_needed(av, 1))
350 5 : gerepileall(av, 6, &M2, &M3, &T4, &C11, &V1, &M6); /* destroy S4, B22 */
351 67291 : A22 = matslice(A, m1 + 1, m, n1 + 1, n);
352 67291 : M7 = ZM_mul_i(A22, T4, m2 + 1, n2 + 1, p1 + 1);
353 67291 : if (gc_needed(av, 1))
354 0 : gerepileall(av, 6, &M2, &M3, &C11, &V1, &M6, &M7); /* destroy A22, T4 */
355 67291 : V3 = add_slices(m1, p2, V1, 0, m1, 0, p2, M3, 0, m2, 0, p2);
356 67291 : C12 = add_slices(m1, p2, V3, 0, m1, 0, p2, M6, 0, m1, 0, p2);
357 67291 : if (gc_needed(av, 1))
358 6 : gerepileall(av, 6, &M2, &M3, &C11, &V1, &M7, &C12); /* destroy V3, M6 */
359 67291 : V2 = add_slices(m2, p1, V1, 0, m2, 0, p1, M2, 0, m2, 0, p2);
360 67291 : if (gc_needed(av, 1))
361 0 : gerepileall(av, 5, &M3, &C11, &M7, &C12, &V2); /* destroy V1, M2 */
362 67291 : C21 = add_slices(m2, p1, V2, 0, m2, 0, p1, M7, 0, m2, 0, p1);
363 67291 : if (gc_needed(av, 1))
364 5 : gerepileall(av, 5, &M3, &C11, &C12, &V2, &C21); /* destroy M7 */
365 67291 : C22 = add_slices(m2, p2, V2, 0, m2, 0, p2, M3, 0, m2, 0, p2);
366 67291 : if (gc_needed(av, 1))
367 0 : gerepileall(av, 4, &C11, &C12, &C21, &C22); /* destroy V2, M3 */
368 67291 : C = shallowconcat(vconcat(C11, C21), vconcat(C12, C22));
369 67291 : return gerepilecopy(av, C);
370 : }
371 :
372 : /* x[i,]*y. Assume lg(x) > 1 and 0 < i < lgcols(x) */
373 : static GEN
374 495951541 : ZMrow_ZC_mul_i(GEN x, GEN y, long i, long lx)
375 : {
376 495951541 : pari_sp av = avma;
377 495951541 : GEN c = mulii(gcoeff(x,i,1), gel(y,1)), ZERO = gen_0;
378 : long k;
379 5755108098 : for (k = 2; k < lx; k++)
380 : {
381 5259156557 : GEN t = mulii(gcoeff(x,i,k), gel(y,k));
382 5259156557 : if (t != ZERO) c = addii(c, t);
383 : }
384 495951541 : return gerepileuptoint(av, c);
385 : }
386 : GEN
387 115914442 : ZMrow_ZC_mul(GEN x, GEN y, long i)
388 115914442 : { 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 64143720 : ZM_ZC_mul_i(GEN x, GEN y, long lx, long l)
393 : {
394 64143720 : GEN z = cgetg(l,t_COL);
395 : long i;
396 444180819 : for (i=1; i<l; i++) gel(z,i) = ZMrow_ZC_mul_i(x,y,i,lx);
397 64143720 : return z;
398 : }
399 :
400 : static GEN
401 11423390 : ZM_mul_classical(GEN x, GEN y, long l, long lx, long ly)
402 : {
403 : long j;
404 11423390 : GEN z = cgetg(ly, t_MAT);
405 54948809 : for (j = 1; j < ly; j++)
406 43525419 : gel(z, j) = ZM_ZC_mul_i(x, gel(y, j), lx, l);
407 11423390 : return z;
408 : }
409 :
410 : static GEN
411 925 : ZM_mul_slice(GEN A, GEN B, GEN P, GEN *mod)
412 : {
413 925 : pari_sp av = avma;
414 925 : long i, n = lg(P)-1;
415 : GEN H, T;
416 925 : 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 925 : T = ZV_producttree(P);
425 925 : A = ZM_nv_mod_tree(A, P, T);
426 925 : B = ZM_nv_mod_tree(B, P, T);
427 925 : H = cgetg(n+1, t_VEC);
428 6581 : for(i=1; i <= n; i++)
429 5656 : gel(H,i) = Flm_mul(gel(A,i),gel(B,i),P[i]);
430 925 : H = nmV_chinese_center_tree_seq(H, P, T, ZV_chinesetree(P, T));
431 925 : *mod = gmael(T, lg(T)-1, 1); return gc_all(av, 2, &H, mod);
432 : }
433 :
434 : GEN
435 925 : ZM_mul_worker(GEN P, GEN A, GEN B)
436 : {
437 925 : GEN V = cgetg(3, t_VEC);
438 925 : gel(V,1) = ZM_mul_slice(A, B, P, &gel(V,2));
439 925 : return V;
440 : }
441 :
442 : static GEN
443 703 : ZM_mul_fast(GEN A, GEN B, long lA, long lB, long sA, long sB)
444 : {
445 703 : pari_sp av = avma;
446 : forprime_t S;
447 : GEN H, worker;
448 : long h;
449 703 : if (sA == 2 || sB == 2) return zeromat(nbrows(A),lB-1);
450 693 : h = 1 + (sA + sB - 4) * BITS_IN_LONG + expu(lA-1);
451 693 : init_modular_big(&S);
452 693 : worker = snm_closure(is_entry("_ZM_mul_worker"), mkvec2(A,B));
453 693 : H = gen_crt("ZM_mul", worker, &S, NULL, h, 0, NULL,
454 : nmV_chinese_center, FpM_center);
455 693 : 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 11490681 : sw_bound(long s)
462 11490681 : { 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 17762267 : 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 14632371 : long Flm_sw_bound = 70;
472 : #else
473 3129896 : long Flm_sw_bound = 120;
474 : #endif
475 17762267 : if (l == 1) return zeromat(0, ly-1);
476 17760628 : if (lx==2 && l==2 && ly==2)
477 313585 : { retmkmat(mkcol(mulii(gcoeff(x,1,1), gcoeff(y,1,1)))); }
478 17447043 : if (lx==3 && l==3 && ly==3) return ZM2_mul(x, y);
479 11471134 : sx = ZM_max_lg_i(x, lx, l);
480 11471134 : 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 11471134 : if (lx > Flm_sw_bound && ly > Flm_sw_bound && l > Flm_sw_bound
484 713 : && sx <= 10 * sy && sy <= 10 * sx) return ZM_mul_fast(x,y, lx,ly, sx,sy);
485 :
486 11470431 : B = sw_bound(minss(sx, sy));
487 11470431 : if (l <= B || lx <= B || ly <= B)
488 11403166 : return ZM_mul_classical(x, y, l, lx, ly);
489 : else
490 67265 : return ZM_mul_sw(x, y, l - 1, lx - 1, ly - 1);
491 : }
492 :
493 : GEN
494 17408464 : ZM_mul(GEN x, GEN y)
495 : {
496 17408464 : long lx = lg(x), ly = lg(y);
497 17408464 : if (ly == 1) return cgetg(1,t_MAT);
498 17292586 : if (lx == 1) return zeromat(0, ly-1);
499 17291230 : 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 478804 : QM_mul(GEN x, GEN y)
550 : {
551 478804 : GEN dx, nx = Q_primitive_part(x, &dx);
552 478804 : GEN dy, ny = Q_primitive_part(y, &dy);
553 478804 : GEN z = ZM_mul(nx, ny);
554 478804 : if (dx || dy)
555 : {
556 406704 : GEN d = dx ? dy ? gmul(dx, dy): dx : dy;
557 406704 : if (!gequal1(d)) z = ZM_Q_mul(z, d);
558 : }
559 478804 : return z;
560 : }
561 :
562 : GEN
563 600 : QM_sqr(GEN x)
564 : {
565 600 : GEN dx, nx = Q_primitive_part(x, &dx);
566 600 : GEN z = ZM_sqr(nx);
567 600 : if (dx)
568 600 : z = ZM_Q_mul(z, gsqr(dx));
569 600 : return z;
570 : }
571 :
572 : GEN
573 108730 : QM_QC_mul(GEN x, GEN y)
574 : {
575 108730 : GEN dx, nx = Q_primitive_part(x, &dx);
576 108730 : GEN dy, ny = Q_primitive_part(y, &dy);
577 108730 : GEN z = ZM_ZC_mul(nx, ny);
578 108730 : if (dx || dy)
579 : {
580 108712 : GEN d = dx ? dy ? gmul(dx, dy): dx : dy;
581 108712 : if (!gequal1(d)) z = ZC_Q_mul(z, d);
582 : }
583 108730 : 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 18 : ZM_mul_diag(GEN m, GEN d)
611 : {
612 : long j, l;
613 18 : GEN y = cgetg_copy(m, &l);
614 66 : for (j=1; j<l; j++)
615 : {
616 48 : GEN c = gel(d,j);
617 48 : gel(y,j) = equali1(c)? gel(m,j): ZC_Z_mul(gel(m,j), c);
618 : }
619 18 : return y;
620 : }
621 : /* compute diagonal(d)*m, assume lg(d) = lg(m~). Shallow */
622 : GEN
623 457872 : ZM_diag_mul(GEN d, GEN m)
624 : {
625 457872 : long i, j, l = lg(d), lm = lg(m);
626 457872 : GEN y = cgetg(lm, t_MAT);
627 1304602 : for (j=1; j<lm; j++) gel(y,j) = cgetg(l, t_COL);
628 1413556 : for (i=1; i<l; i++)
629 : {
630 955684 : GEN c = gel(d,i);
631 955684 : if (equali1(c))
632 203038 : for (j=1; j<lm; j++) gcoeff(y,i,j) = gcoeff(m,i,j);
633 : else
634 2916660 : for (j=1; j<lm; j++) gcoeff(y,i,j) = mulii(gcoeff(m,i,j), c);
635 : }
636 457872 : return y;
637 : }
638 :
639 : /* assume lx > 1 is lg(x) = lg(y) */
640 : static GEN
641 16444497 : ZV_dotproduct_i(GEN x, GEN y, long lx)
642 : {
643 16444497 : pari_sp av = avma;
644 16444497 : GEN c = mulii(gel(x,1), gel(y,1));
645 : long i;
646 123634749 : for (i = 2; i < lx; i++)
647 : {
648 107190252 : GEN t = mulii(gel(x,i), gel(y,i));
649 107190252 : if (t != gen_0) c = addii(c, t);
650 : }
651 16444497 : return gerepileuptoint(av, c);
652 : }
653 :
654 : /* x~ * y, assuming result is symmetric */
655 : GEN
656 456175 : ZM_transmultosym(GEN x, GEN y)
657 : {
658 456175 : long i, j, l, ly = lg(y);
659 : GEN M;
660 456175 : if (ly == 1) return cgetg(1,t_MAT);
661 : /* lg(x) = ly */
662 456175 : l = lgcols(y); /* = lgcols(x) */
663 456175 : M = cgetg(ly, t_MAT);
664 2321698 : for (i=1; i<ly; i++)
665 : {
666 1865523 : GEN xi = gel(x,i), c = cgetg(ly,t_COL);
667 1865523 : gel(M,i) = c;
668 6093269 : for (j=1; j<i; j++)
669 4227746 : gcoeff(M,i,j) = gel(c,j) = ZV_dotproduct_i(xi,gel(y,j),l);
670 1865523 : gel(c,i) = ZV_dotproduct_i(xi,gel(y,i),l);
671 : }
672 456175 : return M;
673 : }
674 : /* x~ * y */
675 : GEN
676 1962 : ZM_transmul(GEN x, GEN y)
677 : {
678 1962 : long i, j, l, lx, ly = lg(y);
679 : GEN M;
680 1962 : if (ly == 1) return cgetg(1,t_MAT);
681 1962 : lx = lg(x);
682 1962 : l = lgcols(y);
683 1962 : if (lgcols(x) != l) pari_err_OP("operation 'ZM_transmul'", x,y);
684 1962 : M = cgetg(ly, t_MAT);
685 5994 : for (i=1; i<ly; i++)
686 : {
687 4032 : GEN yi = gel(y,i), c = cgetg(lx,t_COL);
688 4032 : gel(M,i) = c;
689 10482 : for (j=1; j<lx; j++) gel(c,j) = ZV_dotproduct_i(yi,gel(x,j),l);
690 : }
691 1962 : 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 20916 : ZM_sqr_i(GEN x, long l)
699 : {
700 : long s;
701 20916 : if (l == 2) { retmkmat(mkcol(sqri(gcoeff(x,1,1)))); }
702 20916 : if (l == 3) return ZM2_sqr(x);
703 20250 : s = ZM_max_lg_i(x, l, l);
704 20250 : if (l > 70) return ZM_sqr_fast(x, l, s);
705 20250 : if (l <= sw_bound(s))
706 20224 : return ZM_mul_classical(x, x, l, l, l);
707 : else
708 26 : return ZM_mul_sw(x, x, l - 1, l - 1, l - 1);
709 : }
710 :
711 : GEN
712 20916 : ZM_sqr(GEN x)
713 : {
714 20916 : long lx=lg(x);
715 20916 : if (lx==1) return cgetg(1,t_MAT);
716 20916 : return ZM_sqr_i(x, lx);
717 : }
718 : GEN
719 20692083 : ZM_ZC_mul(GEN x, GEN y)
720 : {
721 20692083 : long lx = lg(x);
722 20692083 : return lx==1? cgetg(1,t_COL): ZM_ZC_mul_i(x, y, lx, lgcols(x));
723 : }
724 :
725 : GEN
726 3150823 : ZC_Z_div(GEN x, GEN c)
727 14928448 : { pari_APPLY_type(t_COL, Qdivii(gel(x,i), c)) }
728 :
729 : GEN
730 16554 : ZM_Z_div(GEN x, GEN c)
731 189560 : { pari_APPLY_same(ZC_Z_div(gel(x, i), c)) }
732 :
733 : GEN
734 2337379 : ZC_Q_mul(GEN A, GEN z)
735 : {
736 2337379 : pari_sp av = avma;
737 2337379 : long i, l = lg(A);
738 : GEN d, n, Ad, B, u;
739 2337379 : if (typ(z)==t_INT) return ZC_Z_mul(A,z);
740 2333629 : n = gel(z, 1); d = gel(z, 2);
741 2333629 : Ad = FpC_red(A, d);
742 2333629 : u = gcdii(d, FpV_factorback(Ad, NULL, d));
743 2333629 : B = cgetg(l, t_COL);
744 2333629 : if (equali1(u))
745 : {
746 357391 : for(i=1; i<l; i++)
747 302515 : gel(B, i) = mkfrac(mulii(n, gel(A,i)), d);
748 : } else
749 : {
750 16130602 : for(i=1; i<l; i++)
751 : {
752 13851849 : GEN di = gcdii(gel(Ad, i), u), ni = mulii(n, diviiexact(gel(A,i), di));
753 13851849 : if (equalii(d, di))
754 9591976 : gel(B, i) = ni;
755 : else
756 4259873 : gel(B, i) = mkfrac(ni, diviiexact(d, di));
757 : }
758 : }
759 2333629 : return gerepilecopy(av, B);
760 : }
761 :
762 : GEN
763 935967 : ZM_Q_mul(GEN x, GEN z)
764 : {
765 935967 : if (typ(z)==t_INT) return ZM_Z_mul(x,z);
766 2793972 : pari_APPLY_same(ZC_Q_mul(gel(x, i), z));
767 : }
768 :
769 : long
770 168342102 : zv_dotproduct(GEN x, GEN y)
771 : {
772 168342102 : long i, lx = lg(x);
773 : ulong c;
774 168342102 : if (lx == 1) return 0;
775 168342102 : c = uel(x,1)*uel(y,1);
776 2612022132 : for (i = 2; i < lx; i++)
777 2443680030 : c += uel(x,i)*uel(y,i);
778 168342102 : return c;
779 : }
780 :
781 : GEN
782 197682 : ZV_ZM_mul(GEN x, GEN y)
783 : {
784 197682 : long i, lx = lg(x), ly = lg(y);
785 : GEN z;
786 197682 : if (lx == 1) return zerovec(ly-1);
787 197580 : z = cgetg(ly, t_VEC);
788 756054 : for (i = 1; i < ly; i++) gel(z,i) = ZV_dotproduct_i(x, gel(y,i), lx);
789 197580 : 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 5710767 : ZV_dotsquare(GEN x)
809 : {
810 : long i, lx;
811 : pari_sp av;
812 : GEN z;
813 5710767 : lx = lg(x);
814 5710767 : if (lx == 1) return gen_0;
815 5710767 : av = avma; z = sqri(gel(x,1));
816 22610498 : for (i=2; i<lx; i++) z = addii(z, sqri(gel(x,i)));
817 5710767 : return gerepileuptoint(av,z);
818 : }
819 :
820 : GEN
821 13885725 : ZV_dotproduct(GEN x,GEN y)
822 : {
823 : long lx;
824 13885725 : if (x == y) return ZV_dotsquare(x);
825 9786304 : lx = lg(x);
826 9786304 : if (lx == 1) return gen_0;
827 9786304 : return ZV_dotproduct_i(x, y, lx);
828 : }
829 :
830 : static GEN
831 240 : _ZM_mul(void *data /*ignored*/, GEN x, GEN y)
832 240 : { (void)data; return ZM_mul(x,y); }
833 : static GEN
834 19998 : _ZM_sqr(void *data /*ignored*/, GEN x)
835 19998 : { (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 19548 : ZM_powu(GEN x, ulong n)
847 : {
848 19548 : pari_sp av = avma;
849 19548 : if (!n) return matid(lg(x)-1);
850 19548 : 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 30919762 : ZC_add_i(GEN x, GEN y, long lx)
859 : {
860 30919762 : GEN A = cgetg(lx, t_COL);
861 : long i;
862 415279790 : for (i=1; i<lx; i++) gel(A,i) = addii(gel(x,i), gel(y,i));
863 30919762 : return A;
864 : }
865 : GEN
866 22320401 : ZC_add(GEN x, GEN y) { return ZC_add_i(x, y, lg(x)); }
867 : GEN
868 313178 : ZC_Z_add(GEN x, GEN y)
869 : {
870 313178 : long k, lx = lg(x);
871 313178 : GEN z = cgetg(lx, t_COL);
872 313178 : if (lx == 1) pari_err_TYPE2("+",x,y);
873 313178 : gel(z,1) = addii(y,gel(x,1));
874 2152404 : for (k = 2; k < lx; k++) gel(z,k) = icopy(gel(x,k));
875 313178 : return z;
876 : }
877 :
878 : static GEN
879 7788500 : ZC_sub_i(GEN x, GEN y, long lx)
880 : {
881 : long i;
882 7788500 : GEN A = cgetg(lx, t_COL);
883 55318981 : for (i=1; i<lx; i++) gel(A,i) = subii(gel(x,i), gel(y,i));
884 7788500 : return A;
885 : }
886 : GEN
887 7507254 : 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 467417 : Z_ZC_sub(GEN a, GEN x)
900 : {
901 467417 : long k, lx = lg(x);
902 467417 : GEN z = cgetg(lx, t_COL);
903 467417 : if (lx == 1) pari_err_TYPE2("-",a,x);
904 467417 : gel(z,1) = subii(a, gel(x,1));
905 1298500 : for (k = 2; k < lx; k++) gel(z,k) = negi(gel(x,k));
906 467417 : return z;
907 : }
908 :
909 : GEN
910 723409 : ZM_add(GEN x, GEN y)
911 : {
912 723409 : long lx = lg(x), l, j;
913 : GEN z;
914 723409 : if (lx == 1) return cgetg(1, t_MAT);
915 689305 : z = cgetg(lx, t_MAT); l = lgcols(x);
916 9288666 : for (j = 1; j < lx; j++) gel(z,j) = ZC_add_i(gel(x,j), gel(y,j), l);
917 689305 : return z;
918 : }
919 : GEN
920 133330 : ZM_sub(GEN x, GEN y)
921 : {
922 133330 : long lx = lg(x), l, j;
923 : GEN z;
924 133330 : if (lx == 1) return cgetg(1, t_MAT);
925 99226 : z = cgetg(lx, t_MAT); l = lgcols(x);
926 380472 : for (j = 1; j < lx; j++) gel(z,j) = ZC_sub_i(gel(x,j), gel(y,j), l);
927 99226 : return z;
928 : }
929 : /********************************************************************/
930 : /** **/
931 : /** LINEAR COMBINATION **/
932 : /** **/
933 : /********************************************************************/
934 : /* return X/c assuming division is exact */
935 : GEN
936 4720199 : ZC_Z_divexact(GEN x, GEN c)
937 67264826 : { pari_APPLY_type(t_COL, diviiexact(gel(x,i), c)) }
938 : GEN
939 1938 : ZC_divexactu(GEN x, ulong c)
940 9750 : { pari_APPLY_type(t_COL, diviuexact(gel(x,i), c)) }
941 :
942 : GEN
943 421162 : ZM_Z_divexact(GEN x, GEN c)
944 3364613 : { pari_APPLY_same(ZC_Z_divexact(gel(x,i), c)) }
945 :
946 : GEN
947 378 : ZM_divexactu(GEN x, ulong c)
948 2316 : { pari_APPLY_same(ZC_divexactu(gel(x,i), c)) }
949 :
950 : GEN
951 29738485 : ZC_Z_mul(GEN x, GEN c)
952 : {
953 29738485 : if (!signe(c)) return zerocol(lg(x)-1);
954 28592098 : if (is_pm1(c)) return (signe(c) > 0)? ZC_copy(x): ZC_neg(x);
955 217550714 : pari_APPLY_type(t_COL, mulii(gel(x,i), c))
956 : }
957 :
958 : GEN
959 53318 : ZC_z_mul(GEN x, long c)
960 : {
961 53318 : if (!c) return zerocol(lg(x)-1);
962 44850 : if (c == 1) return ZC_copy(x);
963 40637 : if (c ==-1) return ZC_neg(x);
964 408278 : pari_APPLY_type(t_COL, mulsi(c, gel(x,i)))
965 : }
966 :
967 : GEN
968 24265 : zv_z_mul(GEN x, long n)
969 369655 : { 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 2805015 : ZM_Z_mul(GEN X, GEN c)
995 : {
996 2805015 : long i, j, h, l = lg(X);
997 : GEN A;
998 2805015 : if (l == 1) return cgetg(1, t_MAT);
999 2805015 : h = lgcols(X);
1000 2805015 : if (!signe(c)) return zeromat(h-1, l-1);
1001 2766140 : if (is_pm1(c)) return (signe(c) > 0)? ZM_copy(X): ZM_neg(X);
1002 2423223 : A = cgetg(l, t_MAT);
1003 13521565 : for (j = 1; j < l; j++)
1004 : {
1005 11098342 : GEN a = cgetg(h, t_COL), x = gel(X, j);
1006 188757305 : for (i = 1; i < h; i++) gel(a,i) = mulii(c, gel(x,i));
1007 11098342 : gel(A,j) = a;
1008 : }
1009 2423223 : return A;
1010 : }
1011 : void
1012 62802264 : ZC_lincomb1_inplace_i(GEN X, GEN Y, GEN v, long n)
1013 : {
1014 : long i;
1015 1391214464 : for (i = n; i; i--) gel(X,i) = addmulii_inplace(gel(X,i), gel(Y,i), v);
1016 62802264 : }
1017 : /* X <- X + v Y (elementary col operation) */
1018 : void
1019 53423498 : ZC_lincomb1_inplace(GEN X, GEN Y, GEN v)
1020 : {
1021 53423498 : if (lgefint(v) != 2) return ZC_lincomb1_inplace_i(X, Y, v, lg(X)-1);
1022 : }
1023 : void
1024 25342998 : Flc_lincomb1_inplace(GEN X, GEN Y, ulong v, ulong q)
1025 : {
1026 : long i;
1027 25342998 : if (!v) return; /* v = 0 */
1028 634932697 : 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 13190696 : ZC_lincomb1(GEN v, GEN x, GEN y)
1034 104542913 : { pari_APPLY_type(t_COL, addmulii(gel(x,i), gel(y,i), v)) }
1035 :
1036 : /* -X + vY */
1037 : static GEN
1038 635055 : ZC_lincomb_1(GEN v, GEN x, GEN y)
1039 3930071 : { 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 28537858 : ZC_lincomb(GEN u, GEN v, GEN X, GEN Y)
1044 : {
1045 : long su, sv;
1046 : GEN A;
1047 :
1048 28537858 : su = signe(u); if (!su) return ZC_Z_mul(Y, v);
1049 28537042 : sv = signe(v); if (!sv) return ZC_Z_mul(X, u);
1050 28168610 : if (is_pm1(v))
1051 : {
1052 9503692 : if (is_pm1(u))
1053 : {
1054 8366483 : if (su != sv) A = ZC_sub(X, Y);
1055 2322085 : else A = ZC_add(X, Y);
1056 8366483 : if (su < 0) ZV_togglesign(A); /* in place but was created above */
1057 : }
1058 : else
1059 : {
1060 1137209 : if (sv > 0) A = ZC_lincomb1 (u, Y, X);
1061 520182 : else A = ZC_lincomb_1(u, Y, X);
1062 : }
1063 : }
1064 18664918 : else if (is_pm1(u))
1065 : {
1066 12688542 : if (su > 0) A = ZC_lincomb1 (v, X, Y);
1067 114873 : else A = ZC_lincomb_1(v, X, Y);
1068 : }
1069 : else
1070 : { /* not cgetg_copy: x may be a t_VEC */
1071 5976376 : long i, lx = lg(X);
1072 5976376 : A = cgetg(lx,t_COL);
1073 38354048 : for (i=1; i<lx; i++) gel(A,i) = lincombii(u,v,gel(X,i),gel(Y,i));
1074 : }
1075 28168610 : return A;
1076 : }
1077 :
1078 : /********************************************************************/
1079 : /** **/
1080 : /** CONVERSIONS **/
1081 : /** **/
1082 : /********************************************************************/
1083 : GEN
1084 399975 : ZV_to_nv(GEN x)
1085 758244 : { pari_APPLY_ulong(itou(gel(x,i))) }
1086 :
1087 : GEN
1088 134200 : zm_to_ZM(GEN x)
1089 781692 : { pari_APPLY_type(t_MAT, zc_to_ZC(gel(x,i))) }
1090 :
1091 : GEN
1092 108 : zmV_to_ZMV(GEN x)
1093 678 : { 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 876 : ZM_to_zm(GEN x)
1098 14976 : { pari_APPLY_same(ZV_to_zv(gel(x,i))) }
1099 :
1100 : GEN
1101 314268 : zv_to_Flv(GEN x, ulong p)
1102 4644696 : { pari_APPLY_ulong(umodsu(x[i], p)) }
1103 :
1104 : GEN
1105 19452 : zm_to_Flm(GEN x, ulong p)
1106 301500 : { pari_APPLY_same(zv_to_Flv(gel(x,i),p)) }
1107 :
1108 : GEN
1109 42 : ZMV_to_zmV(GEN x)
1110 342 : { pari_APPLY_type(t_VEC, ZM_to_zm(gel(x,i))) }
1111 :
1112 : /********************************************************************/
1113 : /** **/
1114 : /** COPY, NEGATION **/
1115 : /** **/
1116 : /********************************************************************/
1117 : GEN
1118 14922967 : ZC_copy(GEN x)
1119 : {
1120 14922967 : long i, lx = lg(x);
1121 14922967 : GEN y = cgetg(lx, t_COL);
1122 115886409 : for (i=1; i<lx; i++)
1123 : {
1124 100963442 : GEN c = gel(x,i);
1125 100963442 : gel(y,i) = lgefint(c) == 2? gen_0: icopy(c);
1126 : }
1127 14922967 : return y;
1128 : }
1129 :
1130 : GEN
1131 560649 : ZM_copy(GEN x)
1132 5128847 : { pari_APPLY_same(ZC_copy(gel(x,i))) }
1133 :
1134 : void
1135 296151 : ZV_neg_inplace(GEN M)
1136 : {
1137 296151 : long l = lg(M);
1138 1116333 : while (--l > 0) gel(M,l) = negi(gel(M,l));
1139 296151 : }
1140 : GEN
1141 5437388 : ZC_neg(GEN x)
1142 31867894 : { pari_APPLY_type(t_COL, negi(gel(x,i))) }
1143 :
1144 : GEN
1145 44308 : zv_neg(GEN x)
1146 567171 : { pari_APPLY_long(-x[i]) }
1147 : GEN
1148 108 : zv_neg_inplace(GEN M)
1149 : {
1150 108 : long l = lg(M);
1151 372 : while (--l > 0) M[l] = -M[l];
1152 108 : return M;
1153 : }
1154 : GEN
1155 66 : zv_abs(GEN x)
1156 4668 : { pari_APPLY_ulong(labs(x[i])) }
1157 : GEN
1158 1476475 : ZM_neg(GEN x)
1159 4420448 : { pari_APPLY_same(ZC_neg(gel(x,i))) }
1160 :
1161 : void
1162 4316688 : ZV_togglesign(GEN M)
1163 : {
1164 4316688 : long l = lg(M);
1165 64711148 : while (--l > 0) togglesign_safe(&gel(M,l));
1166 4316688 : }
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 9214976 : ZC_hnfdivrem_i(GEN x, GEN y, GEN *Q, GEN (*div)(GEN,GEN))
1182 : {
1183 9214976 : long i, l = lg(x);
1184 9214976 : pari_sp av = avma;
1185 :
1186 9214976 : if (Q) *Q = cgetg(l,t_COL);
1187 9214976 : if (l == 1) return cgetg(1,t_COL);
1188 52928416 : for (i = l-1; i>0; i--)
1189 : {
1190 43713440 : GEN q = div(gel(x,i), gcoeff(y,i,i));
1191 43713440 : if (signe(q)) x = ZC_lincomb(gen_1, negi(q), x, gel(y,i));
1192 43713440 : if (Q) gel(*Q, i) = q;
1193 : }
1194 9214976 : if (avma == av) return ZC_copy(x);
1195 6864324 : if (!Q) return gerepileupto(av, x);
1196 52260 : gerepileall(av, 2, &x, Q); return x;
1197 : }
1198 : GEN
1199 8736381 : ZC_hnfdivrem(GEN x, GEN y, GEN *Q)
1200 8736381 : { return ZC_hnfdivrem_i(x, y, Q, diviiround); }
1201 : GEN
1202 24 : ZC_modhnf(GEN x, GEN y, GEN *Q)
1203 24 : { 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 389619 : ZM_hnfdivrem_i(GEN x, GEN y, GEN *Q, GEN (*div)(GEN,GEN))
1208 : {
1209 : long l, i;
1210 389619 : GEN R = cgetg_copy(x, &l);
1211 389619 : if (Q)
1212 : {
1213 109410 : GEN q = cgetg(l, t_MAT); *Q = q;
1214 161658 : for (i = 1; i < l; i++)
1215 52248 : gel(R,i) = ZC_hnfdivrem_i(gel(x,i),y,&gel(q,i),div);
1216 : }
1217 : else
1218 706532 : for (i = 1; i < l; i++)
1219 426323 : gel(R,i) = ZC_hnfdivrem_i(gel(x,i),y,NULL,div);
1220 389619 : return R;
1221 : }
1222 : GEN
1223 389607 : ZM_hnfdivrem(GEN x, GEN y, GEN *Q)
1224 389607 : { return ZM_hnfdivrem_i(x, y, Q, diviiround); }
1225 : GEN
1226 12 : ZM_modhnf(GEN x, GEN y, GEN *Q)
1227 12 : { return ZM_hnfdivrem_i(x, y, Q, truedivii); }
1228 :
1229 : static GEN
1230 6 : ZV_ZV_divrem(GEN x, GEN y, GEN *pQ)
1231 : {
1232 6 : long i, l = lg(x), tx = typ(x);
1233 : GEN Q, R;
1234 :
1235 6 : 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 36 : RgM_issquare(GEN x) { long l = lg(x); return l == 1 || lg(gel(x,1)) == l; }
1248 : static void
1249 84 : matmodhnf_check(GEN x)
1250 : {
1251 84 : switch(typ(x))
1252 : {
1253 36 : case t_VEC: case t_COL:
1254 36 : if (!RgV_is_ZV(x)) pari_err_TYPE("matmodhnf", x);
1255 36 : break;
1256 48 : case t_MAT:
1257 48 : if (!RgM_is_ZM(x)) pari_err_TYPE("matmodhnf", x);
1258 48 : break;
1259 0 : default: pari_err_TYPE("matmodhnf", x);
1260 : }
1261 84 : }
1262 : GEN
1263 42 : matmodhnf(GEN x, GEN y, GEN *Q)
1264 : {
1265 42 : long tx = typ(x), ty = typ(y), ly, lx;
1266 42 : matmodhnf_check(x); lx = lg(x);
1267 42 : matmodhnf_check(y); ly = lg(y);
1268 42 : if (ty == t_MAT && !RgM_issquare(y)) pari_err_TYPE("matmodhnf", y);
1269 42 : 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 42 : if (is_vec_t(ty))
1276 6 : return tx == t_MAT? ZM_ZV_divrem(x, y, Q): ZV_ZV_divrem(x, y, Q);
1277 : /* ty = t_MAT */
1278 36 : if (tx == t_MAT) return ZM_modhnf(x, y, Q);
1279 24 : x = ZC_modhnf(x, y, Q);
1280 24 : if (tx == t_VEC) { settyp(x, tx); if (Q) settyp(*Q, tx); }
1281 24 : return x;
1282 : }
1283 :
1284 : /********************************************************************/
1285 : /** **/
1286 : /** TESTS **/
1287 : /** **/
1288 : /********************************************************************/
1289 : int
1290 19802911 : zv_equal0(GEN V)
1291 : {
1292 19802911 : long l = lg(V);
1293 32211073 : while (--l > 0)
1294 26404564 : if (V[l]) return 0;
1295 5806509 : return 1;
1296 : }
1297 :
1298 : int
1299 12337302 : ZV_equal0(GEN V)
1300 : {
1301 12337302 : long l = lg(V);
1302 21841427 : while (--l > 0)
1303 21475177 : if (signe(gel(V,l))) return 0;
1304 366250 : return 1;
1305 : }
1306 : int
1307 13927 : ZMrow_equal0(GEN V, long i)
1308 : {
1309 13927 : long l = lg(V);
1310 21604 : while (--l > 0)
1311 18598 : if (signe(gcoeff(V,i,l))) return 0;
1312 3006 : return 1;
1313 : }
1314 :
1315 : static int
1316 5340047 : ZV_equal_lg(GEN V, GEN W, long l)
1317 : {
1318 22114856 : while (--l > 0)
1319 17197353 : if (!equalii(gel(V,l), gel(W,l))) return 0;
1320 4917503 : return 1;
1321 : }
1322 : int
1323 251708 : ZV_equal(GEN V, GEN W)
1324 : {
1325 251708 : long l = lg(V);
1326 251708 : if (lg(W) != l) return 0;
1327 251690 : return ZV_equal_lg(V, W, l);
1328 : }
1329 : int
1330 2855134 : ZM_equal(GEN A, GEN B)
1331 : {
1332 2855134 : long i, m, l = lg(A);
1333 2855134 : if (lg(B) != l) return 0;
1334 2855134 : if (l == 1) return 1;
1335 2855134 : m = lgcols(A);
1336 2855134 : if (lgcols(B) != m) return 0;
1337 7708081 : for (i = 1; i < l; i++)
1338 5088357 : if (!ZV_equal_lg(gel(A,i), gel(B,i), m)) return 0;
1339 2619724 : return 1;
1340 : }
1341 : int
1342 62220 : ZM_equal0(GEN A)
1343 : {
1344 62220 : long i, j, m, l = lg(A);
1345 62220 : if (l == 1) return 1;
1346 62220 : m = lgcols(A);
1347 109877 : for (j = 1; j < l; j++)
1348 2367464 : for (i = 1; i < m; i++)
1349 2319807 : if (signe(gcoeff(A,i,j))) return 0;
1350 14698 : return 1;
1351 : }
1352 : int
1353 26514600 : zv_equal(GEN V, GEN W)
1354 : {
1355 26514600 : long l = lg(V);
1356 26514600 : if (lg(W) != l) return 0;
1357 225142966 : while (--l > 0)
1358 199585002 : if (V[l] != W[l]) return 0;
1359 25557964 : return 1;
1360 : }
1361 :
1362 : int
1363 1404 : zvV_equal(GEN V, GEN W)
1364 : {
1365 1404 : long l = lg(V);
1366 1404 : if (lg(W) != l) return 0;
1367 68904 : while (--l > 0)
1368 68496 : if (!zv_equal(gel(V,l),gel(W,l))) return 0;
1369 408 : return 1;
1370 : }
1371 :
1372 : int
1373 644851 : ZM_ishnf(GEN x)
1374 : {
1375 644851 : long i,j, lx = lg(x);
1376 1937441 : for (i=1; i<lx; i++)
1377 : {
1378 1389292 : GEN xii = gcoeff(x,i,i);
1379 1389292 : if (signe(xii) <= 0) return 0;
1380 2702768 : for (j=1; j<i; j++)
1381 1345592 : if (signe(gcoeff(x,i,j))) return 0;
1382 2752928 : for (j=i+1; j<lx; j++)
1383 : {
1384 1460338 : GEN xij = gcoeff(x,i,j);
1385 1460338 : if (signe(xij)<0 || cmpii(xij,xii)>=0) return 0;
1386 : }
1387 : }
1388 548149 : return 1;
1389 : }
1390 : int
1391 563554 : ZM_isidentity(GEN x)
1392 : {
1393 563554 : long i,j, lx = lg(x);
1394 :
1395 563554 : if (lx == 1) return 1;
1396 563548 : if (lx != lgcols(x)) return 0;
1397 2753328 : for (j=1; j<lx; j++)
1398 : {
1399 2189829 : GEN c = gel(x,j);
1400 6869797 : for (i=1; i<j; )
1401 4679987 : if (signe(gel(c,i++))) return 0;
1402 : /* i = j */
1403 2189810 : if (!equali1(gel(c,i++))) return 0;
1404 6869810 : for ( ; i<lx; )
1405 4680018 : if (signe(gel(c,i++))) return 0;
1406 : }
1407 563499 : return 1;
1408 : }
1409 : int
1410 477297 : ZM_isdiagonal(GEN x)
1411 : {
1412 477297 : long i,j, lx = lg(x);
1413 477297 : if (lx == 1) return 1;
1414 477297 : if (lx != lgcols(x)) return 0;
1415 :
1416 1234014 : for (j=1; j<lx; j++)
1417 : {
1418 1036944 : GEN c = gel(x,j);
1419 1455948 : for (i=1; i<j; i++)
1420 699201 : if (signe(gel(c,i))) return 0;
1421 1733750 : for (i++; i<lx; i++)
1422 977033 : if (signe(gel(c,i))) return 0;
1423 : }
1424 197070 : return 1;
1425 : }
1426 : int
1427 137168 : ZM_isscalar(GEN x, GEN s)
1428 : {
1429 137168 : long i, j, lx = lg(x);
1430 :
1431 137168 : if (lx == 1) return 1;
1432 137168 : if (!s) s = gcoeff(x,1,1);
1433 137168 : if (equali1(s)) return ZM_isidentity(x);
1434 135920 : if (lx != lgcols(x)) return 0;
1435 611040 : for (j=1; j<lx; j++)
1436 : {
1437 527862 : GEN c = gel(x,j);
1438 2430251 : for (i=1; i<j; )
1439 1953289 : if (signe(gel(c,i++))) return 0;
1440 : /* i = j */
1441 476962 : if (!equalii(gel(c,i++), s)) return 0;
1442 2454422 : for ( ; i<lx; )
1443 1979302 : if (signe(gel(c,i++))) return 0;
1444 : }
1445 83178 : return 1;
1446 : }
1447 :
1448 : long
1449 149517 : ZC_is_ei(GEN x)
1450 : {
1451 149517 : long i, j = 0, l = lg(x);
1452 1534558 : for (i = 1; i < l; i++)
1453 : {
1454 1385042 : GEN c = gel(x,i);
1455 1385042 : long s = signe(c);
1456 1385042 : if (!s) continue;
1457 149512 : if (s < 0 || !is_pm1(c) || j) return 0;
1458 149511 : j = i;
1459 : }
1460 149516 : return j;
1461 : }
1462 :
1463 : /********************************************************************/
1464 : /** **/
1465 : /** MISCELLANEOUS **/
1466 : /** **/
1467 : /********************************************************************/
1468 : /* assume lg(x) = lg(y), x,y in Z^n */
1469 : int
1470 2697965 : ZV_cmp(GEN x, GEN y)
1471 : {
1472 2697965 : long fl,i, lx = lg(x);
1473 5465610 : for (i=1; i<lx; i++)
1474 4354909 : if (( fl = cmpii(gel(x,i), gel(y,i)) )) return fl;
1475 1110701 : return 0;
1476 : }
1477 : /* assume lg(x) = lg(y), x,y in Z^n */
1478 : int
1479 16920 : ZV_abscmp(GEN x, GEN y)
1480 : {
1481 16920 : long fl,i, lx = lg(x);
1482 46460 : for (i=1; i<lx; i++)
1483 46351 : if (( fl = abscmpii(gel(x,i), gel(y,i)) )) return fl;
1484 109 : return 0;
1485 : }
1486 :
1487 : long
1488 15817127 : zv_content(GEN x)
1489 : {
1490 15817127 : long i, s, l = lg(x);
1491 15817127 : if (l == 1) return 0;
1492 15817121 : s = labs(x[1]);
1493 36002928 : for (i = 2; i < l && s != 1; i++) s = ugcd(s, labs(x[i]));
1494 15817121 : return s;
1495 : }
1496 : GEN
1497 254909 : ZV_content(GEN x)
1498 : {
1499 254909 : long i, l = lg(x);
1500 254909 : pari_sp av = avma;
1501 : GEN c;
1502 254909 : if (l == 1) return gen_0;
1503 254909 : if (l == 2) return absi(gel(x,1));
1504 175601 : c = gel(x,1);
1505 478310 : for (i = 2; i < l; i++)
1506 : {
1507 346017 : c = gcdii(c, gel(x,i));
1508 346017 : if (is_pm1(c)) { set_avma(av); return gen_1; }
1509 : }
1510 132293 : return gerepileuptoint(av, c);
1511 : }
1512 :
1513 : GEN
1514 3388356 : ZM_det_triangular(GEN mat)
1515 : {
1516 : pari_sp av;
1517 3388356 : long i,l = lg(mat);
1518 : GEN s;
1519 :
1520 3388356 : if (l<3) return l<2? gen_1: icopy(gcoeff(mat,1,1));
1521 3019164 : av = avma; s = gcoeff(mat,1,1);
1522 8177287 : for (i=2; i<l; i++) s = mulii(s,gcoeff(mat,i,i));
1523 3019164 : return gerepileuptoint(av,s);
1524 : }
1525 :
1526 : /* assumes no overflow */
1527 : long
1528 815001 : zv_prod(GEN v)
1529 : {
1530 815001 : long n, i, l = lg(v);
1531 815001 : if (l == 1) return 1;
1532 823596 : n = v[1]; for (i = 2; i < l; i++) n *= v[i];
1533 662019 : return n;
1534 : }
1535 :
1536 : static GEN
1537 274369823 : _mulii(void *E, GEN a, GEN b)
1538 274369823 : { (void) E; return mulii(a, b); }
1539 :
1540 : /* product of ulongs */
1541 : GEN
1542 1598231 : zv_prod_Z(GEN v)
1543 : {
1544 : pari_sp av;
1545 1598231 : long k, m, n = lg(v)-1;
1546 1598231 : int stop = 0;
1547 : GEN V;
1548 1598231 : switch(n) {
1549 18030 : case 0: return gen_1;
1550 111510 : case 1: return utoi(v[1]);
1551 837380 : case 2: return muluu(v[1], v[2]);
1552 : }
1553 631311 : av = avma; m = n >> 1;
1554 631311 : V = cgetg(m + (odd(n)? 2: 1), t_VEC);
1555 129655861 : for (k = n; k; k--) /* start from the end: v is usually sorted */
1556 129025616 : if (v[k] & HIGHMASK) { stop = 1; break; }
1557 2072075 : while (!stop)
1558 : { /* HACK: handle V as a t_VECSMALL; gain a few iterations */
1559 69751848 : for (k = 1; k <= m; k++)
1560 : {
1561 67789834 : V[k] = uel(v,k<<1) * uel(v,(k<<1)-1);
1562 67789834 : if (V[k] & HIGHMASK) stop = 1; /* last "free" iteration */
1563 : }
1564 1962014 : if (odd(n))
1565 : {
1566 1157060 : if (n == 1) { set_avma(av); return utoi(v[1]); }
1567 635810 : V[++m] = v[n];
1568 : }
1569 1440764 : v = V; n = m; m = n >> 1;
1570 : }
1571 : /* n > 1; m > 0 */
1572 110061 : if (n == 2) { set_avma(av); return muluu(v[1], v[2]); }
1573 41418206 : for (k = 1; k <= m; k++) gel(V,k) = muluu(v[k<<1], v[(k<<1)-1]);
1574 76584 : if (odd(n)) gel(V, ++m) = utoipos(v[n]);
1575 76584 : setlg(V, m+1); /* HACK: now V is a bona fide t_VEC */
1576 76584 : return gerepileuptoint(av, gen_product(V, NULL, &_mulii));
1577 : }
1578 : GEN
1579 12595194 : vecsmall_prod(GEN v)
1580 : {
1581 12595194 : pari_sp av = avma;
1582 12595194 : long k, m, n = lg(v)-1;
1583 : GEN V;
1584 12595194 : switch (n) {
1585 0 : case 0: return gen_1;
1586 0 : case 1: return stoi(v[1]);
1587 18 : case 2: return mulss(v[1], v[2]);
1588 : }
1589 12595176 : m = n >> 1;
1590 12595176 : V = cgetg(m + (odd(n)? 2: 1), t_VEC);
1591 138477348 : for (k = 1; k <= m; k++) gel(V,k) = mulss(v[k<<1], v[(k<<1)-1]);
1592 12595176 : if (odd(n)) gel(V,k) = stoi(v[n]);
1593 12595176 : return gerepileuptoint(av, gen_product(V, NULL, &_mulii));
1594 : }
1595 :
1596 : GEN
1597 7568620 : ZV_prod(GEN v)
1598 : {
1599 7568620 : pari_sp av = avma;
1600 7568620 : long i, l = lg(v);
1601 : GEN n;
1602 7568620 : if (l == 1) return gen_1;
1603 7410083 : if (l > 7) return gerepileuptoint(av, gen_product(v, NULL, _mulii));
1604 1116103 : n = gel(v,1);
1605 1116103 : if (l == 2) return icopy(n);
1606 1804150 : for (i = 2; i < l; i++) n = mulii(n, gel(v,i));
1607 725807 : return gerepileuptoint(av, n);
1608 : }
1609 : /* assumes no overflow */
1610 : long
1611 14178 : zv_sum(GEN v)
1612 : {
1613 14178 : long n, i, l = lg(v);
1614 14178 : if (l == 1) return 0;
1615 82232 : n = v[1]; for (i = 2; i < l; i++) n += v[i];
1616 14160 : 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 66 : ZV_sum(GEN v)
1629 : {
1630 66 : pari_sp av = avma;
1631 66 : long i, l = lg(v);
1632 : GEN n;
1633 66 : if (l == 1) return gen_0;
1634 66 : n = gel(v,1);
1635 66 : if (l == 2) return icopy(n);
1636 498 : for (i = 2; i < l; i++) n = addii(n, gel(v,i));
1637 66 : 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 308947 : Zupdate_row(long k, long l, GEN q, GEN L, GEN B)
1649 : {
1650 308947 : long i, qq = itos_or_0(q);
1651 308947 : if (!qq)
1652 : {
1653 32081 : for(i=1;i<l;i++) gcoeff(L,k,i) = addii(gcoeff(L,k,i),mulii(q,gcoeff(L,l,i)));
1654 6541 : gcoeff(L,k,l) = addii(gcoeff(L,k,l), mulii(q,B));
1655 6541 : return;
1656 : }
1657 302406 : if (qq == 1) {
1658 127714 : for (i=1;i<l; i++) gcoeff(L,k,i) = addii(gcoeff(L,k,i),gcoeff(L,l,i));
1659 87729 : gcoeff(L,k,l) = addii(gcoeff(L,k,l), B);
1660 214677 : } else if (qq == -1) {
1661 129255 : for (i=1;i<l; i++) gcoeff(L,k,i) = subii(gcoeff(L,k,i),gcoeff(L,l,i));
1662 75690 : gcoeff(L,k,l) = addii(gcoeff(L,k,l), negi(B));
1663 : } else {
1664 245436 : for(i=1;i<l;i++) gcoeff(L,k,i) = addii(gcoeff(L,k,i),mulsi(qq,gcoeff(L,l,i)));
1665 138987 : gcoeff(L,k,l) = addii(gcoeff(L,k,l), mulsi(qq,B));
1666 : }
1667 : }
1668 :
1669 : /* update L[k,] */
1670 : static void
1671 889754 : ZRED(long k, long l, GEN x, GEN L, GEN B)
1672 : {
1673 889754 : GEN q = truedivii(addii(B,shifti(gcoeff(L,k,l),1)), shifti(B,1));
1674 889754 : if (!signe(q)) return;
1675 308947 : q = negi(q);
1676 308947 : Zupdate_row(k,l,q,L,B);
1677 308947 : 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 1020153 : ZincrementalGS(GEN x, GEN L, GEN B, long k)
1683 : {
1684 : long i, j;
1685 3260712 : for (j=1; j<=k; j++)
1686 : {
1687 2240559 : pari_sp av = avma;
1688 2240559 : GEN u = ZV_dotproduct(gel(x,k), gel(x,j));
1689 4849449 : for (i=1; i<j; i++)
1690 : {
1691 2608890 : u = subii(mulii(gel(B,i+1), u), mulii(gcoeff(L,k,i), gcoeff(L,j,i)));
1692 2608890 : u = diviiexact(u, gel(B,i));
1693 : }
1694 2240559 : gcoeff(L,k,j) = gerepileuptoint(av, u);
1695 : }
1696 1020153 : gel(B,k+1) = gcoeff(L,k,k); gcoeff(L,k,k) = gen_1;
1697 1020153 : }
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 96199 : ZC_reducemodmatrix_i(GEN v, GEN y)
1703 : {
1704 96199 : GEN B, L, x = shallowconcat(y, v);
1705 96199 : long k, lx = lg(x), nx = lx-1;
1706 :
1707 96199 : B = scalarcol_shallow(gen_1, lx);
1708 96199 : L = zeromatcopy(nx, nx);
1709 395734 : for (k=1; k <= nx; k++) ZincrementalGS(x, L, B, k);
1710 299535 : for (k = nx-1; k >= 1; k--) ZRED(nx,k, x,L,gel(B,k+1));
1711 96199 : return gel(x,nx);
1712 : }
1713 : GEN
1714 96199 : ZC_reducemodmatrix(GEN v, GEN y) {
1715 96199 : pari_sp av = avma;
1716 96199 : 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 194646 : ZM_reducemodmatrix_i(GEN v, GEN y)
1723 : {
1724 : GEN B, L, V;
1725 194646 : long j, k, lv = lg(v), nx = lg(y), lx = nx+1;
1726 :
1727 194646 : V = cgetg(lv, t_MAT);
1728 194646 : B = scalarcol_shallow(gen_1, lx);
1729 194646 : L = zeromatcopy(nx, nx);
1730 515706 : for (k=1; k < nx; k++) ZincrementalGS(y, L, B, k);
1731 594204 : for (j = 1; j < lg(v); j++)
1732 : {
1733 399558 : GEN x = shallowconcat(y, gel(v,j));
1734 399558 : ZincrementalGS(x, L, B, nx); /* overwrite last */
1735 1085976 : for (k = nx-1; k >= 1; k--) ZRED(nx,k, x,L,gel(B,k+1));
1736 399558 : gel(V,j) = gel(x,nx);
1737 : }
1738 194646 : return V;
1739 : }
1740 : GEN
1741 194646 : ZM_reducemodmatrix(GEN v, GEN y) {
1742 194646 : pari_sp av = avma;
1743 194646 : return gerepilecopy(av, ZM_reducemodmatrix_i(v,y));
1744 : }
1745 :
1746 : GEN
1747 83938 : ZC_reducemodlll(GEN x,GEN y)
1748 : {
1749 83938 : pari_sp av = avma;
1750 83938 : GEN z = ZC_reducemodmatrix(x, ZM_lll(y, 0.75, LLL_INPLACE));
1751 83938 : 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 : }
|