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 1945024 : check_ZV(GEN x, long l)
20 : {
21 : long i;
22 12902363 : for (i=1; i<l; i++)
23 10957395 : if (typ(gel(x,i)) != t_INT) return 0;
24 1944968 : return 1;
25 : }
26 : void
27 1498120 : RgV_check_ZV(GEN A, const char *s)
28 : {
29 1498120 : if (!RgV_is_ZV(A)) pari_err_TYPE(stack_strcat(s," [integer vector]"), A);
30 1498113 : }
31 : void
32 598128 : RgM_check_ZM(GEN A, const char *s)
33 : {
34 598128 : long n = lg(A);
35 598128 : if (n != 1)
36 : {
37 597995 : long j, m = lgcols(A);
38 2542966 : for (j=1; j<n; j++)
39 1945024 : if (!check_ZV(gel(A,j), m))
40 56 : pari_err_TYPE(stack_strcat(s," [integer matrix]"), A);
41 : }
42 598075 : }
43 :
44 : /* assume m > 1 */
45 : static long
46 105722446 : ZV_max_lg_i(GEN x, long m)
47 : {
48 105722446 : long i, l = lgefint(gel(x,1));
49 762613406 : for (i = 2; i < m; i++) l = maxss(l, lgefint(gel(x,i)));
50 105722298 : 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 27330450 : ZM_max_lg_i(GEN x, long n, long m)
62 : {
63 27330450 : long j, l = ZV_max_lg_i(gel(x,1), m);
64 105712542 : for (j = 2; j < n; j++) l = maxss(l, ZV_max_lg_i(gel(x,j), m));
65 27330638 : return l;
66 : }
67 : long
68 22614 : ZM_max_lg(GEN x)
69 : {
70 22614 : long n = lg(x), m;
71 22614 : if (n == 1) return 2;
72 22614 : 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 4260 : ZM_supnorm(GEN x)
108 : {
109 4260 : long i, j, h, lx = lg(x);
110 4260 : GEN s = gen_0;
111 4260 : if (lx == 1) return gen_1;
112 4260 : h = lgcols(x);
113 26221 : for (j=1; j<lx; j++)
114 : {
115 21961 : GEN xj = gel(x,j);
116 295088 : for (i=1; i<h; i++)
117 : {
118 273127 : GEN c = gel(xj,i);
119 273127 : if (abscmpii(c, s) > 0) s = c;
120 : }
121 : }
122 4260 : 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) = gc_INT(av,s);
144 : }
145 0 : return z;
146 : }
147 :
148 : /* x ZV, y a compatible zc. */
149 : GEN
150 2229556 : ZV_zc_mul(GEN x, GEN y)
151 : {
152 2229556 : long j, l = lg(x);
153 2229556 : pari_sp av = avma;
154 2229556 : GEN s = mulis(gel(x,1),y[1]);
155 50293138 : for (j=2; j<l; j++)
156 48063582 : if (y[j]) s = addii(s, mulis(gel(x,j),y[j]));
157 2229556 : return gc_INT(av,s);
158 : }
159 :
160 : /* x nonempty ZM, y a compatible zc (dimension > 0). */
161 : static GEN
162 19423582 : ZM_zc_mul_i(GEN x, GEN y, long c, long l)
163 : {
164 : long i, j;
165 19423582 : GEN z = cgetg(l,t_COL);
166 :
167 121999967 : for (i=1; i<l; i++)
168 : {
169 102577119 : pari_sp av = avma;
170 102577119 : GEN s = mulis(gcoeff(x,i,1),y[1]);
171 1169734721 : for (j=2; j<c; j++)
172 1067148915 : if (y[j]) s = addii(s, mulis(gcoeff(x,i,j),y[j]));
173 102585806 : gel(z,i) = gc_INT(av,s);
174 : }
175 19422848 : return z;
176 : }
177 : GEN
178 17358735 : ZM_zc_mul(GEN x, GEN y) {
179 17358735 : long l = lg(x);
180 17358735 : if (l == 1) return cgetg(1, t_COL);
181 17358735 : 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 2408 : zv_ZM_mul(GEN x, GEN y) {
187 2408 : long i,j, lx = lg(x), ly = lg(y);
188 : GEN z;
189 2408 : if (lx == 1) return zerovec(ly-1);
190 2408 : z = cgetg(ly,t_VEC);
191 6734 : for (j=1; j<ly; j++)
192 : {
193 4326 : pari_sp av = avma;
194 4326 : GEN s = mulsi(x[1], gcoeff(y,1,j));
195 10038 : for (i=2; i<lx; i++)
196 5712 : if (x[i]) s = addii(s, mulsi(x[i], gcoeff(y,i,j)));
197 4326 : gel(z,j) = gc_INT(av,s);
198 : }
199 2408 : return z;
200 : }
201 :
202 : /* x ZM, y a compatible zm (dimension > 0). */
203 : GEN
204 975600 : ZM_zm_mul(GEN x, GEN y)
205 : {
206 975600 : long j, c, l = lg(x), ly = lg(y);
207 975600 : GEN z = cgetg(ly, t_MAT);
208 975600 : if (l == 1) return z;
209 975593 : c = lgcols(x);
210 3040470 : for (j = 1; j < ly; j++) gel(z,j) = ZM_zc_mul_i(x, gel(y,j), l,c);
211 975594 : 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 456728 : 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 456728 : long min_d = minss(da, db), min_e = minss(ea, eb), i, j;
235 456728 : GEN M = cgetg(n + 1, t_MAT), C;
236 :
237 2509012 : for (j = 1; j <= min_e; j++) {
238 2052284 : gel(M, j) = C = cgetg(m + 1, t_COL);
239 32821833 : for (i = 1; i <= min_d; i++)
240 30769549 : gel(C, i) = addii(gcoeff(A, ma + i, na + j),
241 30769549 : gcoeff(B, mb + i, nb + j));
242 2114811 : for (; i <= da; i++) gel(C, i) = gcoeff(A, ma + i, na + j);
243 2052284 : for (; i <= db; i++) gel(C, i) = gcoeff(B, mb + i, nb + j);
244 2052284 : for (; i <= m; i++) gel(C, i) = gen_0;
245 : }
246 512046 : for (; j <= ea; j++) {
247 55318 : gel(M, j) = C = cgetg(m + 1, t_COL);
248 211342 : for (i = 1; i <= da; i++) gel(C, i) = gcoeff(A, ma + i, na + j);
249 55318 : for (; i <= m; i++) gel(C, i) = gen_0;
250 : }
251 456728 : 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 456728 : for (; j <= n; j++) gel(M, j) = zerocol(m);
257 456728 : 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 399637 : 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 399637 : long min_d = minss(da, db), min_e = minss(ea, eb), i, j;
268 399637 : GEN M = cgetg(n + 1, t_MAT), C;
269 :
270 2213735 : for (j = 1; j <= min_e; j++) {
271 1814098 : gel(M, j) = C = cgetg(m + 1, t_COL);
272 30208722 : for (i = 1; i <= min_d; i++)
273 28394624 : gel(C, i) = subii(gcoeff(A, ma + i, na + j),
274 28394624 : gcoeff(B, mb + i, nb + j));
275 1868209 : for (; i <= da; i++) gel(C, i) = gcoeff(A, ma + i, na + j);
276 1897554 : for (; i <= db; i++) gel(C, i) = negi(gcoeff(B, mb + i, nb + j));
277 1814098 : for (; i <= m; i++) gel(C, i) = gen_0;
278 : }
279 399637 : 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 432267 : for (; j <= eb; j++) {
285 32630 : gel(M, j) = C = cgetg(m + 1, t_COL);
286 121290 : for (i = 1; i <= db; i++) gel(C, i) = negi(gcoeff(B, mb + i, nb + j));
287 32630 : for (; i <= m; i++) gel(C, i) = gen_0;
288 : }
289 432267 : for (; j <= n; j++) gel(M, j) = zerocol(m);
290 399637 : 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 57091 : ZM_mul_sw(GEN A, GEN B, long m, long n, long p)
298 : {
299 57091 : pari_sp av = avma;
300 57091 : long m1 = (m + 1)/2, m2 = m/2,
301 57091 : n1 = (n + 1)/2, n2 = n/2,
302 57091 : 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 57091 : T2 = subtract_slices(n1, p2, B, 0, n1, p1, p2, B, n1, n2, p1, p2);
309 57091 : S1 = subtract_slices(m2, n1, A, m1, m2, 0, n1, A, 0, m2, 0, n1);
310 57091 : M2 = ZM_mul_i(S1, T2, m2 + 1, n1 + 1, p2 + 1);
311 57091 : if (gc_needed(av, 1))
312 0 : (void)gc_all(av, 2, &T2, &M2); /* destroy S1 */
313 57091 : T3 = subtract_slices(n1, p1, T2, 0, n1, 0, p2, B, 0, n1, 0, p1);
314 57091 : if (gc_needed(av, 1))
315 0 : (void)gc_all(av, 2, &M2, &T3); /* destroy T2 */
316 57091 : S2 = add_slices(m2, n1, A, m1, m2, 0, n1, A, m1, m2, n1, n2);
317 57091 : T1 = subtract_slices(n1, p1, B, 0, n1, p1, p2, B, 0, n1, 0, p2);
318 57091 : M3 = ZM_mul_i(S2, T1, m2 + 1, n1 + 1, p2 + 1);
319 57091 : if (gc_needed(av, 1))
320 0 : (void)gc_all(av, 4, &M2, &T3, &S2, &M3); /* destroy T1 */
321 57091 : S3 = subtract_slices(m1, n1, S2, 0, m2, 0, n1, A, 0, m1, 0, n1);
322 57091 : if (gc_needed(av, 1))
323 0 : (void)gc_all(av, 4, &M2, &T3, &M3, &S3); /* destroy S2 */
324 57091 : A11 = matslice(A, 1, m1, 1, n1);
325 57091 : B11 = matslice(B, 1, n1, 1, p1);
326 57091 : M1 = ZM_mul_i(A11, B11, m1 + 1, n1 + 1, p1 + 1);
327 57091 : if (gc_needed(av, 1))
328 0 : (void)gc_all(av, 5, &M2, &T3, &M3, &S3, &M1); /* destroy A11, B11 */
329 57091 : A12 = matslice(A, 1, m1, n1 + 1, n);
330 57091 : B21 = matslice(B, n1 + 1, n, 1, p1);
331 57091 : M4 = ZM_mul_i(A12, B21, m1 + 1, n2 + 1, p1 + 1);
332 57091 : if (gc_needed(av, 1))
333 0 : (void)gc_all(av, 6, &M2, &T3, &M3, &S3, &M1, &M4); /* destroy A12, B21 */
334 57091 : C11 = add_slices(m1, p1, M1, 0, m1, 0, p1, M4, 0, m1, 0, p1);
335 57091 : if (gc_needed(av, 1))
336 0 : (void)gc_all(av, 6, &M2, &T3, &M3, &S3, &M1, &C11); /* destroy M4 */
337 57091 : M5 = ZM_mul_i(S3, T3, m1 + 1, n1 + 1, p1 + 1);
338 57091 : S4 = subtract_slices(m1, n2, A, 0, m1, n1, n2, S3, 0, m1, 0, n2);
339 57091 : if (gc_needed(av, 1))
340 5 : (void)gc_all(av, 7, &M2, &T3, &M3, &M1, &C11, &M5, &S4); /* destroy S3 */
341 57091 : T4 = add_slices(n2, p1, B, n1, n2, 0, p1, T3, 0, n2, 0, p1);
342 57091 : if (gc_needed(av, 1))
343 0 : (void)gc_all(av, 7, &M2, &M3, &M1, &C11, &M5, &S4, &T4); /* destroy T3 */
344 57091 : V1 = subtract_slices(m1, p1, M1, 0, m1, 0, p1, M5, 0, m1, 0, p1);
345 57091 : if (gc_needed(av, 1))
346 1 : (void)gc_all(av, 6, &M2, &M3, &S4, &T4, &C11, &V1); /* destroy M1, M5 */
347 57091 : B22 = matslice(B, n1 + 1, n, p1 + 1, p);
348 57091 : M6 = ZM_mul_i(S4, B22, m1 + 1, n2 + 1, p2 + 1);
349 57091 : if (gc_needed(av, 1))
350 6 : (void)gc_all(av, 6, &M2, &M3, &T4, &C11, &V1, &M6); /* destroy S4, B22 */
351 57091 : A22 = matslice(A, m1 + 1, m, n1 + 1, n);
352 57091 : M7 = ZM_mul_i(A22, T4, m2 + 1, n2 + 1, p1 + 1);
353 57091 : if (gc_needed(av, 1))
354 0 : (void)gc_all(av, 6, &M2, &M3, &C11, &V1, &M6, &M7); /* destroy A22, T4 */
355 57091 : V3 = add_slices(m1, p2, V1, 0, m1, 0, p2, M3, 0, m2, 0, p2);
356 57091 : C12 = add_slices(m1, p2, V3, 0, m1, 0, p2, M6, 0, m1, 0, p2);
357 57091 : if (gc_needed(av, 1))
358 6 : (void)gc_all(av, 6, &M2, &M3, &C11, &V1, &M7, &C12); /* destroy V3, M6 */
359 57091 : V2 = add_slices(m2, p1, V1, 0, m2, 0, p1, M2, 0, m2, 0, p2);
360 57091 : if (gc_needed(av, 1))
361 0 : (void)gc_all(av, 5, &M3, &C11, &M7, &C12, &V2); /* destroy V1, M2 */
362 57091 : C21 = add_slices(m2, p1, V2, 0, m2, 0, p1, M7, 0, m2, 0, p1);
363 57091 : if (gc_needed(av, 1))
364 6 : (void)gc_all(av, 5, &M3, &C11, &C12, &V2, &C21); /* destroy M7 */
365 57091 : C22 = add_slices(m2, p2, V2, 0, m2, 0, p2, M3, 0, m2, 0, p2);
366 57091 : if (gc_needed(av, 1))
367 0 : (void)gc_all(av, 4, &C11, &C12, &C21, &C22); /* destroy V2, M3 */
368 57091 : C = shallowconcat(vconcat(C11, C21), vconcat(C12, C22));
369 57091 : return gc_GEN(av, C);
370 : }
371 :
372 : /* x[i,]*y. Assume lg(x) > 1 and 0 < i < lgcols(x) */
373 : static GEN
374 551125520 : ZMrow_ZC_mul_i(GEN x, GEN y, long i, long lx)
375 : {
376 551125520 : pari_sp av = avma;
377 551125520 : GEN c = mulii(gcoeff(x,i,1), gel(y,1)), ZERO = gen_0;
378 : long k;
379 5574409245 : for (k = 2; k < lx; k++)
380 : {
381 5023662002 : GEN t = mulii(gcoeff(x,i,k), gel(y,k));
382 5022216895 : if (t != ZERO) c = addii(c, t);
383 : }
384 550747243 : return gc_INT(av, c);
385 : }
386 : GEN
387 135604271 : ZMrow_ZC_mul(GEN x, GEN y, long i)
388 135604271 : { 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 75496325 : ZM_ZC_mul_i(GEN x, GEN y, long lx, long l)
393 : {
394 75496325 : GEN z = cgetg(l,t_COL);
395 : long i;
396 491045460 : for (i=1; i<l; i++) gel(z,i) = ZMrow_ZC_mul_i(x,y,i,lx);
397 75455798 : return z;
398 : }
399 :
400 : static GEN
401 13608395 : ZM_mul_classical(GEN x, GEN y, long l, long lx, long ly)
402 : {
403 : long j;
404 13608395 : GEN z = cgetg(ly, t_MAT);
405 63210361 : for (j = 1; j < ly; j++)
406 49604638 : gel(z, j) = ZM_ZC_mul_i(x, gel(y, j), lx, l);
407 13605723 : return z;
408 : }
409 :
410 : static GEN
411 1101 : ZM_mul_slice(GEN A, GEN B, GEN P, GEN *mod)
412 : {
413 1101 : pari_sp av = avma;
414 1101 : long i, n = lg(P)-1;
415 : GEN H, T;
416 1101 : 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 = gc_upto(av, Flm_to_ZM(Flm_mul(a, b, p)));
422 0 : *mod = utoi(p); return Hp;
423 : }
424 1101 : 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 1101 : ZM_mul_worker(GEN P, GEN A, GEN B)
436 : {
437 1101 : GEN V = cgetg(3, t_VEC);
438 1101 : 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 gc_upto(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 13665454 : sw_bound(long s)
462 13665454 : { 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 21138364 : 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 17922777 : long Flm_sw_bound = 70;
472 : #else
473 3215587 : long Flm_sw_bound = 120;
474 : #endif
475 21138364 : if (l == 1) return zeromat(0, ly-1);
476 21136702 : if (lx==2 && l==2 && ly==2)
477 350114 : { retmkmat(mkcol(mulii(gcoeff(x,1,1), gcoeff(y,1,1)))); }
478 20786588 : if (lx==3 && l==3 && ly==3) return ZM2_mul(x, y);
479 13641217 : sx = ZM_max_lg_i(x, lx, l);
480 13641915 : 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 13641937 : 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 13641098 : B = sw_bound(minss(sx, sy));
487 13641109 : if (l <= B || lx <= B || ly <= B)
488 13584034 : return ZM_mul_classical(x, y, l, lx, ly);
489 : else
490 57075 : return ZM_mul_sw(x, y, l - 1, lx - 1, ly - 1);
491 : }
492 :
493 : GEN
494 20873946 : ZM_mul(GEN x, GEN y)
495 : {
496 20873946 : long lx = lg(x), ly = lg(y);
497 20873946 : if (ly == 1) return cgetg(1,t_MAT);
498 20740221 : if (lx == 1) return zeromat(0, ly-1);
499 20738639 : 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 = gc_upto(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 gc_upto(av, H);
546 : }
547 :
548 : GEN
549 558748 : QM_mul(GEN x, GEN y)
550 : {
551 558748 : GEN dx, nx = Q_primitive_part(x, &dx);
552 558748 : GEN dy, ny = Q_primitive_part(y, &dy);
553 558748 : GEN z = ZM_mul(nx, ny);
554 558748 : if (dx || dy)
555 : {
556 474633 : GEN d = dx ? dy ? gmul(dx, dy): dx : dy;
557 474633 : if (!gequal1(d)) z = ZM_Q_mul(z, d);
558 : }
559 558748 : 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 125303 : QM_QC_mul(GEN x, GEN y)
574 : {
575 125303 : GEN dx, nx = Q_primitive_part(x, &dx);
576 125303 : GEN dy, ny = Q_primitive_part(y, &dy);
577 125303 : GEN z = ZM_ZC_mul(nx, ny);
578 125303 : if (dx || dy)
579 : {
580 125282 : GEN d = dx ? dy ? gmul(dx, dy): dx : dy;
581 125282 : if (!gequal1(d)) z = ZC_Q_mul(z, d);
582 : }
583 125303 : 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 576866 : ZM_diag_mul(GEN d, GEN m)
624 : {
625 576866 : long i, j, l = lg(d), lm = lg(m);
626 576866 : GEN y = cgetg(lm, t_MAT);
627 1669171 : for (j=1; j<lm; j++) gel(y,j) = cgetg(l, t_COL);
628 1814376 : for (i=1; i<l; i++)
629 : {
630 1237842 : GEN c = gel(d,i);
631 1237842 : if (equali1(c))
632 274682 : for (j=1; j<lm; j++) gcoeff(y,i,j) = gcoeff(m,i,j);
633 : else
634 3915023 : for (j=1; j<lm; j++) gcoeff(y,i,j) = mulii(gcoeff(m,i,j), c);
635 : }
636 576534 : return y;
637 : }
638 :
639 : /* assume lx > 1 is lg(x) = lg(y) */
640 : static GEN
641 19321950 : ZV_dotproduct_i(GEN x, GEN y, long lx)
642 : {
643 19321950 : pari_sp av = avma;
644 19321950 : GEN c = mulii(gel(x,1), gel(y,1));
645 : long i;
646 139721061 : for (i = 2; i < lx; i++)
647 : {
648 120400299 : GEN t = mulii(gel(x,i), gel(y,i));
649 120399753 : if (t != gen_0) c = addii(c, t);
650 : }
651 19320762 : return gc_INT(av, c);
652 : }
653 :
654 : /* x~ * y, assuming result is symmetric */
655 : GEN
656 532172 : ZM_transmultosym(GEN x, GEN y)
657 : {
658 532172 : long i, j, l, ly = lg(y);
659 : GEN M;
660 532172 : if (ly == 1) return cgetg(1,t_MAT);
661 : /* lg(x) = ly */
662 532172 : l = lgcols(y); /* = lgcols(x) */
663 532172 : M = cgetg(ly, t_MAT);
664 2708044 : for (i=1; i<ly; i++)
665 : {
666 2175872 : GEN xi = gel(x,i), c = cgetg(ly,t_COL);
667 2175872 : gel(M,i) = c;
668 7104247 : for (j=1; j<i; j++)
669 4928375 : gcoeff(M,i,j) = gel(c,j) = ZV_dotproduct_i(xi,gel(y,j),l);
670 2175872 : gel(c,i) = ZV_dotproduct_i(xi,gel(y,i),l);
671 : }
672 532172 : 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 25137 : ZM_sqr_i(GEN x, long l)
699 : {
700 : long s;
701 25137 : if (l == 2) { retmkmat(mkcol(sqri(gcoeff(x,1,1)))); }
702 25137 : if (l == 3) return ZM2_sqr(x);
703 24360 : s = ZM_max_lg_i(x, l, l);
704 24360 : if (l > 70) return ZM_sqr_fast(x, l, s);
705 24360 : if (l <= sw_bound(s))
706 24344 : return ZM_mul_classical(x, x, l, l, l);
707 : else
708 16 : return ZM_mul_sw(x, x, l - 1, l - 1, l - 1);
709 : }
710 :
711 : GEN
712 25137 : ZM_sqr(GEN x)
713 : {
714 25137 : long lx=lg(x);
715 25137 : if (lx==1) return cgetg(1,t_MAT);
716 25137 : return ZM_sqr_i(x, lx);
717 : }
718 : GEN
719 25983283 : ZM_ZC_mul(GEN x, GEN y)
720 : {
721 25983283 : long lx = lg(x);
722 25983283 : return lx==1? cgetg(1,t_COL): ZM_ZC_mul_i(x, y, lx, lgcols(x));
723 : }
724 :
725 : GEN
726 3487937 : ZC_Z_div(GEN x, GEN c)
727 16349740 : { pari_APPLY_type(t_COL, Qdivii(gel(x,i), c)) }
728 :
729 : GEN
730 18712 : ZM_Z_div(GEN x, GEN c)
731 205699 : { pari_APPLY_same(ZC_Z_div(gel(x, i), c)) }
732 :
733 : GEN
734 2735108 : ZC_Q_mul(GEN A, GEN z)
735 : {
736 2735108 : pari_sp av = avma;
737 2735108 : long i, l = lg(A);
738 : GEN d, n, Ad, B, u;
739 2735108 : if (typ(z)==t_INT) return ZC_Z_mul(A,z);
740 2730719 : n = gel(z, 1); d = gel(z, 2);
741 2730719 : Ad = FpC_red(A, d);
742 2730656 : u = gcdii(d, FpV_factorback(Ad, NULL, d));
743 2730645 : B = cgetg(l, t_COL);
744 2730647 : if (equali1(u))
745 : {
746 417283 : for(i=1; i<l; i++)
747 351909 : gel(B, i) = mkfrac(mulii(n, gel(A,i)), d);
748 : } else
749 : {
750 18484730 : for(i=1; i<l; i++)
751 : {
752 15819548 : GEN di = gcdii(gel(Ad, i), u), ni = mulii(n, diviiexact(gel(A,i), di));
753 15819452 : if (equalii(d, di))
754 10847510 : gel(B, i) = ni;
755 : else
756 4971922 : gel(B, i) = mkfrac(ni, diviiexact(d, di));
757 : }
758 : }
759 2730556 : return gc_GEN(av, B);
760 : }
761 :
762 : GEN
763 1121951 : ZM_Q_mul(GEN x, GEN z)
764 : {
765 1121951 : if (typ(z)==t_INT) return ZM_Z_mul(x,z);
766 3276756 : 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 231238 : ZV_ZM_mul(GEN x, GEN y)
783 : {
784 231238 : long i, lx = lg(x), ly = lg(y);
785 : GEN z;
786 231238 : if (lx == 1) return zerovec(ly-1);
787 231119 : z = cgetg(ly, t_VEC);
788 883743 : for (i = 1; i < ly; i++) gel(z,i) = ZV_dotproduct_i(x, gel(y,i), lx);
789 231119 : 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 6831844 : ZV_dotsquare(GEN x)
809 : {
810 : long i, lx;
811 : pari_sp av;
812 : GEN z;
813 6831844 : lx = lg(x);
814 6831844 : if (lx == 1) return gen_0;
815 6831844 : av = avma; z = sqri(gel(x,1));
816 26626153 : for (i=2; i<lx; i++) z = addii(z, sqri(gel(x,i)));
817 6830278 : return gc_INT(av,z);
818 : }
819 :
820 : GEN
821 16476760 : ZV_dotproduct(GEN x,GEN y)
822 : {
823 : long lx;
824 16476760 : if (x == y) return ZV_dotsquare(x);
825 11556769 : lx = lg(x);
826 11556769 : if (lx == 1) return gen_0;
827 11556769 : 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 24059 : _ZM_sqr(void *data /*ignored*/, GEN x)
835 24059 : { (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 gc_GEN(av, gen_pow_i(x, n, NULL, &_ZM_sqr, &_ZM_mul));
844 : }
845 : GEN
846 23534 : ZM_powu(GEN x, ulong n)
847 : {
848 23534 : pari_sp av = avma;
849 23534 : if (!n) return matid(lg(x)-1);
850 23534 : return gc_GEN(av, gen_powu_i(x, n, NULL, &_ZM_sqr, &_ZM_mul));
851 : }
852 : /********************************************************************/
853 : /** **/
854 : /** ADD, SUB **/
855 : /** **/
856 : /********************************************************************/
857 : static GEN
858 35034231 : ZC_add_i(GEN x, GEN y, long lx)
859 : {
860 35034231 : GEN A = cgetg(lx, t_COL);
861 : long i;
862 424196532 : for (i=1; i<lx; i++) gel(A,i) = addii(gel(x,i), gel(y,i));
863 35026280 : return A;
864 : }
865 : GEN
866 26781214 : ZC_add(GEN x, GEN y) { return ZC_add_i(x, y, lg(x)); }
867 : GEN
868 358504 : ZC_Z_add(GEN x, GEN y)
869 : {
870 358504 : long k, lx = lg(x);
871 358504 : GEN z = cgetg(lx, t_COL);
872 358500 : if (lx == 1) pari_err_TYPE2("+",x,y);
873 358500 : gel(z,1) = addii(y,gel(x,1));
874 2447101 : for (k = 2; k < lx; k++) gel(z,k) = icopy(gel(x,k));
875 358529 : return z;
876 : }
877 :
878 : static GEN
879 9397835 : ZC_sub_i(GEN x, GEN y, long lx)
880 : {
881 : long i;
882 9397835 : GEN A = cgetg(lx, t_COL);
883 65384301 : for (i=1; i<lx; i++) gel(A,i) = subii(gel(x,i), gel(y,i));
884 9396880 : return A;
885 : }
886 : GEN
887 9069467 : 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 583874 : Z_ZC_sub(GEN a, GEN x)
900 : {
901 583874 : long k, lx = lg(x);
902 583874 : GEN z = cgetg(lx, t_COL);
903 583876 : if (lx == 1) pari_err_TYPE2("-",a,x);
904 583876 : gel(z,1) = subii(a, gel(x,1));
905 1593030 : for (k = 2; k < lx; k++) gel(z,k) = negi(gel(x,k));
906 583878 : return z;
907 : }
908 :
909 : GEN
910 788688 : ZM_add(GEN x, GEN y)
911 : {
912 788688 : long lx = lg(x), l, j;
913 : GEN z;
914 788688 : if (lx == 1) return cgetg(1, t_MAT);
915 749229 : z = cgetg(lx, t_MAT); l = lgcols(x);
916 9002077 : for (j = 1; j < lx; j++) gel(z,j) = ZC_add_i(gel(x,j), gel(y,j), l);
917 749229 : return z;
918 : }
919 : GEN
920 154606 : ZM_sub(GEN x, GEN y)
921 : {
922 154606 : long lx = lg(x), l, j;
923 : GEN z;
924 154606 : if (lx == 1) return cgetg(1, t_MAT);
925 115147 : z = cgetg(lx, t_MAT); l = lgcols(x);
926 443502 : for (j = 1; j < lx; j++) gel(z,j) = ZC_sub_i(gel(x,j), gel(y,j), l);
927 115148 : return z;
928 : }
929 : /********************************************************************/
930 : /** **/
931 : /** LINEAR COMBINATION **/
932 : /** **/
933 : /********************************************************************/
934 : /* return X/c assuming division is exact */
935 : GEN
936 4905865 : ZC_Z_divexact(GEN x, GEN c)
937 45546465 : { 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 502228 : ZM_Z_divexact(GEN x, GEN c)
944 3146641 : { 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 35662770 : ZC_Z_mul(GEN x, GEN c)
952 : {
953 35662770 : if (!signe(c)) return zerocol(lg(x)-1);
954 34155172 : if (is_pm1(c)) return (signe(c) > 0)? ZC_copy(x): ZC_neg(x);
955 252397996 : pari_APPLY_type(t_COL, mulii(gel(x,i), c))
956 : }
957 :
958 : GEN
959 62856 : ZC_z_mul(GEN x, long c)
960 : {
961 62856 : if (!c) return zerocol(lg(x)-1);
962 53186 : if (c == 1) return ZC_copy(x);
963 48278 : if (c ==-1) return ZC_neg(x);
964 476662 : pari_APPLY_type(t_COL, mulsi(c, gel(x,i)))
965 : }
966 :
967 : GEN
968 28570 : zv_z_mul(GEN x, long n)
969 436677 : { 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 3284449 : ZM_Z_mul(GEN X, GEN c)
995 : {
996 3284449 : long i, j, h, l = lg(X);
997 : GEN A;
998 3284449 : if (l == 1) return cgetg(1, t_MAT);
999 3284449 : h = lgcols(X);
1000 3284434 : if (!signe(c)) return zeromat(h-1, l-1);
1001 3240121 : if (is_pm1(c)) return (signe(c) > 0)? ZM_copy(X): ZM_neg(X);
1002 2831749 : A = cgetg(l, t_MAT);
1003 13862524 : for (j = 1; j < l; j++)
1004 : {
1005 11031070 : GEN a = cgetg(h, t_COL), x = gel(X, j);
1006 147721848 : for (i = 1; i < h; i++) gel(a,i) = mulii(c, gel(x,i));
1007 11030772 : gel(A,j) = a;
1008 : }
1009 2831454 : return A;
1010 : }
1011 : void
1012 73939730 : ZC_lincomb1_inplace_i(GEN X, GEN Y, GEN v, long n)
1013 : {
1014 : long i;
1015 1591533655 : for (i = n; i; i--) gel(X,i) = addmulii_inplace(gel(X,i), gel(Y,i), v);
1016 73893893 : }
1017 : /* X <- X + v Y (elementary col operation) */
1018 : void
1019 63257397 : ZC_lincomb1_inplace(GEN X, GEN Y, GEN v)
1020 : {
1021 63257397 : if (lgefint(v) != 2) return ZC_lincomb1_inplace_i(X, Y, v, lg(X)-1);
1022 : }
1023 : void
1024 29741659 : Flc_lincomb1_inplace(GEN X, GEN Y, ulong v, ulong q)
1025 : {
1026 : long i;
1027 29741659 : if (!v) return; /* v = 0 */
1028 741825589 : 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 15942742 : ZC_lincomb1(GEN v, GEN x, GEN y)
1034 123053138 : { pari_APPLY_type(t_COL, addmulii(gel(x,i), gel(y,i), v)) }
1035 :
1036 : /* -X + vY */
1037 : static GEN
1038 752416 : ZC_lincomb_1(GEN v, GEN x, GEN y)
1039 4531631 : { 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 34611917 : ZC_lincomb(GEN u, GEN v, GEN X, GEN Y)
1044 : {
1045 : long su, sv;
1046 : GEN A;
1047 :
1048 34611917 : su = signe(u); if (!su) return ZC_Z_mul(Y, v);
1049 34610965 : sv = signe(v); if (!sv) return ZC_Z_mul(X, u);
1050 33895074 : if (is_pm1(v))
1051 : {
1052 11455532 : if (is_pm1(u))
1053 : {
1054 10117346 : if (su != sv) A = ZC_sub(X, Y);
1055 2821885 : else A = ZC_add(X, Y);
1056 10116632 : if (su < 0) ZV_togglesign(A); /* in place but was created above */
1057 : }
1058 : else
1059 : {
1060 1338161 : if (sv > 0) A = ZC_lincomb1 (u, Y, X);
1061 606990 : else A = ZC_lincomb_1(u, Y, X);
1062 : }
1063 : }
1064 22440568 : else if (is_pm1(u))
1065 : {
1066 15356538 : if (su > 0) A = ZC_lincomb1 (v, X, Y);
1067 144909 : else A = ZC_lincomb_1(v, X, Y);
1068 : }
1069 : else
1070 : { /* not cgetg_copy: x may be a t_VEC */
1071 7087787 : long i, lx = lg(X);
1072 7087787 : A = cgetg(lx,t_COL);
1073 44556685 : for (i=1; i<lx; i++) gel(A,i) = lincombii(u,v,gel(X,i),gel(Y,i));
1074 : }
1075 33886041 : return A;
1076 : }
1077 :
1078 : /********************************************************************/
1079 : /** **/
1080 : /** CONVERSIONS **/
1081 : /** **/
1082 : /********************************************************************/
1083 : GEN
1084 495162 : ZV_to_nv(GEN x)
1085 927528 : { pari_APPLY_ulong(itou(gel(x,i))) }
1086 :
1087 : GEN
1088 156685 : zm_to_ZM(GEN x)
1089 914246 : { 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 17617736 : ZC_copy(GEN x)
1119 : {
1120 17617736 : long i, lx = lg(x);
1121 17617736 : GEN y = cgetg(lx, t_COL);
1122 112651641 : for (i=1; i<lx; i++)
1123 : {
1124 95033193 : GEN c = gel(x,i);
1125 95033193 : gel(y,i) = lgefint(c) == 2? gen_0: icopy(c);
1126 : }
1127 17618448 : return y;
1128 : }
1129 :
1130 : GEN
1131 661196 : ZM_copy(GEN x)
1132 5594198 : { pari_APPLY_same(ZC_copy(gel(x,i))) }
1133 :
1134 : void
1135 390916 : ZV_neg_inplace(GEN M)
1136 : {
1137 390916 : long l = lg(M);
1138 1608343 : while (--l > 0) gel(M,l) = negi(gel(M,l));
1139 390951 : }
1140 : GEN
1141 6327895 : ZC_neg(GEN x)
1142 35992853 : { pari_APPLY_type(t_COL, negi(gel(x,i))) }
1143 :
1144 : GEN
1145 51665 : zv_neg(GEN x)
1146 662439 : { pari_APPLY_long(-x[i]) }
1147 : GEN
1148 126 : zv_neg_inplace(GEN M)
1149 : {
1150 126 : long l = lg(M);
1151 432 : 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 1730039 : ZM_neg(GEN x)
1159 5176104 : { pari_APPLY_same(ZC_neg(gel(x,i))) }
1160 : int
1161 22921521 : zv_canon_inplace(GEN x)
1162 : {
1163 22921521 : long l = lg(x), j, k;
1164 41338311 : for (j = 1; j < l && x[j] == 0; ++j);
1165 22921521 : if (j < l && x[j] < 0)
1166 : {
1167 52277022 : for (k = j; k < l; ++k) x[k] = -x[k];
1168 10255819 : return -1;
1169 : }
1170 12665702 : return 1;
1171 : }
1172 :
1173 : void
1174 5120732 : ZV_togglesign(GEN M)
1175 : {
1176 5120732 : long l = lg(M);
1177 76147866 : while (--l > 0) togglesign_safe(&gel(M,l));
1178 5120805 : }
1179 : void
1180 0 : ZM_togglesign(GEN M)
1181 : {
1182 0 : long l = lg(M);
1183 0 : while (--l > 0) ZV_togglesign(gel(M,l));
1184 0 : }
1185 :
1186 : /********************************************************************/
1187 : /** **/
1188 : /** "DIVISION" mod HNF **/
1189 : /** **/
1190 : /********************************************************************/
1191 : /* Reduce ZC x modulo ZM y in HNF */
1192 : static GEN
1193 11292236 : ZC_hnfdivrem_i(GEN x, GEN y, GEN *Q, GEN (*div)(GEN,GEN))
1194 : {
1195 11292236 : long i, l = lg(x);
1196 11292236 : pari_sp av = avma;
1197 :
1198 11292236 : if (Q) *Q = cgetg(l,t_COL);
1199 11292236 : if (l == 1) return cgetg(1,t_COL);
1200 63651415 : for (i = l-1; i>0; i--)
1201 : {
1202 52362387 : GEN q = div(gel(x,i), gcoeff(y,i,i));
1203 52361347 : if (signe(q)) x = ZC_lincomb(gen_1, negi(q), x, gel(y,i));
1204 52359179 : if (Q) gel(*Q, i) = q;
1205 : }
1206 11289028 : if (avma == av) return ZC_copy(x);
1207 8473526 : if (!Q) return gc_upto(av, x);
1208 56841 : (void)gc_all(av, 2, &x, Q); return x;
1209 : }
1210 : GEN
1211 10667161 : ZC_hnfdivrem(GEN x, GEN y, GEN *Q)
1212 10667161 : { return ZC_hnfdivrem_i(x, y, Q, diviiround); }
1213 : GEN
1214 280 : ZC_modhnf(GEN x, GEN y, GEN *Q)
1215 280 : { return ZC_hnfdivrem_i(x, y, Q, truedivii); }
1216 :
1217 : /* Return R such that x = y Q + R, y integral HNF */
1218 : static GEN
1219 489470 : ZM_hnfdivrem_i(GEN x, GEN y, GEN *Q, GEN (*div)(GEN,GEN))
1220 : {
1221 : long l, i;
1222 489470 : GEN R = cgetg_copy(x, &l);
1223 489504 : if (Q)
1224 : {
1225 126074 : GEN q = cgetg(l, t_MAT); *Q = q;
1226 186052 : for (i = 1; i < l; i++)
1227 59978 : gel(R,i) = ZC_hnfdivrem_i(gel(x,i),y,&gel(q,i),div);
1228 : }
1229 : else
1230 928424 : for (i = 1; i < l; i++)
1231 564983 : gel(R,i) = ZC_hnfdivrem_i(gel(x,i),y,NULL,div);
1232 489515 : return R;
1233 : }
1234 : GEN
1235 489456 : ZM_hnfdivrem(GEN x, GEN y, GEN *Q)
1236 489456 : { return ZM_hnfdivrem_i(x, y, Q, diviiround); }
1237 : GEN
1238 14 : ZM_modhnf(GEN x, GEN y, GEN *Q)
1239 14 : { return ZM_hnfdivrem_i(x, y, Q, truedivii); }
1240 :
1241 : static GEN
1242 7 : ZV_ZV_divrem(GEN x, GEN y, GEN *pQ)
1243 : {
1244 7 : long i, l = lg(x), tx = typ(x);
1245 : GEN Q, R;
1246 :
1247 7 : if (!pQ) return ZV_ZV_mod(x, y);
1248 0 : Q = cgetg(l,tx);
1249 0 : R = cgetg(l,tx);
1250 0 : for (i = 1; i < l; i++) gel(Q,i) = truedvmdii(gel(x,i), gel(y,i), &gel(R,i));
1251 0 : *pQ = Q; return R;
1252 : }
1253 : static GEN
1254 0 : ZM_ZV_divrem(GEN x, GEN y, GEN *Q)
1255 0 : { if (!Q) return ZM_ZV_mod(x, y);
1256 0 : pari_APPLY_same(ZV_ZV_divrem(gel(x,i), y, Q)); }
1257 :
1258 : static int
1259 42 : RgM_issquare(GEN x) { long l = lg(x); return l == 1 || lg(gel(x,1)) == l; }
1260 : static void
1261 98 : matmodhnf_check(GEN x)
1262 : {
1263 98 : switch(typ(x))
1264 : {
1265 42 : case t_VEC: case t_COL:
1266 42 : if (!RgV_is_ZV(x)) pari_err_TYPE("matmodhnf", x);
1267 42 : break;
1268 56 : case t_MAT:
1269 56 : if (!RgM_is_ZM(x)) pari_err_TYPE("matmodhnf", x);
1270 56 : break;
1271 0 : default: pari_err_TYPE("matmodhnf", x);
1272 : }
1273 98 : }
1274 : GEN
1275 49 : matmodhnf(GEN x, GEN y, GEN *Q)
1276 : {
1277 49 : long tx = typ(x), ty = typ(y), ly, lx;
1278 49 : matmodhnf_check(x); lx = lg(x);
1279 49 : matmodhnf_check(y); ly = lg(y);
1280 49 : if (ty == t_MAT && !RgM_issquare(y)) pari_err_TYPE("matmodhnf", y);
1281 49 : if (tx == t_MAT && lx == 1)
1282 : {
1283 0 : if (ly != 1) pari_err_DIM("matmodhnf");
1284 0 : if (!Q) *Q = cgetg(1, t_MAT);
1285 0 : return cgetg(1, t_MAT);
1286 : }
1287 49 : if (is_vec_t(ty))
1288 7 : return tx == t_MAT? ZM_ZV_divrem(x, y, Q): ZV_ZV_divrem(x, y, Q);
1289 : /* ty = t_MAT */
1290 42 : if (tx == t_MAT) return ZM_modhnf(x, y, Q);
1291 28 : x = ZC_modhnf(x, y, Q);
1292 28 : if (tx == t_VEC) { settyp(x, tx); if (Q) settyp(*Q, tx); }
1293 28 : return x;
1294 : }
1295 :
1296 : /********************************************************************/
1297 : /** **/
1298 : /** TESTS **/
1299 : /** **/
1300 : /********************************************************************/
1301 : int
1302 23113165 : zv_equal0(GEN V)
1303 : {
1304 23113165 : long l = lg(V);
1305 37597586 : while (--l > 0)
1306 30820883 : if (V[l]) return 0;
1307 6776703 : return 1;
1308 : }
1309 :
1310 : int
1311 14613955 : ZV_equal0(GEN V)
1312 : {
1313 14613955 : long l = lg(V);
1314 25441234 : while (--l > 0)
1315 24996900 : if (signe(gel(V,l))) return 0;
1316 444334 : return 1;
1317 : }
1318 : int
1319 16231 : ZMrow_equal0(GEN V, long i)
1320 : {
1321 16231 : long l = lg(V);
1322 25183 : while (--l > 0)
1323 21679 : if (signe(gcoeff(V,i,l))) return 0;
1324 3504 : return 1;
1325 : }
1326 :
1327 : static int
1328 6482940 : ZV_equal_lg(GEN V, GEN W, long l)
1329 : {
1330 26563308 : while (--l > 0)
1331 20579503 : if (!equalii(gel(V,l), gel(W,l))) return 0;
1332 5983805 : return 1;
1333 : }
1334 : int
1335 292665 : ZV_equal(GEN V, GEN W)
1336 : {
1337 292665 : long l = lg(V);
1338 292665 : if (lg(W) != l) return 0;
1339 292644 : return ZV_equal_lg(V, W, l);
1340 : }
1341 : int
1342 3474918 : ZM_equal(GEN A, GEN B)
1343 : {
1344 3474918 : long i, m, l = lg(A);
1345 3474918 : if (lg(B) != l) return 0;
1346 3474864 : if (l == 1) return 1;
1347 3474864 : m = lgcols(A);
1348 3474861 : if (lgcols(B) != m) return 0;
1349 9383628 : for (i = 1; i < l; i++)
1350 6190277 : if (!ZV_equal_lg(gel(A,i), gel(B,i), m)) return 0;
1351 3193351 : return 1;
1352 : }
1353 : int
1354 32292 : ZM_equal0(GEN A)
1355 : {
1356 32292 : long i, j, m, l = lg(A);
1357 32292 : if (l == 1) return 1;
1358 32292 : m = lgcols(A);
1359 87954 : for (j = 1; j < l; j++)
1360 2721073 : for (i = 1; i < m; i++)
1361 2665411 : if (signe(gcoeff(A,i,j))) return 0;
1362 17181 : return 1;
1363 : }
1364 : int
1365 30816848 : zv_equal(GEN V, GEN W)
1366 : {
1367 30816848 : long l = lg(V);
1368 30816848 : if (lg(W) != l) return 0;
1369 262053871 : while (--l > 0)
1370 232353614 : if (V[l] != W[l]) return 0;
1371 29700257 : return 1;
1372 : }
1373 :
1374 : int
1375 1638 : zvV_equal(GEN V, GEN W)
1376 : {
1377 1638 : long l = lg(V);
1378 1638 : if (lg(W) != l) return 0;
1379 80388 : while (--l > 0)
1380 79912 : if (!zv_equal(gel(V,l),gel(W,l))) return 0;
1381 476 : return 1;
1382 : }
1383 :
1384 : int
1385 819030 : ZM_ishnf(GEN x)
1386 : {
1387 819030 : long i,j, lx = lg(x);
1388 2530509 : for (i=1; i<lx; i++)
1389 : {
1390 1824095 : GEN xii = gcoeff(x,i,i);
1391 1824095 : if (signe(xii) <= 0) return 0;
1392 3582696 : for (j=1; j<i; j++)
1393 1796107 : if (signe(gcoeff(x,i,j))) return 0;
1394 3640893 : for (j=i+1; j<lx; j++)
1395 : {
1396 1929414 : GEN xij = gcoeff(x,i,j);
1397 1929414 : if (signe(xij)<0 || cmpii(xij,xii)>=0) return 0;
1398 : }
1399 : }
1400 706414 : return 1;
1401 : }
1402 : int
1403 656741 : ZM_isidentity(GEN x)
1404 : {
1405 656741 : long i,j, lx = lg(x);
1406 :
1407 656741 : if (lx == 1) return 1;
1408 656734 : if (lx != lgcols(x)) return 0;
1409 3206592 : for (j=1; j<lx; j++)
1410 : {
1411 2549936 : GEN c = gel(x,j);
1412 7926915 : for (i=1; i<j; )
1413 5377001 : if (signe(gel(c,i++))) return 0;
1414 : /* i = j */
1415 2549914 : if (!equali1(gel(c,i++))) return 0;
1416 7926930 : for ( ; i<lx; )
1417 5377037 : if (signe(gel(c,i++))) return 0;
1418 : }
1419 656656 : return 1;
1420 : }
1421 : int
1422 624359 : ZM_isdiagonal(GEN x)
1423 : {
1424 624359 : long i,j, lx = lg(x);
1425 624359 : if (lx == 1) return 1;
1426 624359 : if (lx != lgcols(x)) return 0;
1427 :
1428 1629382 : for (j=1; j<lx; j++)
1429 : {
1430 1368471 : GEN c = gel(x,j);
1431 1931009 : for (i=1; i<j; i++)
1432 925943 : if (signe(gel(c,i))) return 0;
1433 2330615 : for (i++; i<lx; i++)
1434 1325584 : if (signe(gel(c,i))) return 0;
1435 : }
1436 260911 : return 1;
1437 : }
1438 : int
1439 157638 : ZM_isscalar(GEN x, GEN s)
1440 : {
1441 157638 : long i, j, lx = lg(x);
1442 :
1443 157638 : if (lx == 1) return 1;
1444 157638 : if (!s) s = gcoeff(x,1,1);
1445 157638 : if (equali1(s)) return ZM_isidentity(x);
1446 156182 : if (lx != lgcols(x)) return 0;
1447 688224 : for (j=1; j<lx; j++)
1448 : {
1449 591676 : GEN c = gel(x,j);
1450 2490135 : for (i=1; i<j; )
1451 1955976 : if (signe(gel(c,i++))) return 0;
1452 : /* i = j */
1453 534159 : if (!equalii(gel(c,i++), s)) return 0;
1454 2518620 : for ( ; i<lx; )
1455 1986578 : if (signe(gel(c,i++))) return 0;
1456 : }
1457 96548 : return 1;
1458 : }
1459 :
1460 : long
1461 171608 : ZC_is_ei(GEN x)
1462 : {
1463 171608 : long i, j = 0, l = lg(x);
1464 1625987 : for (i = 1; i < l; i++)
1465 : {
1466 1454380 : GEN c = gel(x,i);
1467 1454380 : long s = signe(c);
1468 1454380 : if (!s) continue;
1469 171602 : if (s < 0 || !is_pm1(c) || j) return 0;
1470 171601 : j = i;
1471 : }
1472 171607 : return j;
1473 : }
1474 :
1475 : /********************************************************************/
1476 : /** **/
1477 : /** MISCELLANEOUS **/
1478 : /** **/
1479 : /********************************************************************/
1480 : /* assume lg(x) = lg(y), x,y in Z^n */
1481 : int
1482 3233786 : ZV_cmp(GEN x, GEN y)
1483 : {
1484 3233786 : long fl,i, lx = lg(x);
1485 6486616 : for (i=1; i<lx; i++)
1486 5177640 : if (( fl = cmpii(gel(x,i), gel(y,i)) )) return fl;
1487 1308976 : return 0;
1488 : }
1489 : /* assume lg(x) = lg(y), x,y in Z^n */
1490 : int
1491 19530 : ZV_abscmp(GEN x, GEN y)
1492 : {
1493 19530 : long fl,i, lx = lg(x);
1494 53046 : for (i=1; i<lx; i++)
1495 52919 : if (( fl = abscmpii(gel(x,i), gel(y,i)) )) return fl;
1496 127 : return 0;
1497 : }
1498 :
1499 : long
1500 18583744 : zv_content(GEN x)
1501 : {
1502 18583744 : long i, s, l = lg(x);
1503 18583744 : if (l == 1) return 0;
1504 18583737 : s = labs(x[1]);
1505 42319184 : for (i = 2; i < l && s != 1; i++) s = ugcd(s, labs(x[i]));
1506 18583765 : return s;
1507 : }
1508 : GEN
1509 300436 : ZV_content(GEN x)
1510 : {
1511 300436 : long i, l = lg(x);
1512 300436 : pari_sp av = avma;
1513 : GEN c;
1514 300436 : if (l == 1) return gen_0;
1515 300436 : if (l == 2) return absi(gel(x,1));
1516 207483 : c = gel(x,1);
1517 562390 : for (i = 2; i < l; i++)
1518 : {
1519 407107 : c = gcdii(c, gel(x,i));
1520 407107 : if (is_pm1(c)) { set_avma(av); return gen_1; }
1521 : }
1522 155283 : return gc_INT(av, c);
1523 : }
1524 :
1525 : GEN
1526 4133045 : ZM_det_triangular(GEN mat)
1527 : {
1528 : pari_sp av;
1529 4133045 : long i,l = lg(mat);
1530 : GEN s;
1531 :
1532 4133045 : if (l<3) return l<2? gen_1: icopy(gcoeff(mat,1,1));
1533 3731541 : av = avma; s = gcoeff(mat,1,1);
1534 10076779 : for (i=2; i<l; i++) s = mulii(s,gcoeff(mat,i,i));
1535 3731113 : return gc_INT(av,s);
1536 : }
1537 :
1538 : /* assumes no overflow */
1539 : long
1540 949442 : zv_prod(GEN v)
1541 : {
1542 949442 : long n, i, l = lg(v);
1543 949442 : if (l == 1) return 1;
1544 959729 : n = v[1]; for (i = 2; i < l; i++) n *= v[i];
1545 771061 : return n;
1546 : }
1547 :
1548 : static GEN
1549 319217430 : _mulii(void *E, GEN a, GEN b)
1550 319217430 : { (void) E; return mulii(a, b); }
1551 :
1552 : /* product of ulongs */
1553 : GEN
1554 1919666 : zv_prod_Z(GEN v)
1555 : {
1556 : pari_sp av;
1557 1919666 : long k, m, n = lg(v)-1;
1558 1919666 : int stop = 0;
1559 : GEN V;
1560 1919666 : switch(n) {
1561 21014 : case 0: return gen_1;
1562 128919 : case 1: return utoi(v[1]);
1563 1035800 : case 2: return muluu(v[1], v[2]);
1564 : }
1565 733933 : av = avma; m = n >> 1;
1566 733933 : V = cgetg(m + (odd(n)? 2: 1), t_VEC);
1567 155050078 : for (k = n; k; k--) /* start from the end: v is usually sorted */
1568 154317215 : if (v[k] & HIGHMASK) { stop = 1; break; }
1569 2412496 : while (!stop)
1570 : { /* HACK: handle V as a t_VECSMALL; gain a few iterations */
1571 83347514 : for (k = 1; k <= m; k++)
1572 : {
1573 81063337 : V[k] = uel(v,k<<1) * uel(v,(k<<1)-1);
1574 81063337 : if (V[k] & HIGHMASK) stop = 1; /* last "free" iteration */
1575 : }
1576 2284177 : if (odd(n))
1577 : {
1578 1347953 : if (n == 1) { set_avma(av); return utoi(v[1]); }
1579 742339 : V[++m] = v[n];
1580 : }
1581 1678563 : v = V; n = m; m = n >> 1;
1582 : }
1583 : /* n > 1; m > 0 */
1584 128319 : if (n == 2) { set_avma(av); return muluu(v[1], v[2]); }
1585 47423534 : for (k = 1; k <= m; k++) gel(V,k) = muluu(v[k<<1], v[(k<<1)-1]);
1586 88181 : if (odd(n)) gel(V, ++m) = utoipos(v[n]);
1587 88180 : setlg(V, m+1); /* HACK: now V is a bona fide t_VEC */
1588 88180 : return gc_INT(av, gen_product(V, NULL, &_mulii));
1589 : }
1590 : GEN
1591 14694393 : vecsmall_prod(GEN v)
1592 : {
1593 14694393 : pari_sp av = avma;
1594 14694393 : long k, m, n = lg(v)-1;
1595 : GEN V;
1596 14694393 : switch (n) {
1597 0 : case 0: return gen_1;
1598 0 : case 1: return stoi(v[1]);
1599 21 : case 2: return mulss(v[1], v[2]);
1600 : }
1601 14694372 : m = n >> 1;
1602 14694372 : V = cgetg(m + (odd(n)? 2: 1), t_VEC);
1603 161556906 : for (k = 1; k <= m; k++) gel(V,k) = mulss(v[k<<1], v[(k<<1)-1]);
1604 14694372 : if (odd(n)) gel(V,k) = stoi(v[n]);
1605 14694372 : return gc_INT(av, gen_product(V, NULL, &_mulii));
1606 : }
1607 :
1608 : GEN
1609 8850600 : ZV_prod(GEN v)
1610 : {
1611 8850600 : pari_sp av = avma;
1612 8850600 : long i, l = lg(v);
1613 : GEN n;
1614 8850600 : if (l == 1) return gen_1;
1615 8668722 : if (l > 7) return gc_INT(av, gen_product(v, NULL, _mulii));
1616 1323171 : n = gel(v,1);
1617 1323171 : if (l == 2) return icopy(n);
1618 2182604 : for (i = 2; i < l; i++) n = mulii(n, gel(v,i));
1619 866560 : return gc_INT(av, n);
1620 : }
1621 : /* assumes no overflow */
1622 : long
1623 16145 : zv_sum(GEN v)
1624 : {
1625 16145 : long n, i, l = lg(v);
1626 16145 : if (l == 1) return 0;
1627 92998 : n = v[1]; for (i = 2; i < l; i++) n += v[i];
1628 16124 : return n;
1629 : }
1630 : /* assumes no overflow and 0 <= n <= #v */
1631 : long
1632 0 : zv_sumpart(GEN v, long n)
1633 : {
1634 : long i, p;
1635 0 : if (!n) return 0;
1636 0 : p = v[1]; for (i = 2; i <= n; i++) p += v[i];
1637 0 : return p;
1638 : }
1639 : GEN
1640 77 : ZV_sum(GEN v)
1641 : {
1642 77 : pari_sp av = avma;
1643 77 : long i, l = lg(v);
1644 : GEN n;
1645 77 : if (l == 1) return gen_0;
1646 77 : n = gel(v,1);
1647 77 : if (l == 2) return icopy(n);
1648 581 : for (i = 2; i < l; i++) n = addii(n, gel(v,i));
1649 77 : return gc_INT(av, n);
1650 : }
1651 :
1652 : /********************************************************************/
1653 : /** **/
1654 : /** GRAM SCHMIDT REDUCTION (integer matrices) **/
1655 : /** **/
1656 : /********************************************************************/
1657 :
1658 : /* L[k,] += q * L[l,], l < k. Inefficient if q = 0 */
1659 : static void
1660 400336 : Zupdate_row(long k, long l, GEN q, GEN L, GEN B)
1661 : {
1662 400336 : long i, qq = itos_or_0(q);
1663 400339 : if (!qq)
1664 : {
1665 32306 : for(i=1;i<l;i++) gcoeff(L,k,i) = addii(gcoeff(L,k,i),mulii(q,gcoeff(L,l,i)));
1666 7179 : gcoeff(L,k,l) = addii(gcoeff(L,k,l), mulii(q,B));
1667 7179 : return;
1668 : }
1669 393160 : if (qq == 1) {
1670 182281 : for (i=1;i<l; i++) gcoeff(L,k,i) = addii(gcoeff(L,k,i),gcoeff(L,l,i));
1671 123579 : gcoeff(L,k,l) = addii(gcoeff(L,k,l), B);
1672 269581 : } else if (qq == -1) {
1673 174778 : for (i=1;i<l; i++) gcoeff(L,k,i) = subii(gcoeff(L,k,i),gcoeff(L,l,i));
1674 106814 : gcoeff(L,k,l) = addii(gcoeff(L,k,l), negi(B));
1675 : } else {
1676 265939 : for(i=1;i<l;i++) gcoeff(L,k,i) = addii(gcoeff(L,k,i),mulsi(qq,gcoeff(L,l,i)));
1677 162778 : gcoeff(L,k,l) = addii(gcoeff(L,k,l), mulsi(qq,B));
1678 : }
1679 : }
1680 :
1681 : /* update L[k,] */
1682 : static void
1683 1165028 : ZRED(long k, long l, GEN x, GEN L, GEN B)
1684 : {
1685 1165028 : GEN q = truedivii(addii(B,shifti(gcoeff(L,k,l),1)), shifti(B,1));
1686 1165136 : if (!signe(q)) return;
1687 400307 : q = negi(q);
1688 400327 : Zupdate_row(k,l,q,L,B);
1689 400310 : gel(x,k) = ZC_lincomb(gen_1, q, gel(x,k), gel(x,l));
1690 : }
1691 :
1692 : /* Gram-Schmidt reduction, x a ZM */
1693 : static void
1694 1271913 : ZincrementalGS(GEN x, GEN L, GEN B, long k)
1695 : {
1696 : long i, j;
1697 4064136 : for (j=1; j<=k; j++)
1698 : {
1699 2793151 : pari_sp av = avma;
1700 2793151 : GEN u = ZV_dotproduct(gel(x,k), gel(x,j));
1701 5170730 : for (i=1; i<j; i++)
1702 : {
1703 2379294 : u = subii(mulii(gel(B,i+1), u), mulii(gcoeff(L,k,i), gcoeff(L,j,i)));
1704 2378691 : u = diviiexact(u, gel(B,i));
1705 : }
1706 2791436 : gcoeff(L,k,j) = gc_INT(av, u);
1707 : }
1708 1270985 : gel(B,k+1) = gcoeff(L,k,k); gcoeff(L,k,k) = gen_1;
1709 1270985 : }
1710 :
1711 : /* Variant reducemodinvertible(ZC v, ZM y), when y singular.
1712 : * Very inefficient if y is not LLL-reduced of maximal rank */
1713 : static GEN
1714 111435 : ZC_reducemodmatrix_i(GEN v, GEN y)
1715 : {
1716 111435 : GEN B, L, x = shallowconcat(y, v);
1717 111435 : long k, lx = lg(x), nx = lx-1;
1718 :
1719 111435 : B = scalarcol_shallow(gen_1, lx);
1720 111434 : L = zeromatcopy(nx, nx);
1721 456471 : for (k=1; k <= nx; k++) ZincrementalGS(x, L, B, k);
1722 345037 : for (k = nx-1; k >= 1; k--) ZRED(nx,k, x,L,gel(B,k+1));
1723 111417 : return gel(x,nx);
1724 : }
1725 : GEN
1726 111435 : ZC_reducemodmatrix(GEN v, GEN y) {
1727 111435 : pari_sp av = avma;
1728 111435 : return gc_GEN(av, ZC_reducemodmatrix_i(v,y));
1729 : }
1730 :
1731 : /* Variant reducemodinvertible(ZM v, ZM y), when y singular.
1732 : * Very inefficient if y is not LLL-reduced of maximal rank */
1733 : static GEN
1734 238751 : ZM_reducemodmatrix_i(GEN v, GEN y)
1735 : {
1736 : GEN B, L, V;
1737 238751 : long j, k, lv = lg(v), nx = lg(y), lx = nx+1;
1738 :
1739 238751 : V = cgetg(lv, t_MAT);
1740 238767 : B = scalarcol_shallow(gen_1, lx);
1741 238769 : L = zeromatcopy(nx, nx);
1742 638346 : for (k=1; k < nx; k++) ZincrementalGS(y, L, B, k);
1743 766019 : for (j = 1; j < lg(v); j++)
1744 : {
1745 527302 : GEN x = shallowconcat(y, gel(v,j));
1746 527342 : ZincrementalGS(x, L, B, nx); /* overwrite last */
1747 1458767 : for (k = nx-1; k >= 1; k--) ZRED(nx,k, x,L,gel(B,k+1));
1748 527265 : gel(V,j) = gel(x,nx);
1749 : }
1750 238717 : return V;
1751 : }
1752 : GEN
1753 238751 : ZM_reducemodmatrix(GEN v, GEN y) {
1754 238751 : pari_sp av = avma;
1755 238751 : return gc_GEN(av, ZM_reducemodmatrix_i(v,y));
1756 : }
1757 :
1758 : GEN
1759 97602 : ZC_reducemodlll(GEN x,GEN y)
1760 : {
1761 97602 : pari_sp av = avma;
1762 97602 : GEN z = ZC_reducemodmatrix(x, ZM_lll(y, 0.75, LLL_INPLACE));
1763 97606 : return gc_GEN(av, z);
1764 : }
1765 : GEN
1766 0 : ZM_reducemodlll(GEN x,GEN y)
1767 : {
1768 0 : pari_sp av = avma;
1769 0 : GEN z = ZM_reducemodmatrix(x, ZM_lll(y, 0.75, LLL_INPLACE));
1770 0 : return gc_GEN(av, z);
1771 : }
|