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 : /*******************************************************************/
19 : /** **/
20 : /** SPECIAL POLYNOMIALS **/
21 : /** **/
22 : /*******************************************************************/
23 : /* Tchebichev polynomial: T0=1; T1=X; T(n)=2*X*T(n-1)-T(n-2)
24 : * T(n) = (n/2) sum_{k=0}^{n/2} a_k x^(n-2k)
25 : * where a_k = (-1)^k 2^(n-2k) (n-k-1)! / k!(n-2k)! is an integer
26 : * and a_0 = 2^(n-1), a_k / a_{k-1} = - (n-2k+2)(n-2k+1) / 4k(n-k) */
27 : GEN
28 2156 : polchebyshev1(long n, long v) /* Assume 4*n < LONG_MAX */
29 : {
30 : long k, l;
31 : pari_sp av;
32 : GEN q,a,r;
33 :
34 2156 : if (v<0) v = 0;
35 : /* polchebyshev(-n,1) = polchebyshev(n,1) */
36 2156 : if (n < 0) n = -n;
37 2156 : if (n==0) return pol_1(v);
38 2135 : if (n==1) return pol_x(v);
39 :
40 2093 : q = cgetg(n+3, t_POL); r = q + n+2;
41 2093 : a = int2n(n-1);
42 2093 : gel(r--,0) = a;
43 2093 : gel(r--,0) = gen_0;
44 31955 : for (k=1,l=n; l>1; k++,l-=2)
45 : {
46 29862 : av = avma;
47 29862 : a = diviuuexact(muluui(l, l-1, a), 4*k, n-k);
48 29862 : togglesign(a); a = gerepileuptoint(av, a);
49 29862 : gel(r--,0) = a;
50 29862 : gel(r--,0) = gen_0;
51 : }
52 2093 : q[1] = evalsigne(1) | evalvarn(v);
53 2093 : return q;
54 : }
55 : static void
56 70 : polchebyshev1_eval_aux(long n, GEN x, GEN *pt1, GEN *pt2)
57 : {
58 : GEN t1, t2, b;
59 70 : if (n == 1) { *pt1 = gen_1; *pt2 = x; return; }
60 56 : if (n == 0) { *pt1 = x; *pt2 = gen_1; return; }
61 56 : polchebyshev1_eval_aux((n+1) >> 1, x, &t1, &t2);
62 56 : b = gsub(gmul(gmul2n(t1,1), t2), x);
63 56 : if (odd(n)) { *pt1 = gadd(gmul2n(gsqr(t1), 1), gen_m1); *pt2 = b; }
64 42 : else { *pt1 = b; *pt2 = gadd(gmul2n(gsqr(t2), 1), gen_m1); }
65 : }
66 : static GEN
67 14 : polchebyshev1_eval(long n, GEN x)
68 : {
69 : GEN t1, t2;
70 : long i, v;
71 : pari_sp av;
72 :
73 14 : if (n < 0) n = -n;
74 14 : if (n==0) return gen_1;
75 14 : if (n==1) return gcopy(x);
76 14 : av = avma;
77 14 : v = u_lvalrem(n, 2, (ulong*)&n);
78 14 : polchebyshev1_eval_aux((n+1)>>1, x, &t1, &t2);
79 14 : if (n != 1) t2 = gsub(gmul(gmul2n(t1,1), t2), x);
80 35 : for (i = 1; i <= v; i++) t2 = gadd(gmul2n(gsqr(t2), 1), gen_m1);
81 14 : return gerepileupto(av, t2);
82 : }
83 :
84 : /* Chebychev polynomial of the second kind U(n,x): the coefficient in front of
85 : * x^(n-2*m) is (-1)^m * 2^(n-2m)*(n-m)!/m!/(n-2m)! for m=0,1,...,n/2 */
86 : GEN
87 2135 : polchebyshev2(long n, long v)
88 : {
89 : pari_sp av;
90 : GEN q, a, r;
91 : long m;
92 2135 : int neg = 0;
93 :
94 2135 : if (v<0) v = 0;
95 : /* polchebyshev(-n,2) = -polchebyshev(n-2,2) */
96 2135 : if (n < 0) {
97 1050 : if (n == -1) return zeropol(v);
98 1029 : neg = 1; n = -n-2;
99 : }
100 2114 : if (n==0) return neg ? scalar_ZX_shallow(gen_m1, v): pol_1(v);
101 :
102 2072 : q = cgetg(n+3, t_POL); r = q + n+2;
103 2072 : a = int2n(n);
104 2072 : if (neg) togglesign(a);
105 2072 : gel(r--,0) = a;
106 2072 : gel(r--,0) = gen_0;
107 30807 : for (m=1; 2*m<= n; m++)
108 : {
109 28735 : av = avma;
110 28735 : a = diviuuexact(muluui(n-2*m+2, n-2*m+1, a), 4*m, n-m+1);
111 28735 : togglesign(a); a = gerepileuptoint(av, a);
112 28735 : gel(r--,0) = a;
113 28735 : gel(r--,0) = gen_0;
114 : }
115 2072 : q[1] = evalsigne(1) | evalvarn(v);
116 2072 : return q;
117 : }
118 : static void
119 91 : polchebyshev2_eval_aux(long n, GEN x, GEN *pu1, GEN *pu2)
120 : {
121 : GEN u1, u2, u, mu1;
122 91 : if (n == 1) { *pu1 = gen_1; *pu2 = gmul2n(x,1); return; }
123 70 : if (n == 0) { *pu1 = gen_0; *pu2 = gen_1; return; }
124 70 : polchebyshev2_eval_aux(n >> 1, x, &u1, &u2);
125 70 : mu1 = gneg(u1);
126 70 : u = gmul(gadd(u2,u1), gadd(u2,mu1));
127 70 : if (odd(n)) { *pu1 = u; *pu2 = gmul(gmul2n(u2,1), gadd(gmul(x,u2), mu1)); }
128 35 : else { *pu2 = u; *pu1 = gmul(gmul2n(u1,1), gadd(u2, gmul(x,mu1))); }
129 : }
130 : static GEN
131 35 : polchebyshev2_eval(long n, GEN x)
132 : {
133 : GEN u1, u2, mu1;
134 35 : long neg = 0;
135 : pari_sp av;
136 :
137 35 : if (n < 0) {
138 14 : if (n == -1) return gen_0;
139 7 : neg = 1; n = -n-2;
140 : }
141 28 : if (n==0) return neg ? gen_m1: gen_1;
142 21 : av = avma;
143 21 : polchebyshev2_eval_aux(n>>1, x, &u1, &u2);
144 21 : mu1 = gneg(u1);
145 21 : if (odd(n)) u2 = gmul(gmul2n(u2,1), gadd(gmul(x,u2), mu1));
146 14 : else u2 = gmul(gadd(u2,u1), gadd(u2,mu1));
147 21 : if (neg) u2 = gneg(u2);
148 21 : return gerepileupto(av, u2);
149 : }
150 :
151 : GEN
152 4284 : polchebyshev(long n, long kind, long v)
153 : {
154 4284 : switch (kind)
155 : {
156 2149 : case 1: return polchebyshev1(n, v);
157 2135 : case 2: return polchebyshev2(n, v);
158 0 : default: pari_err_FLAG("polchebyshev");
159 : }
160 : return NULL; /* LCOV_EXCL_LINE */
161 : }
162 : GEN
163 4333 : polchebyshev_eval(long n, long kind, GEN x)
164 : {
165 4333 : if (!x) return polchebyshev(n, kind, 0);
166 63 : if (gequalX(x)) return polchebyshev(n, kind, varn(x));
167 49 : switch (kind)
168 : {
169 14 : case 1: return polchebyshev1_eval(n, x);
170 35 : case 2: return polchebyshev2_eval(n, x);
171 0 : default: pari_err_FLAG("polchebyshev");
172 : }
173 : return NULL; /* LCOV_EXCL_LINE */
174 : }
175 :
176 : /* Hermite polynomial H(n,x): H(n+1) = 2x H(n) - 2n H(n-1)
177 : * The coefficient in front of x^(n-2*m) is
178 : * (-1)^m * n! * 2^(n-2m)/m!/(n-2m)! for m=0,1,...,n/2.. */
179 : GEN
180 1442 : polhermite(long n, long v)
181 : {
182 : long m;
183 : pari_sp av;
184 : GEN q,a,r;
185 :
186 1442 : if (v<0) v = 0;
187 1442 : if (n==0) return pol_1(v);
188 :
189 1435 : q = cgetg(n+3, t_POL); r = q + n+2;
190 1435 : a = int2n(n);
191 1435 : gel(r--,0) = a;
192 1435 : gel(r--,0) = gen_0;
193 40327 : for (m=1; 2*m<= n; m++)
194 : {
195 38892 : av = avma;
196 38892 : a = diviuexact(muluui(n-2*m+2, n-2*m+1, a), 4*m);
197 38892 : togglesign(a);
198 38892 : gel(r--,0) = a = gerepileuptoint(av, a);
199 38892 : gel(r--,0) = gen_0;
200 : }
201 1435 : q[1] = evalsigne(1) | evalvarn(v);
202 1435 : return q;
203 : }
204 : static void
205 21 : err_hermite(long n)
206 21 : { pari_err_DOMAIN("polhermite", "degree", "<", gen_0, stoi(n)); }
207 : GEN
208 1477 : polhermite_eval0(long n, GEN x, long flag)
209 : {
210 : long i;
211 : pari_sp av, av2;
212 : GEN x2, u, v;
213 :
214 1477 : if (n < 0) err_hermite(n);
215 1470 : if (!x || gequalX(x))
216 : {
217 1442 : long v = x? varn(x): 0;
218 1442 : if (flag)
219 : {
220 14 : if (!n) err_hermite(-1);
221 7 : retmkvec2(polhermite(n-1,v),polhermite(n,v));
222 : }
223 1428 : return polhermite(n, v);
224 : }
225 28 : if (n==0)
226 : {
227 7 : if (flag) err_hermite(-1);
228 0 : return gen_1;
229 : }
230 21 : if (n==1)
231 : {
232 0 : if (flag) retmkvec2(gen_1, gmul2n(x,1));
233 0 : return gmul2n(x,1);
234 : }
235 21 : av = avma; x2 = gmul2n(x,1); v = gen_1; u = x2;
236 21 : av2= avma;
237 7070 : for (i=1; i<n; i++)
238 : { /* u = H_i(x), v = H_{i-1}(x), compute t = H_{i+1}(x) */
239 : GEN t;
240 7049 : if ((i & 0xff) == 0) gerepileall(av2,2,&u, &v);
241 7049 : t = gsub(gmul(x2, u), gmulsg(2*i,v));
242 7049 : v = u; u = t;
243 : }
244 21 : if (flag) return gerepilecopy(av, mkvec2(v, u));
245 14 : return gerepileupto(av, u);
246 : }
247 : GEN
248 0 : polhermite_eval(long n, GEN x) { return polhermite_eval0(n, x, 0); }
249 :
250 : /* Legendre polynomial
251 : * L0=1; L1=X; (n+1)*L(n+1)=(2*n+1)*X*L(n)-n*L(n-1)
252 : * L(n) = 2^-n sum_{k=0}^{n/2} a_k x^(n-2k)
253 : * where a_k = (-1)^k (2n-2k)! / k! (n-k)! (n-2k)! is an integer
254 : * and a_0 = binom(2n,n), a_k / a_{k-1} = - (n-2k+1)(n-2k+2) / 2k (2n-2k+1) */
255 : GEN
256 2163 : pollegendre(long n, long v)
257 : {
258 : long k, l;
259 : pari_sp av;
260 : GEN a, r, q;
261 :
262 2163 : if (v<0) v = 0;
263 : /* pollegendre(-n) = pollegendre(n-1) */
264 2163 : if (n < 0) n = -n-1;
265 2163 : if (n==0) return pol_1(v);
266 2121 : if (n==1) return pol_x(v);
267 :
268 2079 : av = avma;
269 2079 : q = cgetg(n+3, t_POL); r = q + n+2;
270 2079 : gel(r--,0) = a = binomialuu(n<<1,n);
271 2079 : gel(r--,0) = gen_0;
272 31423 : for (k=1,l=n; l>1; k++,l-=2)
273 : { /* l = n-2*k+2 */
274 29344 : av = avma;
275 29344 : a = diviuuexact(muluui(l, l-1, a), 2*k, n+l-1);
276 29344 : togglesign(a); a = gerepileuptoint(av, a);
277 29344 : gel(r--,0) = a;
278 29344 : gel(r--,0) = gen_0;
279 : }
280 2079 : q[1] = evalsigne(1) | evalvarn(v);
281 2079 : return gerepileupto(av, gmul2n(q,-n));
282 : }
283 : /* q such that Ln * 2^n = q(x^2) [n even] or x q(x^2) [n odd] */
284 : GEN
285 0 : pollegendre_reduced(long n, long v)
286 : {
287 : long k, l, N;
288 : pari_sp av;
289 : GEN a, r, q;
290 :
291 0 : if (v<0) v = 0;
292 : /* pollegendre(-n) = pollegendre(n-1) */
293 0 : if (n < 0) n = -n-1;
294 0 : if (n<=1) return n? scalarpol_shallow(gen_2,v): pol_1(v);
295 :
296 0 : N = n >> 1;
297 0 : q = cgetg(N+3, t_POL); r = q + N+2;
298 0 : gel(r--,0) = a = binomialuu(n<<1,n);
299 0 : for (k=1,l=n; l>1; k++,l-=2)
300 : { /* l = n-2*k+2 */
301 0 : av = avma;
302 0 : a = diviuuexact(muluui(l, l-1, a), 2*k, n+l-1);
303 0 : togglesign(a);
304 0 : gel(r--,0) = a = gerepileuptoint(av, a);
305 : }
306 0 : q[1] = evalsigne(1) | evalvarn(v);
307 0 : return q;
308 : }
309 :
310 : GEN
311 2177 : pollegendre_eval0(long n, GEN x, long flag)
312 : {
313 : pari_sp av;
314 : GEN u, v;
315 : long i;
316 :
317 2177 : if (n < 0) n = -n-1; /* L(-n) = L(n-1) */
318 : /* n >= 0 */
319 2177 : if (flag && flag != 1) pari_err_FLAG("pollegendre");
320 2177 : if (!x || gequalX(x))
321 : {
322 2156 : long v = x? varn(x): 0;
323 2156 : if (flag) retmkvec2(pollegendre(n-1,v), pollegendre(n,v));
324 2149 : return pollegendre(n, v);
325 : }
326 21 : if (n==0)
327 : {
328 0 : if (flag) retmkvec2(gen_1, gcopy(x));
329 0 : return gen_1;
330 : }
331 21 : if (n==1)
332 : {
333 0 : if (flag) retmkvec2(gcopy(x), gen_1);
334 0 : return gcopy(x);
335 : }
336 21 : av = avma; v = gen_1; u = x;
337 7070 : for (i=1; i<n; i++)
338 : { /* u = P_i(x), v = P_{i-1}(x), compute t = P_{i+1}(x) */
339 : GEN t;
340 7049 : if ((i & 0xff) == 0) gerepileall(av,2,&u, &v);
341 7049 : t = gdivgu(gsub(gmul(gmulsg(2*i+1,x), u), gmulsg(i,v)), i+1);
342 7049 : v = u; u = t;
343 : }
344 21 : if (flag) return gerepilecopy(av, mkvec2(v, u));
345 14 : return gerepileupto(av, u);
346 : }
347 : GEN
348 0 : pollegendre_eval(long n, GEN x) { return pollegendre_eval0(n, x, 0); }
349 :
350 : /* Laguerre polynomial
351 : * L0^a = 1; L1^a = -X+a+1;
352 : * (n+1)*L^a(n+1) = (-X+(2*n+a+1))*L^a(n) - (n+a)*L^a(n-1)
353 : * L^a(n) = sum_{k=0}^n (-1)^k * binom(n+a,n-k) * x^k/k! */
354 : GEN
355 2128 : pollaguerre(long n, GEN a, long v)
356 : {
357 2128 : pari_sp av = avma;
358 2128 : GEN L = cgetg(n+3, t_POL), c1 = gen_1, c2 = mpfact(n);
359 : long i;
360 :
361 2128 : L[1] = evalsigne(1) | evalvarn(v);
362 2128 : if (odd(n)) togglesign_safe(&c2);
363 117404 : for (i = n; i >= 0; i--)
364 : {
365 115276 : gel(L, i+2) = gdiv(c1, c2);
366 115276 : if (i)
367 : {
368 113148 : c2 = divis(c2,-i);
369 113148 : c1 = gdivgu(gmul(c1, gaddsg(i,a)), n+1-i);
370 : }
371 : }
372 2128 : return gerepilecopy(av, L);
373 : }
374 : static void
375 21 : err_lag(long n)
376 21 : { pari_err_DOMAIN("pollaguerre", "degree", "<", gen_0, stoi(n)); }
377 : GEN
378 2163 : pollaguerre_eval0(long n, GEN a, GEN x, long flag)
379 : {
380 2163 : pari_sp av = avma;
381 : long i;
382 : GEN v, u;
383 :
384 2163 : if (n < 0) err_lag(n);
385 2156 : if (flag && flag != 1) pari_err_FLAG("pollaguerre");
386 2156 : if (!a) a = gen_0;
387 2156 : if (!x || gequalX(x))
388 : {
389 2128 : long v = x? varn(x): 0;
390 2128 : if (flag)
391 : {
392 14 : if (!n) err_lag(-1);
393 7 : retmkvec2(pollaguerre(n-1,a,v), pollaguerre(n,a,v));
394 : }
395 2114 : return pollaguerre(n,a,v);
396 : }
397 28 : if (n==0)
398 : {
399 7 : if (flag) err_lag(-1);
400 0 : return gen_1;
401 : }
402 21 : if (n==1)
403 : {
404 0 : if (flag) retmkvec2(gsub(gaddgs(a,1),x), gen_1);
405 0 : return gsub(gaddgs(a,1),x);
406 : }
407 21 : av = avma; v = gen_1; u = gsub(gaddgs(a,1),x);
408 7070 : for (i=1; i<n; i++)
409 : { /* u = P_i(x), v = P_{i-1}(x), compute t = P_{i+1}(x) */
410 : GEN t;
411 7049 : if ((i & 0xff) == 0) gerepileall(av,2,&u, &v);
412 7049 : t = gdivgu(gsub(gmul(gsub(gaddsg(2*i+1,a),x), u), gmul(gaddsg(i,a),v)), i+1);
413 7049 : v = u; u = t;
414 : }
415 21 : if (flag) return gerepilecopy(av, mkvec2(v, u));
416 14 : return gerepileupto(av, u);
417 : }
418 : GEN
419 0 : pollaguerre_eval(long n, GEN x, GEN a) { return pollaguerre_eval0(n, x, a, 0); }
420 :
421 : /* polcyclo(p) = X^(p-1) + ... + 1 */
422 : static GEN
423 504774 : polcyclo_prime(long p, long v)
424 : {
425 504774 : GEN T = cgetg(p+2, t_POL);
426 : long i;
427 504774 : T[1] = evalsigne(1) | evalvarn(v);
428 3432807 : for (i = 2; i < p+2; i++) gel(T,i) = gen_1;
429 504774 : return T;
430 : }
431 :
432 : /* cyclotomic polynomial */
433 : GEN
434 629968 : polcyclo(long n, long v)
435 : {
436 : long s, q, i, l;
437 629968 : pari_sp av=avma;
438 : GEN T, P;
439 :
440 629968 : if (v<0) v = 0;
441 629968 : if (n < 3)
442 125194 : switch(n)
443 : {
444 33124 : case 1: return deg1pol_shallow(gen_1, gen_m1, v);
445 92070 : case 2: return deg1pol_shallow(gen_1, gen_1, v);
446 0 : default: pari_err_DOMAIN("polcyclo", "index", "<=", gen_0, stoi(n));
447 : }
448 504774 : P = gel(factoru(n), 1); l = lg(P);
449 504774 : s = P[1]; T = polcyclo_prime(s, v);
450 793699 : for (i = 2; i < l; i++)
451 : { /* Phi_{np}(X) = Phi_n(X^p) / Phi_n(X) */
452 288925 : s *= P[i];
453 288925 : T = RgX_div(RgX_inflate(T, P[i]), T);
454 : }
455 : /* s = squarefree part of n */
456 504774 : q = n / s;
457 504774 : if (q == 1) return gerepileupto(av, T);
458 243161 : return gerepilecopy(av, RgX_inflate(T,q));
459 : }
460 :
461 : /* cyclotomic polynomial */
462 : GEN
463 100239 : polcyclo_eval(long n, GEN x)
464 : {
465 100239 : pari_sp av= avma;
466 : GEN P, md, xd, yneg, ypos;
467 100239 : long vpx, l, s, i, j, q, tx, root_of_1 = 0;
468 :
469 100239 : if (!x) return polcyclo(n, 0);
470 15133 : tx = typ(x);
471 15133 : if (gequalX(x)) return polcyclo(n, varn(x));
472 14538 : if (n <= 0) pari_err_DOMAIN("polcyclo", "index", "<=", gen_0, stoi(n));
473 14538 : if (n == 1) return gsubgs(x, 1);
474 14538 : if (tx == t_INT && !signe(x)) return gen_1;
475 15714 : while ((n & 3) == 0) { n >>= 1; x = gsqr(x); } /* Phi_4n(x) = Phi_2n(x^2) */
476 : /* n not divisible by 4 */
477 14538 : if (n == 2) return gerepileupto(av, gaddgs(x,1));
478 6411 : if (!odd(n)) { n >>= 1; x = gneg(x); } /* Phi_2n(x) = Phi_n(-x) for n>1 odd */
479 : /* n odd > 2. s largest squarefree divisor of n */
480 6411 : P = gel(factoru(n), 1); s = zv_prod(P);
481 : /* replace n by largest squarefree divisor */
482 6411 : q = n/s; if (q != 1) { x = gpowgs(x, q); n = s; }
483 6411 : l = lg(P)-1;
484 : /* n squarefree odd > 2, l distinct prime divisors. Now handle x = 1 or -1 */
485 6411 : if (tx == t_INT) { /* shortcut */
486 1715 : if (is_pm1(x))
487 : {
488 56 : set_avma(av);
489 56 : if (signe(x) > 0 && l == 1) return utoipos(P[1]);
490 35 : return gen_1;
491 : }
492 : } else {
493 4696 : if (gequal1(x))
494 : { /* n is prime, return n; multiply by x to keep the type */
495 14 : if (l == 1) return gerepileupto(av, gmulgu(x,n));
496 7 : return gerepilecopy(av, x); /* else 1 */
497 : }
498 4682 : if (gequalm1(x)) return gerepileupto(av, gneg(x)); /* -1 */
499 : }
500 : /* Heuristic: evaluation will probably not improve things */
501 6334 : if (tx == t_POL || tx == t_MAT || lg(x) > n)
502 24 : return gerepileupto(av, poleval(polcyclo(n,0), x));
503 :
504 6310 : xd = cgetg((1L<<l) + 1, t_VEC); /* the x^d, where d | n */
505 6310 : md = cgetg((1L<<l) + 1, t_VECSMALL); /* the mu(d), where d | n */
506 6310 : gel(xd, 1) = x;
507 6310 : md[1] = 1;
508 : /* Use Phi_n(x) = Prod_{d|n} (x^d-1)^mu(n/d).
509 : * If x has exact order D, n = Dq, then the result is 0 if q = 1. Otherwise
510 : * the factors with x^d-1, D|d are omitted and we multiply at the end by
511 : * prod_{d | q} d^mu(q/d) = q if prime, 1 otherwise */
512 : /* We store the factors with mu(d)= 1 (resp.-1) in ypos (resp yneg).
513 : * At the end we return ypos/yneg if mu(n)=1 and yneg/ypos if mu(n)=-1 */
514 6310 : ypos = gsubgs(x,1);
515 6310 : yneg = gen_1;
516 6310 : vpx = (typ(x) == t_PADIC)? valp(x): 0;
517 13961 : for (i = 1; i <= l; i++)
518 : {
519 7651 : long ti = 1L<<(i-1), p = P[i];
520 16741 : for (j = 1; j <= ti; j++) {
521 9090 : GEN X = gel(xd,j), t;
522 9090 : if (vpx > 0)
523 : { /* ypos, X t_PADIC */
524 98 : ulong a = umuluu_or_0(p, valp(X)), b = precp(ypos) - 1;
525 98 : long e = (a && a < b) ? b - a : 0;
526 98 : if (precp(X) > e) X = cvtop(X, gel(ypos,2), e);
527 98 : if (e > 0) X = gpowgs(X, p); /* avoid valp overflow of p-adic 0*/
528 : }
529 : else
530 8992 : X = gpowgs(X, p);
531 9090 : md[ti+j] = -md[j];
532 9090 : gel(xd,ti+j) = X;
533 : /* avoid precp overflow */
534 9090 : t = (vpx > 0 && gequal0(X))? gen_m1: gsubgs(X,1);
535 9090 : if (gequal0(t))
536 : { /* x^d = 1; root_of_1 := the smallest index ti+j such that X == 1
537 : * (whose bits code d: bit i-1 is set iff P[i] | d). If no such index
538 : * exists, then root_of_1 remains 0. Do not multiply with X-1 if X = 1,
539 : * we handle these factors at the end */
540 28 : if (!root_of_1) root_of_1 = ti+j;
541 : }
542 : else
543 : {
544 9062 : if (md[ti+j] == 1) ypos = gmul(ypos, t);
545 7686 : else yneg = gmul(yneg, t);
546 : }
547 : }
548 : }
549 6310 : ypos = odd(l)? gdiv(yneg,ypos): gdiv(ypos,yneg);
550 6310 : if (root_of_1)
551 : {
552 21 : GEN X = gel(xd,(1<<l)); /* = x^n = 1 */
553 21 : long bitmask_q = (1<<l) - root_of_1;
554 : /* bitmask_q encodes q = n/d: bit (i-1) is 1 iff P[i] | q */
555 :
556 : /* x is a root of unity. If bitmask_q = 0, then x was a primitive n-th
557 : * root of 1 and the result is zero. Return X - 1 to preserve type. */
558 21 : if (!bitmask_q) return gerepileupto(av, gsubgs(X, 1));
559 : /* x is a primitive d-th root of unity, where d|n and d<n: we
560 : * must multiply ypos by if(isprime(n/d), n/d, 1) */
561 7 : ypos = gmul(ypos, X); /* multiply by X = 1 to preserve type */
562 : /* If bitmask_q = 1<<(i-1) for some i <= l, then q == P[i] and we multiply
563 : * by P[i]; otherwise q is composite and nothing more needs to be done */
564 7 : if (!(bitmask_q & (bitmask_q-1))) /* detects power of 2, since bitmask!=0 */
565 : {
566 7 : i = vals(bitmask_q)+1; /* q = P[i] */
567 7 : ypos = gmulgu(ypos, P[i]);
568 : }
569 : }
570 6296 : return gerepileupto(av, ypos);
571 : }
572 : /********************************************************************/
573 : /** **/
574 : /** HILBERT & PASCAL MATRICES **/
575 : /** **/
576 : /********************************************************************/
577 : GEN
578 133 : mathilbert(long n) /* Hilbert matrix of order n */
579 : {
580 : long i,j;
581 : GEN p;
582 :
583 133 : if (n < 0) pari_err_DOMAIN("mathilbert", "dimension", "<", gen_0, stoi(n));
584 133 : p = cgetg(n+1,t_MAT);
585 1120 : for (j=1; j<=n; j++)
586 : {
587 987 : gel(p,j) = cgetg(n+1,t_COL);
588 16583 : for (i=1+(j==1); i<=n; i++)
589 15596 : gcoeff(p,i,j) = mkfrac(gen_1, utoipos(i+j-1));
590 : }
591 133 : if (n) gcoeff(p,1,1) = gen_1;
592 133 : return p;
593 : }
594 :
595 : /* q-Pascal triangle = (choose(i,j)_q) (ordinary binomial if q = NULL) */
596 : GEN
597 4004 : matqpascal(long n, GEN q)
598 : {
599 : long i, j, I;
600 4004 : pari_sp av = avma;
601 4004 : GEN m, qpow = NULL; /* gcc -Wall */
602 :
603 4004 : if (n < -1) pari_err_DOMAIN("matpascal", "n", "<", gen_m1, stoi(n));
604 4004 : n++; m = cgetg(n+1,t_MAT);
605 41174 : for (j=1; j<=n; j++) gel(m,j) = cgetg(n+1,t_COL);
606 4004 : if (q)
607 : {
608 42 : I = (n+1)/2;
609 42 : if (I > 1) { qpow = new_chunk(I+1); gel(qpow,2)=q; }
610 84 : for (j=3; j<=I; j++) gel(qpow,j) = gmul(q, gel(qpow,j-1));
611 : }
612 41174 : for (i=1; i<=n; i++)
613 : {
614 37170 : I = (i+1)/2; gcoeff(m,i,1)= gen_1;
615 37170 : if (q)
616 : {
617 483 : for (j=2; j<=I; j++)
618 238 : gcoeff(m,i,j) = gadd(gmul(gel(qpow,j),gcoeff(m,i-1,j)),
619 238 : gcoeff(m,i-1,j-1));
620 : }
621 : else
622 : {
623 1000335 : for (j=2; j<=I; j++)
624 963410 : gcoeff(m,i,j) = addii(gcoeff(m,i-1,j), gcoeff(m,i-1,j-1));
625 : }
626 1018374 : for ( ; j<=i; j++) gcoeff(m,i,j) = gcoeff(m,i,i+1-j);
627 1982022 : for ( ; j<=n; j++) gcoeff(m,i,j) = gen_0;
628 : }
629 4004 : return gerepilecopy(av, m);
630 : }
631 :
632 : GEN
633 77 : eulerianpol(long N, long v)
634 : {
635 77 : pari_sp av = avma;
636 77 : long n, n2, k = 0;
637 : GEN A;
638 77 : if (v < 0) v = 0;
639 77 : if (N < 0) pari_err_DOMAIN("eulerianpol", "index", "<", gen_0, stoi(N));
640 70 : if (N <= 1) return pol_1(v);
641 42 : if (N == 2) return deg1pol_shallow(gen_1, gen_1, v);
642 35 : A = cgetg(N+1, t_VEC);
643 35 : gel(A,1) = gen_1; gel(A,2) = gen_1; /* A_2 = x+1 */
644 567 : for (n = 3; n <= N; n++)
645 : { /* A(n,k) = (n-k)A(n-1,k-1) + (k+1)A(n-1,k) */
646 532 : n2 = n >> 1;
647 532 : if (odd(n)) gel(A,n2+1) = mului(n+1, gel(A,n2));
648 8652 : for (k = n2-1; k; k--)
649 8120 : gel(A,k+1) = addii(mului(n-k, gel(A,k)), mului(k+1, gel(A,k+1)));
650 532 : if (gc_needed(av,1))
651 : {
652 0 : if (DEBUGMEM>1) pari_warn(warnmem,"eulerianpol, %ld/%ld",n,N);
653 0 : for (k = odd(n)? n2+1: n2; k < N; k++) gel(A,k+1) = gen_0;
654 0 : A = gerepilecopy(av, A);
655 : }
656 : }
657 35 : k = N >> 1; if (odd(N)) k++;
658 329 : for (; k < N; k++) gel(A,k+1) = gel(A, N-k);
659 35 : return gerepilecopy(av, RgV_to_RgX(A, v));
660 : }
661 :
662 : /******************************************************************/
663 : /** **/
664 : /** PRECISION CHANGES **/
665 : /** **/
666 : /******************************************************************/
667 :
668 : GEN
669 98 : gprec(GEN x, long d)
670 : {
671 98 : pari_sp av = avma;
672 98 : if (d <= 0) pari_err_DOMAIN("gprec", "precision", "<=", gen_0, stoi(d));
673 98 : return gerepilecopy(av, gprec_w(x, ndec2prec(d)));
674 : }
675 :
676 : /* not GC-safe */
677 : GEN
678 11427624 : gprec_w(GEN x, long pr)
679 : {
680 11427624 : switch(typ(x))
681 : {
682 7750655 : case t_REAL:
683 7750655 : if (signe(x)) return realprec(x) != pr? rtor(x,pr): x;
684 67632 : return real_0_bit(minss(-prec2nbits(pr), expo(x)));
685 1907811 : case t_COMPLEX:
686 1907811 : retmkcomplex(gprec_w(gel(x,1),pr), gprec_w(gel(x,2),pr));
687 3867506 : case t_POL: pari_APPLY_pol_normalized(gprec_w(gel(x,i),pr));
688 30037 : case t_SER: pari_APPLY_ser_normalized(gprec_w(gel(x,i),pr));
689 444610 : case t_POLMOD: case t_RFRAC: case t_VEC: case t_COL: case t_MAT:
690 1977933 : pari_APPLY_same(gprec_w(gel(x,i), pr));
691 : }
692 711879 : return x;
693 : }
694 : /* not GC-safe */
695 : GEN
696 6290450 : gprec_wensure(GEN x, long pr)
697 : {
698 6290450 : switch(typ(x))
699 : {
700 5411529 : case t_REAL:
701 5411529 : if (signe(x)) return realprec(x) < pr? rtor(x,pr): x;
702 17179 : return real_0_bit(minss(-prec2nbits(pr), expo(x)));
703 331544 : case t_COMPLEX:
704 331544 : retmkcomplex(gprec_wensure(gel(x,1),pr), gprec_wensure(gel(x,2),pr));
705 854 : case t_POL: pari_APPLY_pol_normalized(gprec_wensure(gel(x,i),pr));
706 867482 : case t_SER: pari_APPLY_ser_normalized(gprec_wensure(gel(x,i),pr));
707 83529 : case t_POLMOD: case t_RFRAC: case t_VEC: case t_COL: case t_MAT:
708 1300863 : pari_APPLY_same(gprec_wensure(gel(x,i), pr));
709 : }
710 414064 : return x;
711 : }
712 :
713 : /* not GC-safe; truncate mantissa to precision 'pr' but never increase it */
714 : GEN
715 3683927 : gprec_wtrunc(GEN x, long pr)
716 : {
717 3683927 : switch(typ(x))
718 : {
719 3205466 : case t_REAL:
720 3205466 : return (signe(x) && realprec(x) > pr)? rtor(x,pr): x;
721 311020 : case t_COMPLEX:
722 311020 : retmkcomplex(gprec_wtrunc(gel(x,1),pr), gprec_wtrunc(gel(x,2),pr));
723 18844 : case t_POL: pari_APPLY_pol_normalized(gprec_wtrunc(gel(x,i),pr));
724 8554 : case t_SER: pari_APPLY_ser_normalized(gprec_wtrunc(gel(x,i),pr));
725 88834 : case t_POLMOD: case t_RFRAC: case t_VEC: case t_COL: case t_MAT:
726 403342 : pari_APPLY_same(gprec_wtrunc(gel(x,i), pr));
727 : }
728 74008 : return x;
729 : }
730 :
731 : /********************************************************************/
732 : /** **/
733 : /** SERIES TRANSFORMS **/
734 : /** **/
735 : /********************************************************************/
736 : /** LAPLACE TRANSFORM (OF A SERIES) **/
737 : /********************************************************************/
738 : static GEN
739 14 : serlaplace(GEN x)
740 : {
741 14 : long i, l = lg(x), e = valser(x);
742 14 : GEN t, y = cgetg(l,t_SER);
743 14 : if (e < 0) pari_err_DOMAIN("laplace","valuation","<",gen_0,stoi(e));
744 14 : t = mpfact(e); y[1] = x[1];
745 154 : for (i=2; i<l; i++)
746 : {
747 140 : gel(y,i) = gmul(t, gel(x,i));
748 140 : e++; t = mului(e,t);
749 : }
750 14 : return y;
751 : }
752 : static GEN
753 14 : pollaplace(GEN x)
754 : {
755 14 : long i, e = 0, l = lg(x);
756 14 : GEN t = gen_1, y = cgetg(l,t_POL);
757 14 : y[1] = x[1];
758 63 : for (i=2; i<l; i++)
759 : {
760 49 : gel(y,i) = gmul(t, gel(x,i));
761 49 : e++; t = mului(e,t);
762 : }
763 14 : return y;
764 : }
765 : GEN
766 35 : laplace(GEN x)
767 : {
768 35 : pari_sp av = avma;
769 35 : switch(typ(x))
770 : {
771 14 : case t_POL: x = pollaplace(x); break;
772 14 : case t_SER: x = serlaplace(x); break;
773 7 : default: if (is_scalar_t(typ(x))) return gcopy(x);
774 0 : pari_err_TYPE("laplace",x);
775 : }
776 28 : return gerepilecopy(av, x);
777 : }
778 :
779 : /********************************************************************/
780 : /** CONVOLUTION PRODUCT (OF TWO SERIES) **/
781 : /********************************************************************/
782 : GEN
783 14 : convol(GEN x, GEN y)
784 : {
785 14 : long j, lx, ly, ex, ey, vx = varn(x);
786 : GEN z;
787 :
788 14 : if (typ(x) != t_SER) pari_err_TYPE("convol",x);
789 14 : if (typ(y) != t_SER) pari_err_TYPE("convol",y);
790 14 : if (varn(y) != vx) pari_err_VAR("convol", x,y);
791 14 : ex = valser(x);
792 14 : ey = valser(y);
793 14 : if (ser_isexactzero(x))
794 : {
795 7 : z = scalarser(gadd(Rg_get_0(x), Rg_get_0(y)), varn(x), 1);
796 7 : setvalser(z, maxss(ex,ey)); return z;
797 : }
798 7 : lx = lg(x) + ex; x -= ex;
799 7 : ly = lg(y) + ey; y -= ey;
800 : /* inputs shifted: x[i] and y[i] now correspond to monomials of same degree */
801 7 : if (ly < lx) lx = ly; /* min length */
802 7 : if (ex < ey) ex = ey; /* max valuation */
803 7 : if (lx - ex < 3) return zeroser(vx, lx-2);
804 :
805 7 : z = cgetg(lx - ex, t_SER);
806 7 : z[1] = evalvalser(ex) | evalvarn(vx);
807 119 : for (j = ex+2; j<lx; j++) gel(z,j-ex) = gmul(gel(x,j),gel(y,j));
808 7 : return normalizeser(z);
809 : }
810 :
811 : /***********************************************************************/
812 : /* OPERATIONS ON DIRICHLET SERIES: *, / */
813 : /* (+, -, scalar multiplication are done on the corresponding vectors) */
814 : /***********************************************************************/
815 : static long
816 869316 : dirval(GEN x)
817 : {
818 869316 : long i = 1, lx = lg(x);
819 869337 : while (i < lx && gequal0(gel(x,i))) i++;
820 869316 : return i;
821 : }
822 :
823 : GEN
824 336 : dirmul(GEN x, GEN y)
825 : {
826 336 : pari_sp av = avma, av2;
827 : long nx, ny, nz, dx, dy, i, j, k;
828 : GEN z;
829 :
830 336 : if (typ(x)!=t_VEC) pari_err_TYPE("dirmul",x);
831 336 : if (typ(y)!=t_VEC) pari_err_TYPE("dirmul",y);
832 336 : dx = dirval(x); nx = lg(x)-1;
833 336 : dy = dirval(y); ny = lg(y)-1;
834 336 : if (ny-dy < nx-dx) { swap(x,y); lswap(nx,ny); lswap(dx,dy); }
835 336 : nz = minss(nx*dy,ny*dx);
836 336 : y = RgV_kill0(y);
837 336 : av2 = avma;
838 336 : z = zerovec(nz);
839 39095 : for (j=dx; j<=nx; j++)
840 : {
841 38759 : GEN c = gel(x,j);
842 38759 : if (gequal0(c)) continue;
843 17031 : if (gequal1(c))
844 : {
845 94199 : for (k=dy,i=j*dy; i<=nz; i+=j,k++)
846 88550 : if (gel(y,k)) gel(z,i) = gadd(gel(z,i),gel(y,k));
847 : }
848 11382 : else if (gequalm1(c))
849 : {
850 5649 : for (k=dy,i=j*dy; i<=nz; i+=j,k++)
851 4298 : if (gel(y,k)) gel(z,i) = gsub(gel(z,i),gel(y,k));
852 : }
853 : else
854 : {
855 46508 : for (k=dy,i=j*dy; i<=nz; i+=j,k++)
856 36477 : if (gel(y,k)) gel(z,i) = gadd(gel(z,i),gmul(c,gel(y,k)));
857 : }
858 17031 : if (gc_needed(av2,3))
859 : {
860 0 : if (DEBUGMEM>1) pari_warn(warnmem,"dirmul, %ld/%ld",j,nx);
861 0 : z = gerepilecopy(av2,z);
862 : }
863 : }
864 336 : return gerepilecopy(av,z);
865 : }
866 :
867 : GEN
868 434322 : dirdiv(GEN x, GEN y)
869 : {
870 434322 : pari_sp av = avma, av2;
871 : long nx,ny,nz, dx,dy, i,j,k;
872 : GEN p1;
873 :
874 434322 : if (typ(x)!=t_VEC) pari_err_TYPE("dirdiv",x);
875 434322 : if (typ(y)!=t_VEC) pari_err_TYPE("dirdiv",y);
876 434322 : dx = dirval(x); nx = lg(x)-1;
877 434322 : dy = dirval(y); ny = lg(y)-1;
878 434322 : if (dy != 1 || !ny) pari_err_INV("dirdiv",y);
879 434322 : nz = minss(nx,ny*dx);
880 434322 : p1 = gel(y,1);
881 434322 : if (gequal1(p1)) p1 = NULL; else y = gdiv(y,p1);
882 434322 : y = RgV_kill0(y);
883 434322 : av2 = avma;
884 434322 : x = p1 ? gdiv(x,p1): leafcopy(x);
885 434329 : for (j=1; j<dx; j++) gel(x,j) = gen_0;
886 434322 : setlg(x,nz+1);
887 109807992 : for (j=dx; j<=nz; j++)
888 : {
889 109373670 : GEN c = gel(x,j);
890 109373670 : if (gequal0(c)) continue;
891 75821501 : if (gequal1(c))
892 : {
893 133758387 : for (i=j+j,k=2; i<=nz; i+=j,k++)
894 131988864 : if (gel(y,k)) gel(x,i) = gsub(gel(x,i),gel(y,k));
895 : }
896 74051978 : else if (gequalm1(c))
897 : {
898 28856261 : for (i=j+j,k=2; i<=nz; i+=j,k++)
899 27302821 : if (gel(y,k)) gel(x,i) = gadd(gel(x,i),gel(y,k));
900 : }
901 : else
902 : {
903 331182936 : for (i=j+j,k=2; i<=nz; i+=j,k++)
904 258684398 : if (gel(y,k)) gel(x,i) = gsub(gel(x,i),gmul(c,gel(y,k)));
905 : }
906 75821501 : if (gc_needed(av2,3))
907 : {
908 0 : if (DEBUGMEM>1) pari_warn(warnmem,"dirdiv, %ld/%ld",j,nz);
909 0 : x = gerepilecopy(av2,x);
910 : }
911 : }
912 434322 : return gerepilecopy(av,x);
913 : }
914 :
915 : /*******************************************************************/
916 : /** **/
917 : /** COMBINATORICS **/
918 : /** **/
919 : /*******************************************************************/
920 : /** BINOMIAL COEFFICIENTS **/
921 : /*******************************************************************/
922 : /* Lucas's formula for v_p(\binom{n}{k}), used in the tough case p <= sqrt(n) */
923 : static long
924 3206 : binomial_lval(ulong n, ulong k, ulong p)
925 : {
926 3206 : ulong r = 0, e = 0;
927 : do
928 : {
929 10290 : ulong a = n % p, b = k % p + r;
930 10290 : n /= p; k /= p;
931 10290 : if (a < b) { e++; r = 1; } else r = 0;
932 10290 : } while (n);
933 3206 : return e;
934 : }
935 : GEN
936 81267 : binomialuu(ulong n, ulong k)
937 : {
938 81267 : pari_sp av = avma;
939 : ulong p, nk, sn;
940 : long c, l;
941 : forprime_t T;
942 : GEN v, z;
943 81267 : if (k > n) return gen_0;
944 81260 : nk = n-k; if (k > nk) lswap(nk, k);
945 81260 : if (!k) return gen_1;
946 79636 : if (k == 1) return utoipos(n);
947 74288 : if (k == 2) return muluu(odd(n)? n: n-1, n>>1);
948 54254 : if (k < 1000 || ((double)k/ n) * log((double)n) < 0.5)
949 : { /* k "small" */
950 54240 : z = diviiexact(mulu_interval(n-k+1, n), mulu_interval(2UL, k));
951 54241 : return gerepileuptoint(av, z);
952 : }
953 14 : sn = usqrt(n);
954 : /* use Lucas's formula, k <= n/2 */
955 14 : l = minuu(1UL << 20, n); v = cgetg(l+1, t_VECSMALL); c = 1;
956 14 : u_forprime_init(&T, nk+1, n);
957 1553958 : while ((p = u_forprime_next(&T))) /* all primes n-k < p <= n occur, v_p = 1 */
958 : {
959 1553944 : if (c == l) { ulong L = l << 1; v = vecsmall_lengthen(v, L); l = L; }
960 1553944 : v[c++] = p;
961 : }
962 14 : u_forprime_init(&T, sn+1, n >> 1);
963 2437785 : while ((p = u_forprime_next(&T))) /* p^2 > n, v_p <= 1 */
964 2437771 : if (n % p < k % p)
965 : {
966 1428679 : if (c == l) { ulong L = l << 1; v = vecsmall_lengthen(v, L); l = L; }
967 1428679 : v[c++] = p;
968 : }
969 14 : setlg(v, c); z = zv_prod_Z(v);
970 14 : u_forprime_init(&T, 3, sn);
971 14 : l = minuu(1UL << 20, sn); v = cgetg(l + 1, t_VEC); c = 1;
972 3220 : while ((p = u_forprime_next(&T))) /* p <= sqrt(n) */
973 : {
974 3206 : ulong e = binomial_lval(n, k, p);
975 3206 : if (e)
976 : {
977 2541 : if (c == l) { ulong L = l << 1; v = vec_lengthen(v, L); l = L; }
978 2541 : gel(v, c++) = powuu(p, e);
979 : }
980 : }
981 14 : setlg(v, c); z = mulii(z, ZV_prod(v));
982 : { /* p = 2 */
983 14 : ulong e = hammingl(k);
984 14 : e += (k == nk)? e: hammingl(nk);
985 14 : e -= hammingl(n); if (e) z = shifti(z, e);
986 : }
987 14 : return gerepileuptoint(av, z);
988 : }
989 :
990 : GEN
991 102722 : binomial(GEN n, long k)
992 : {
993 102722 : long i, prec, tn = typ(n);
994 : pari_sp av;
995 : GEN y;
996 :
997 102722 : av = avma;
998 102722 : if (tn == t_INT)
999 : {
1000 : long sn;
1001 : GEN z;
1002 102547 : if (k == 0) return gen_1;
1003 69157 : sn = signe(n);
1004 69157 : if (sn == 0) return gen_0; /* k != 0 */
1005 69157 : if (sn > 0)
1006 : { /* n > 0 */
1007 68709 : if (k < 0) return gen_0;
1008 68709 : if (k == 1) return icopy(n);
1009 42697 : z = subiu(n, k);
1010 42698 : if (cmpiu(z, k) < 0)
1011 : {
1012 1456 : switch(signe(z))
1013 : {
1014 7 : case -1: return gc_const(av, gen_0);
1015 63 : case 0: return gc_const(av, gen_1);
1016 : }
1017 1386 : k = z[2];
1018 1386 : if (k == 1) { set_avma(av); return icopy(n); }
1019 : }
1020 42152 : set_avma(av);
1021 42152 : if (lgefint(n) == 3) return binomialuu(n[2],(ulong)k);
1022 : }
1023 : else
1024 : { /* n < 0, k != 0; use Kronenburg's definition */
1025 448 : if (k > 0)
1026 427 : z = binomial(subsi(k - 1, n), k);
1027 : else
1028 : {
1029 21 : z = subis(n, k); if (signe(z) < 0) return gen_0;
1030 14 : n = stoi(-k-1); k = itos(z);
1031 14 : z = binomial(n, k);
1032 : }
1033 441 : if (odd(k)) togglesign_safe(&z);
1034 441 : return gerepileuptoint(av, z);
1035 : }
1036 : /* n >= 0 and huge, k != 0 */
1037 8 : if (k < 0) return gen_0;
1038 8 : if (k == 1) return icopy(n);
1039 : /* k > 1 */
1040 8 : y = cgetg(k+1,t_VEC); gel(y,1) = n;
1041 18 : for (i = 2; i <= k; i++) gel(y,i) = subiu(n,i-1);
1042 8 : y = diviiexact(ZV_prod(y), mpfact(k));
1043 8 : return gerepileuptoint(av, y);
1044 : }
1045 175 : if (is_noncalc_t(tn)) pari_err_TYPE("binomial",n);
1046 175 : if (k <= 1)
1047 : {
1048 14 : if (k < 0) return Rg_get_0(n);
1049 7 : if (k == 0) return Rg_get_1(n);
1050 0 : return gcopy(n);
1051 : }
1052 161 : prec = precision(n);
1053 161 : if (prec && k > 200 + 0.8*prec2nbits(prec)) {
1054 7 : GEN A = mpfactr(k, prec), B = ggamma(gsubgs(n,k-1), prec);
1055 7 : return gerepileupto(av, gdiv(ggamma(gaddgs(n,1), prec), gmul(A,B)));
1056 : }
1057 :
1058 154 : y = cgetg(k+1,t_VEC);
1059 12236 : for (i=1; i<=k; i++) gel(y,i) = gsubgs(n,i-1);
1060 154 : return gerepileupto(av, gdiv(RgV_prod(y), mpfact(k)));
1061 : }
1062 :
1063 : GEN
1064 1851 : binomial0(GEN x, GEN k)
1065 : {
1066 1851 : if (!k)
1067 : {
1068 21 : if (typ(x) != t_INT || signe(x) < 0) pari_err_TYPE("binomial", x);
1069 7 : return vecbinomial(itos(x));
1070 : }
1071 1830 : if (typ(k) != t_INT) pari_err_TYPE("binomial", k);
1072 1823 : return binomial(x, itos(k));
1073 : }
1074 :
1075 : /* Assume n >= 0, return bin, bin[k+1] = binomial(n, k) */
1076 : GEN
1077 185725 : vecbinomial(long n)
1078 : {
1079 : long d, k;
1080 : GEN C;
1081 185725 : if (!n) return mkvec(gen_1);
1082 185354 : C = cgetg(n+2, t_VEC) + 1; /* C[k] = binomial(n, k) */
1083 185354 : gel(C,0) = gen_1;
1084 185354 : gel(C,1) = utoipos(n); d = (n + 1) >> 1;
1085 747923 : for (k=2; k <= d; k++)
1086 : {
1087 562569 : pari_sp av = avma;
1088 562569 : gel(C,k) = gerepileuptoint(av, diviuexact(mului(n-k+1, gel(C,k-1)), k));
1089 : }
1090 823601 : for ( ; k <= n; k++) gel(C,k) = gel(C,n-k);
1091 185354 : return C - 1;
1092 : }
1093 :
1094 : /********************************************************************/
1095 : /** STIRLING NUMBERS **/
1096 : /********************************************************************/
1097 : /* Stirling number of the 2nd kind. The number of ways of partitioning
1098 : a set of n elements into m nonempty subsets. */
1099 : GEN
1100 1694 : stirling2(ulong n, ulong m)
1101 : {
1102 1694 : pari_sp av = avma;
1103 : GEN s, bmk;
1104 : ulong k;
1105 1694 : if (n==0) return (m == 0)? gen_1: gen_0;
1106 1694 : if (m > n || m == 0) return gen_0;
1107 1694 : if (m==n) return gen_1;
1108 : /* k = 0 */
1109 1694 : bmk = gen_1; s = powuu(m, n);
1110 20314 : for (k = 1; k <= ((m-1)>>1); ++k)
1111 : { /* bmk = binomial(m, k) */
1112 : GEN c, kn, mkn;
1113 18620 : bmk = diviuexact(mului(m-k+1, bmk), k);
1114 18620 : kn = powuu(k, n); mkn = powuu(m-k, n);
1115 18620 : c = odd(m)? subii(mkn,kn): addii(mkn,kn);
1116 18620 : c = mulii(bmk, c);
1117 18620 : s = odd(k)? subii(s, c): addii(s, c);
1118 18620 : if (gc_needed(av,2))
1119 : {
1120 0 : if(DEBUGMEM>1) pari_warn(warnmem,"stirling2");
1121 0 : gerepileall(av, 2, &s, &bmk);
1122 : }
1123 : }
1124 : /* k = m/2 */
1125 1694 : if (!odd(m))
1126 : {
1127 : GEN c;
1128 805 : bmk = diviuexact(mului(k+1, bmk), k);
1129 805 : c = mulii(bmk, powuu(k,n));
1130 805 : s = odd(k)? subii(s, c): addii(s, c);
1131 : }
1132 1694 : return gerepileuptoint(av, diviiexact(s, mpfact(m)));
1133 : }
1134 :
1135 : /* Stirling number of the first kind. Up to the sign, the number of
1136 : permutations of n symbols which have exactly m cycles. */
1137 : GEN
1138 154 : stirling1(ulong n, ulong m)
1139 : {
1140 154 : pari_sp ltop=avma;
1141 : ulong k;
1142 : GEN s, t;
1143 154 : if (n < m) return gen_0;
1144 154 : else if (n==m) return gen_1;
1145 : /* t = binomial(n-1+k, m-1) * binomial(2n-m, n-m-k) */
1146 : /* k = n-m > 0 */
1147 154 : t = binomialuu(2*n-m-1, m-1);
1148 154 : s = mulii(t, stirling2(2*(n-m), n-m));
1149 154 : if (odd(n-m)) togglesign(s);
1150 1547 : for (k = n-m-1; k > 0; --k)
1151 : {
1152 : GEN c;
1153 1393 : t = diviuuexact(muluui(n-m+k+1, n+k+1, t), n+k, n-m-k);
1154 1393 : c = mulii(t, stirling2(n-m+k, k));
1155 1393 : s = odd(k)? subii(s, c): addii(s, c);
1156 1393 : if ((k & 0x1f) == 0) {
1157 21 : t = gerepileuptoint(ltop, t);
1158 21 : s = gerepileuptoint(avma, s);
1159 : }
1160 : }
1161 154 : return gerepileuptoint(ltop, s);
1162 : }
1163 :
1164 : GEN
1165 301 : stirling(long n, long m, long flag)
1166 : {
1167 301 : if (n < 0) pari_err_DOMAIN("stirling", "n", "<", gen_0, stoi(n));
1168 301 : if (m < 0) pari_err_DOMAIN("stirling", "m", "<", gen_0, stoi(m));
1169 301 : switch (flag)
1170 : {
1171 154 : case 1: return stirling1((ulong)n,(ulong)m);
1172 147 : case 2: return stirling2((ulong)n,(ulong)m);
1173 0 : default: pari_err_FLAG("stirling");
1174 : }
1175 : return NULL; /*LCOV_EXCL_LINE*/
1176 : }
1177 :
1178 : /*******************************************************************/
1179 : /** **/
1180 : /** RECIPROCAL POLYNOMIAL **/
1181 : /** **/
1182 : /*******************************************************************/
1183 : /* return coefficients s.t x = x_0 X^n + ... + x_n */
1184 : GEN
1185 161 : polrecip(GEN x)
1186 : {
1187 161 : long tx = typ(x);
1188 161 : if (is_scalar_t(tx)) return gcopy(x);
1189 154 : if (tx != t_POL) pari_err_TYPE("polrecip",x);
1190 154 : return RgX_recip(x);
1191 : }
1192 :
1193 : /********************************************************************/
1194 : /** **/
1195 : /** POLYNOMIAL INTERPOLATION **/
1196 : /** **/
1197 : /********************************************************************/
1198 : /* given complex roots L[i], i <= n of some monic T in C[X], return
1199 : * the T'(L[i]), computed stably via products of differences */
1200 : GEN
1201 84300 : vandermondeinverseinit(GEN L)
1202 : {
1203 84300 : long i, j, l = lg(L);
1204 84300 : GEN V = cgetg(l, t_VEC);
1205 472545 : for (i = 1; i < l; i++)
1206 : {
1207 388245 : pari_sp av = avma;
1208 388245 : GEN W = cgetg(l-1,t_VEC);
1209 388247 : long k = 1;
1210 4236912 : for (j = 1; j < l; j++)
1211 3848750 : if (i != j) gel(W, k++) = gsub(gel(L,i), gel(L,j));
1212 388162 : gel(V,i) = gerepileupto(av, RgV_prod(W));
1213 : }
1214 84300 : return V;
1215 : }
1216 :
1217 : /* Compute the inverse of the van der Monde matrix of T multiplied by den */
1218 : GEN
1219 53408 : vandermondeinverse(GEN L, GEN T, GEN den, GEN V)
1220 : {
1221 53408 : pari_sp av = avma;
1222 53408 : long i, n = lg(L)-1;
1223 53408 : GEN M = cgetg(n+1, t_MAT);
1224 :
1225 53408 : if (!V) V = vandermondeinverseinit(L);
1226 53408 : if (den && equali1(den)) den = NULL;
1227 291596 : for (i = 1; i <= n; i++)
1228 : {
1229 476367 : GEN d = gel(V,i), P = RgX_Rg_mul(RgX_div_by_X_x(T, gel(L,i), NULL),
1230 238187 : den? gdiv(den,d): ginv(d));
1231 238184 : gel(M,i) = RgX_to_RgC(P, n);
1232 : }
1233 53409 : return gerepilecopy(av, M);
1234 : }
1235 :
1236 : static GEN
1237 224 : RgV_polint_fast(GEN X, GEN Y, long v)
1238 : {
1239 : GEN p, pol;
1240 : long t, pa;
1241 224 : if (X) t = RgV_type2(X,Y, &p, &pol, &pa);
1242 21 : else t = Rg_type(Y, &p, &pol, &pa);
1243 224 : if (t != t_INTMOD) return NULL;
1244 7 : Y = RgC_to_FpC(Y, p);
1245 7 : X = X? RgC_to_FpC(X, p): identity_ZV(lg(Y)-1);
1246 7 : return FpX_to_mod(FpV_polint(X, Y, p, v), p);
1247 : }
1248 : /* allow X = NULL for [1,...,n] */
1249 : GEN
1250 224 : RgV_polint(GEN X, GEN Y, long v)
1251 : {
1252 224 : pari_sp av0 = avma, av;
1253 224 : GEN Q, L, P = NULL;
1254 224 : long i, l = lg(Y);
1255 224 : if ((Q = RgV_polint_fast(X,Y,v))) return Q;
1256 217 : if (!X) X = identity_ZV(l-1);
1257 217 : L = vandermondeinverseinit(X);
1258 217 : Q = roots_to_pol(X, v); av = avma;
1259 553 : for (i=1; i<l; i++)
1260 : {
1261 : GEN T, dP;
1262 336 : if (gequal0(gel(Y,i))) continue;
1263 238 : T = RgX_div_by_X_x(Q, gel(X,i), NULL);
1264 238 : dP = RgX_Rg_mul(T, gdiv(gel(Y,i), gel(L,i)));
1265 238 : P = P? RgX_add(P, dP): dP;
1266 238 : if (gc_needed(av,2))
1267 : {
1268 0 : if (DEBUGMEM>1) pari_warn(warnmem,"RgV_polint i = %ld/%ld", i, l-1);
1269 0 : P = gerepileupto(av, P);
1270 : }
1271 : }
1272 217 : if (!P) { set_avma(av); return zeropol(v); }
1273 147 : return gerepileupto(av0, P);
1274 : }
1275 : static int
1276 17357 : inC(GEN x)
1277 : {
1278 17357 : switch(typ(x)) {
1279 1365 : case t_INT: case t_REAL: case t_FRAC: case t_COMPLEX: case t_QUAD: return 1;
1280 15992 : default: return 0;
1281 : }
1282 : }
1283 : static long
1284 16188 : check_dy(GEN X, GEN x, long n)
1285 : {
1286 16188 : GEN D = NULL;
1287 16188 : long i, ns = 0;
1288 16188 : if (!inC(x)) return -1;
1289 1176 : for (i = 0; i < n; i++)
1290 : {
1291 966 : GEN t = gsub(x, gel(X,i));
1292 966 : if (!inC(t)) return -1;
1293 952 : t = gabs(t, DEFAULTPREC);
1294 952 : if (!D || gcmp(t,D) < 0) { ns = i; D = t; }
1295 : }
1296 : /* X[ns] is closest to x */
1297 210 : return ns;
1298 : }
1299 : /* X,Y are "spec" GEN vectors with n > 0 components ( at X[0], ... X[n-1] ) */
1300 : GEN
1301 16223 : polintspec(GEN X, GEN Y, GEN x, long n, long *pe)
1302 : {
1303 : long i, m, ns;
1304 16223 : pari_sp av = avma, av2;
1305 16223 : GEN y, c, d, dy = NULL; /* gcc -Wall */
1306 :
1307 16223 : if (pe) *pe = -HIGHEXPOBIT;
1308 16223 : if (n == 1) return gmul(gel(Y,0), Rg_get_1(x));
1309 16188 : if (!X) X = identity_ZV(n) + 1;
1310 16188 : av2 = avma;
1311 16188 : ns = check_dy(X, x, n); if (ns < 0) { pe = NULL; ns = 0; }
1312 16188 : c = cgetg(n+1, t_VEC);
1313 81031 : d = cgetg(n+1, t_VEC); for (i=0; i<n; i++) gel(c,i+1) = gel(d,i+1) = gel(Y,i);
1314 16188 : y = gel(d,ns+1);
1315 : /* divided differences */
1316 64836 : for (m = 1; m < n; m++)
1317 : {
1318 146238 : for (i = 0; i < n-m; i++)
1319 : {
1320 97590 : GEN ho = gsub(gel(X,i),x), hp = gsub(gel(X,i+m),x), den = gsub(ho,hp);
1321 97590 : if (gequal0(den))
1322 : {
1323 7 : char *x1 = stack_sprintf("X[%ld]", i+1);
1324 7 : char *x2 = stack_sprintf("X[%ld]", i+m+1);
1325 7 : pari_err_DOMAIN("polinterpolate",x1,"=",strtoGENstr(x2), X);
1326 : }
1327 97583 : den = gdiv(gsub(gel(c,i+2),gel(d,i+1)), den);
1328 97583 : gel(c,i+1) = gmul(ho,den);
1329 97583 : gel(d,i+1) = gmul(hp,den);
1330 : }
1331 48648 : dy = (2*ns < n-m)? gel(c,ns+1): gel(d,ns--);
1332 48648 : y = gadd(y,dy);
1333 48648 : if (gc_needed(av2,2))
1334 : {
1335 0 : if (DEBUGMEM>1) pari_warn(warnmem,"polint, %ld/%ld",m,n-1);
1336 0 : gerepileall(av2, 4, &y, &c, &d, &dy);
1337 : }
1338 : }
1339 16181 : if (pe && inC(dy)) *pe = gexpo(dy);
1340 16181 : return gerepileupto(av, y);
1341 : }
1342 :
1343 : GEN
1344 329 : polint_i(GEN X, GEN Y, GEN t, long *pe)
1345 : {
1346 329 : long lx = lg(X), vt;
1347 :
1348 329 : if (! is_vec_t(typ(X))) pari_err_TYPE("polinterpolate",X);
1349 329 : if (Y)
1350 : {
1351 301 : if (! is_vec_t(typ(Y))) pari_err_TYPE("polinterpolate",Y);
1352 301 : if (lx != lg(Y)) pari_err_DIM("polinterpolate");
1353 : }
1354 : else
1355 : {
1356 28 : Y = X;
1357 28 : X = NULL;
1358 : }
1359 329 : if (pe) *pe = -HIGHEXPOBIT;
1360 329 : vt = t? gvar(t): 0;
1361 329 : if (vt != NO_VARIABLE)
1362 : { /* formal interpolation */
1363 : pari_sp av;
1364 224 : long v0, vY = gvar(Y);
1365 : GEN P;
1366 224 : if (X) vY = varnmax(vY, gvar(X));
1367 : /* shortcut */
1368 224 : if (varncmp(vY, vt) > 0 && (!t || gequalX(t))) return RgV_polint(X, Y, vt);
1369 84 : av = avma;
1370 : /* first interpolate in high priority variable, then substitute t */
1371 84 : v0 = fetch_var_higher();
1372 84 : P = RgV_polint(X, Y, v0);
1373 84 : P = gsubst(P, v0, t? t: pol_x(0));
1374 84 : (void)delete_var();
1375 84 : return gerepileupto(av, P);
1376 : }
1377 : /* numerical interpolation */
1378 105 : if (lx == 1) return Rg_get_0(t);
1379 91 : return polintspec(X? X+1: NULL,Y+1,t,lx-1, pe);
1380 : }
1381 : GEN
1382 329 : polint(GEN X, GEN Y, GEN t, GEN *pe)
1383 : {
1384 : long e;
1385 329 : GEN p = polint_i(X, Y, t, &e);
1386 322 : if (pe) *pe = stoi(e);
1387 322 : return p;
1388 : }
1389 :
1390 : /********************************************************************/
1391 : /** **/
1392 : /** MODREVERSE **/
1393 : /** **/
1394 : /********************************************************************/
1395 : static void
1396 7 : err_reverse(GEN x, GEN T)
1397 : {
1398 7 : pari_err_DOMAIN("modreverse","deg(minpoly(z))", "<", stoi(degpol(T)),
1399 : mkpolmod(x,T));
1400 0 : }
1401 :
1402 : /* return y such that Mod(y, charpoly(Mod(a,T)) = Mod(a,T) */
1403 : GEN
1404 175 : RgXQ_reverse(GEN a, GEN T)
1405 : {
1406 175 : pari_sp av = avma;
1407 175 : long n = degpol(T);
1408 : GEN y;
1409 :
1410 175 : if (n <= 1) {
1411 7 : if (n <= 0) return gcopy(a);
1412 7 : return gerepileupto(av, gneg(gdiv(gel(T,2), gel(T,3))));
1413 : }
1414 168 : if (typ(a) != t_POL || !signe(a)) err_reverse(a,T);
1415 168 : y = RgXV_to_RgM(RgXQ_powers(a,n-1,T), n);
1416 168 : y = RgM_solve(y, col_ei(n, 2));
1417 168 : if (!y) err_reverse(a,T);
1418 161 : return gerepilecopy(av, RgV_to_RgX(y, varn(T)));
1419 : }
1420 : GEN
1421 6104 : QXQ_reverse(GEN a, GEN T)
1422 : {
1423 6104 : pari_sp av = avma;
1424 6104 : long n = degpol(T);
1425 : GEN y;
1426 :
1427 6104 : if (n <= 1) {
1428 14 : if (n <= 0) return gcopy(a);
1429 14 : return gerepileupto(av, gneg(gdiv(gel(T,2), gel(T,3))));
1430 : }
1431 6090 : if (typ(a) != t_POL || !signe(a)) err_reverse(a,T);
1432 6090 : if (gequalX(a)) return gcopy(a);
1433 5936 : y = RgXV_to_RgM(QXQ_powers(a,n-1,T), n);
1434 5936 : y = QM_gauss(y, col_ei(n, 2));
1435 5936 : if (!y) err_reverse(a,T);
1436 5936 : return gerepilecopy(av, RgV_to_RgX(y, varn(T)));
1437 : }
1438 :
1439 : GEN
1440 28 : modreverse(GEN x)
1441 : {
1442 : long v, n;
1443 : GEN T, a;
1444 :
1445 28 : if (typ(x)!=t_POLMOD) pari_err_TYPE("modreverse",x);
1446 28 : T = gel(x,1); n = degpol(T); if (n <= 0) return gcopy(x);
1447 21 : a = gel(x,2);
1448 21 : v = varn(T);
1449 21 : retmkpolmod(RgXQ_reverse(a, T),
1450 : (n==1)? gsub(pol_x(v), a): RgXQ_charpoly(a, T, v));
1451 : }
1452 :
1453 : /********************************************************************/
1454 : /** **/
1455 : /** MERGESORT **/
1456 : /** **/
1457 : /********************************************************************/
1458 : static int
1459 77 : cmp_small(GEN x, GEN y) {
1460 77 : long a = (long)x, b = (long)y;
1461 77 : return a>b? 1: (a<b? -1: 0);
1462 : }
1463 :
1464 : static int
1465 295015 : veccmp(void *data, GEN x, GEN y)
1466 : {
1467 295015 : GEN k = (GEN)data;
1468 295015 : long i, s, lk = lg(k), lx = minss(lg(x), lg(y));
1469 :
1470 295015 : if (!is_vec_t(typ(x))) pari_err_TYPE("lexicographic vecsort",x);
1471 295015 : if (!is_vec_t(typ(y))) pari_err_TYPE("lexicographic vecsort",y);
1472 306684 : for (i=1; i<lk; i++)
1473 : {
1474 295043 : long c = k[i];
1475 295043 : if (c >= lx)
1476 14 : pari_err_TYPE("lexicographic vecsort, index too large", stoi(c));
1477 295029 : s = lexcmp(gel(x,c), gel(y,c));
1478 295029 : if (s) return s;
1479 : }
1480 11641 : return 0;
1481 : }
1482 :
1483 : /* return permutation sorting v[1..n], removing duplicates. Assume n > 0 */
1484 : static GEN
1485 2275457 : gen_sortspec_uniq(GEN v, long n, void *E, int (*cmp)(void*,GEN,GEN))
1486 : {
1487 : pari_sp av;
1488 : long NX, nx, ny, m, ix, iy, i;
1489 : GEN x, y, w, W;
1490 : int s;
1491 2275457 : switch(n)
1492 : {
1493 99301 : case 1: return mkvecsmall(1);
1494 958749 : case 2:
1495 958749 : s = cmp(E,gel(v,1),gel(v,2));
1496 958754 : if (s < 0) return mkvecsmall2(1,2);
1497 328644 : else if (s > 0) return mkvecsmall2(2,1);
1498 31661 : return mkvecsmall(1);
1499 285187 : case 3:
1500 285187 : s = cmp(E,gel(v,1),gel(v,2));
1501 285187 : if (s < 0) {
1502 182162 : s = cmp(E,gel(v,2),gel(v,3));
1503 182162 : if (s < 0) return mkvecsmall3(1,2,3);
1504 67278 : else if (s == 0) return mkvecsmall2(1,2);
1505 66549 : s = cmp(E,gel(v,1),gel(v,3));
1506 66549 : if (s < 0) return mkvecsmall3(1,3,2);
1507 34776 : else if (s > 0) return mkvecsmall3(3,1,2);
1508 2408 : return mkvecsmall2(1,2);
1509 103025 : } else if (s > 0) {
1510 98846 : s = cmp(E,gel(v,1),gel(v,3));
1511 98846 : if (s < 0) return mkvecsmall3(2,1,3);
1512 67298 : else if (s == 0) return mkvecsmall2(2,1);
1513 65129 : s = cmp(E,gel(v,2),gel(v,3));
1514 65129 : if (s < 0) return mkvecsmall3(2,3,1);
1515 31850 : else if (s > 0) return mkvecsmall3(3,2,1);
1516 721 : return mkvecsmall2(2,1);
1517 : } else {
1518 4179 : s = cmp(E,gel(v,1),gel(v,3));
1519 4179 : if (s < 0) return mkvecsmall2(1,3);
1520 1855 : else if (s == 0) return mkvecsmall(1);
1521 1001 : return mkvecsmall2(3,1);
1522 : }
1523 : }
1524 932220 : NX = nx = n>>1; ny = n-nx;
1525 932220 : av = avma;
1526 932220 : x = gen_sortspec_uniq(v, nx,E,cmp); nx = lg(x)-1;
1527 932263 : y = gen_sortspec_uniq(v+NX,ny,E,cmp); ny = lg(y)-1;
1528 932265 : w = cgetg(n+1, t_VECSMALL);
1529 932263 : m = ix = iy = 1;
1530 10267219 : while (ix<=nx && iy<=ny)
1531 : {
1532 9334959 : s = cmp(E, gel(v,x[ix]), gel(v,y[iy]+NX));
1533 9334956 : if (s < 0)
1534 4425538 : w[m++] = x[ix++];
1535 4909418 : else if (s > 0)
1536 3760576 : w[m++] = y[iy++]+NX;
1537 : else {
1538 1148842 : w[m++] = x[ix++];
1539 1148842 : iy++;
1540 : }
1541 : }
1542 1421383 : while (ix<=nx) w[m++] = x[ix++];
1543 2333927 : while (iy<=ny) w[m++] = y[iy++]+NX;
1544 932260 : set_avma(av);
1545 932258 : W = cgetg(m, t_VECSMALL);
1546 12157995 : for (i = 1; i < m; i++) W[i] = w[i];
1547 932252 : return W;
1548 : }
1549 :
1550 : /* return permutation sorting v[1..n]. Assume n > 0 */
1551 : static GEN
1552 190744836 : gen_sortspec(GEN v, long n, void *E, int (*cmp)(void*,GEN,GEN))
1553 : {
1554 : long nx, ny, m, ix, iy;
1555 : GEN x, y, w;
1556 190744836 : switch(n)
1557 : {
1558 5783660 : case 1:
1559 5783660 : (void)cmp(E,gel(v,1),gel(v,1)); /* check for type error */
1560 5783777 : return mkvecsmall(1);
1561 78874615 : case 2:
1562 135372210 : return cmp(E,gel(v,1),gel(v,2)) <= 0? mkvecsmall2(1,2)
1563 135372112 : : mkvecsmall2(2,1);
1564 37359377 : case 3:
1565 37359377 : if (cmp(E,gel(v,1),gel(v,2)) <= 0) {
1566 27135342 : if (cmp(E,gel(v,2),gel(v,3)) <= 0) return mkvecsmall3(1,2,3);
1567 12085783 : return (cmp(E,gel(v,1),gel(v,3)) <= 0)? mkvecsmall3(1,3,2)
1568 12085795 : : mkvecsmall3(3,1,2);
1569 : } else {
1570 10224041 : if (cmp(E,gel(v,1),gel(v,3)) <= 0) return mkvecsmall3(2,1,3);
1571 10723082 : return (cmp(E,gel(v,2),gel(v,3)) <= 0)? mkvecsmall3(2,3,1)
1572 10723084 : : mkvecsmall3(3,2,1);
1573 : }
1574 : }
1575 68727184 : nx = n>>1; ny = n-nx;
1576 68727184 : w = cgetg(n+1,t_VECSMALL);
1577 68729479 : x = gen_sortspec(v, nx,E,cmp);
1578 68729452 : y = gen_sortspec(v+nx,ny,E,cmp);
1579 68729454 : m = ix = iy = 1;
1580 459788855 : while (ix<=nx && iy<=ny)
1581 391059406 : if (cmp(E, gel(v,x[ix]), gel(v,y[iy]+nx))<=0)
1582 216887448 : w[m++] = x[ix++];
1583 : else
1584 174171953 : w[m++] = y[iy++]+nx;
1585 104683495 : while (ix<=nx) w[m++] = x[ix++];
1586 174713997 : while (iy<=ny) w[m++] = y[iy++]+nx;
1587 68729449 : set_avma((pari_sp)w); return w;
1588 : }
1589 :
1590 : static void
1591 46514851 : init_sort(GEN *x, long *tx, long *lx)
1592 : {
1593 46514851 : *tx = typ(*x);
1594 46514851 : if (*tx == t_LIST)
1595 : {
1596 35 : if (list_typ(*x)!=t_LIST_RAW) pari_err_TYPE("sort",*x);
1597 35 : *x = list_data(*x);
1598 35 : *lx = *x? lg(*x): 1;
1599 : } else {
1600 46514816 : if (!is_matvec_t(*tx) && *tx != t_VECSMALL) pari_err_TYPE("gen_sort",*x);
1601 46514811 : *lx = lg(*x);
1602 : }
1603 46514846 : }
1604 :
1605 : /* (x o y)[1..lx-1], destroy y */
1606 : INLINE GEN
1607 3159566 : sort_extract(GEN x, GEN y, long tx, long lx)
1608 : {
1609 : long i;
1610 3159566 : switch(tx)
1611 : {
1612 7 : case t_VECSMALL:
1613 35 : for (i=1; i<lx; i++) y[i] = x[y[i]];
1614 7 : break;
1615 7 : case t_LIST:
1616 7 : settyp(y,t_VEC);
1617 35 : for (i=1; i<lx; i++) gel(y,i) = gel(x,y[i]);
1618 7 : return gtolist(y);
1619 3159552 : default:
1620 3159552 : settyp(y,tx);
1621 9716455 : for (i=1; i<lx; i++) gel(y,i) = gcopy(gel(x,y[i]));
1622 : }
1623 3159629 : return y;
1624 : }
1625 :
1626 : static GEN
1627 1946569 : triv_sort(long tx) { return tx == t_LIST? mklist(): cgetg(1, tx); }
1628 : /* Sort the vector x, using cmp to compare entries. */
1629 : GEN
1630 348495 : gen_sort_uniq(GEN x, void *E, int (*cmp)(void*,GEN,GEN))
1631 : {
1632 : long tx, lx;
1633 : GEN y;
1634 :
1635 348495 : init_sort(&x, &tx, &lx);
1636 348495 : if (lx==1) return triv_sort(tx);
1637 343798 : y = gen_sortspec_uniq(x,lx-1,E,cmp);
1638 343798 : return sort_extract(x, y, tx, lg(y)); /* lg(y) <= lx */
1639 : }
1640 : /* Sort the vector x, using cmp to compare entries. */
1641 : GEN
1642 4757608 : gen_sort(GEN x, void *E, int (*cmp)(void*,GEN,GEN))
1643 : {
1644 : long tx, lx;
1645 : GEN y;
1646 :
1647 4757608 : init_sort(&x, &tx, &lx);
1648 4757604 : if (lx==1) return triv_sort(tx);
1649 2815732 : y = gen_sortspec(x,lx-1,E,cmp);
1650 2815773 : return sort_extract(x, y, tx, lx);
1651 : }
1652 : /* indirect sort: return the permutation that would sort x */
1653 : GEN
1654 73998 : gen_indexsort_uniq(GEN x, void *E, int (*cmp)(void*,GEN,GEN))
1655 : {
1656 : long tx, lx;
1657 73998 : init_sort(&x, &tx, &lx);
1658 74000 : if (lx==1) return cgetg(1, t_VECSMALL);
1659 67202 : return gen_sortspec_uniq(x,lx-1,E,cmp);
1660 : }
1661 : /* indirect sort: return the permutation that would sort x */
1662 : GEN
1663 830005 : gen_indexsort(GEN x, void *E, int (*cmp)(void*,GEN,GEN))
1664 : {
1665 : long tx, lx;
1666 830005 : init_sort(&x, &tx, &lx);
1667 830004 : if (lx==1) return cgetg(1, t_VECSMALL);
1668 829696 : return gen_sortspec(x,lx-1,E,cmp);
1669 : }
1670 :
1671 : /* Sort the vector x in place, using cmp to compare entries */
1672 : void
1673 40102444 : gen_sort_inplace(GEN x, void *E, int (*cmp)(void*,GEN,GEN), GEN *perm)
1674 : {
1675 : long tx, lx, i;
1676 40102444 : pari_sp av = avma;
1677 : GEN y;
1678 :
1679 40102444 : init_sort(&x, &tx, &lx);
1680 40102443 : if (lx<=2)
1681 : {
1682 566517 : if (perm) *perm = lx == 1? cgetg(1, t_VECSMALL): mkvecsmall(1);
1683 566517 : return;
1684 : }
1685 39535926 : y = gen_sortspec(x,lx-1, E, cmp);
1686 39535917 : if (perm)
1687 : {
1688 15351 : GEN z = new_chunk(lx);
1689 121324 : for (i=1; i<lx; i++) gel(z,i) = gel(x,y[i]);
1690 121324 : for (i=1; i<lx; i++) gel(x,i) = gel(z,i);
1691 15351 : *perm = y;
1692 15351 : set_avma((pari_sp)y);
1693 : } else {
1694 283451559 : for (i=1; i<lx; i++) gel(y,i) = gel(x,y[i]);
1695 283451569 : for (i=1; i<lx; i++) gel(x,i) = gel(y,i);
1696 39520566 : set_avma(av);
1697 : }
1698 : }
1699 : GEN
1700 402339 : gen_sort_shallow(GEN x, void *E, int (*cmp)(void*,GEN,GEN))
1701 : {
1702 : long tx, lx, i;
1703 : pari_sp av;
1704 : GEN y, z;
1705 :
1706 402339 : init_sort(&x, &tx, &lx);
1707 402339 : if (lx<=2) return x;
1708 240968 : z = cgetg(lx, tx); av = avma;
1709 240968 : y = gen_sortspec(x,lx-1, E, cmp);
1710 1246616 : for (i=1; i<lx; i++) gel(z,i) = gel(x,y[i]);
1711 240968 : return gc_const(av, z);
1712 : }
1713 :
1714 : static int
1715 7909 : closurecmp(void *data, GEN x, GEN y)
1716 : {
1717 7909 : pari_sp av = avma;
1718 7909 : long s = gsigne(closure_callgen2((GEN)data, x,y));
1719 7909 : set_avma(av); return s;
1720 : }
1721 : static void
1722 133 : check_positive_entries(GEN k)
1723 : {
1724 133 : long i, l = lg(k);
1725 301 : for (i=1; i<l; i++)
1726 168 : if (k[i] <= 0) pari_err_DOMAIN("sort_function", "index", "<", gen_0, stoi(k[i]));
1727 133 : }
1728 :
1729 : typedef int (*CMP_FUN)(void*,GEN,GEN);
1730 : /* return NULL if t_CLOSURE k is a "key" (arity 1) and not a sorting func */
1731 : static CMP_FUN
1732 126861 : sort_function(void **E, GEN x, GEN k)
1733 : {
1734 126861 : int (*cmp)(GEN,GEN) = &lexcmp;
1735 126861 : long tx = typ(x);
1736 126861 : if (!k)
1737 : {
1738 126175 : *E = (void*)((typ(x) == t_VECSMALL)? cmp_small: cmp);
1739 126175 : return &cmp_nodata;
1740 : }
1741 686 : if (tx == t_VECSMALL) pari_err_TYPE("sort_function", x);
1742 672 : switch(typ(k))
1743 : {
1744 98 : case t_INT: k = mkvecsmall(itos(k)); break;
1745 35 : case t_VEC: case t_COL: k = ZV_to_zv(k); break;
1746 0 : case t_VECSMALL: break;
1747 539 : case t_CLOSURE:
1748 539 : if (closure_is_variadic(k))
1749 0 : pari_err_TYPE("sort_function, variadic cmpf",k);
1750 539 : *E = (void*)k;
1751 539 : switch(closure_arity(k))
1752 : {
1753 35 : case 1: return NULL; /* wrt key */
1754 504 : case 2: return &closurecmp;
1755 0 : default: pari_err_TYPE("sort_function, cmpf arity != 1, 2",k);
1756 : }
1757 0 : default: pari_err_TYPE("sort_function",k);
1758 : }
1759 133 : check_positive_entries(k);
1760 133 : *E = (void*)k; return &veccmp;
1761 : }
1762 :
1763 : #define cmp_IND 1
1764 : #define cmp_LEX 2 /* FIXME: backward compatibility, ignored */
1765 : #define cmp_REV 4
1766 : #define cmp_UNIQ 8
1767 : GEN
1768 735 : vecsort0(GEN x, GEN k, long flag)
1769 : {
1770 : void *E;
1771 735 : int (*CMP)(void*,GEN,GEN) = sort_function(&E, x, k);
1772 :
1773 728 : if (flag < 0 || flag > (cmp_REV|cmp_LEX|cmp_IND|cmp_UNIQ))
1774 0 : pari_err_FLAG("vecsort");
1775 728 : if (!CMP)
1776 : { /* wrt key: precompute all values, O(n) calls instead of O(n log n) */
1777 28 : pari_sp av = avma;
1778 : GEN v, y;
1779 : long i, tx, lx;
1780 28 : init_sort(&x, &tx, &lx);
1781 28 : if (lx == 1) return flag&cmp_IND? cgetg(1,t_VECSMALL): triv_sort(tx);
1782 28 : v = cgetg(lx, t_VEC);
1783 140 : for (i = 1; i < lx; i++) gel(v,i) = closure_callgen1(k, gel(x,i));
1784 28 : y = vecsort0(v, NULL, flag | cmp_IND);
1785 28 : y = flag&cmp_IND? y: sort_extract(x, y, tx, lg(y));
1786 28 : return gerepileupto(av, y);
1787 : }
1788 700 : if (flag&cmp_UNIQ)
1789 35 : x = flag&cmp_IND? gen_indexsort_uniq(x, E, CMP): gen_sort_uniq(x, E, CMP);
1790 : else
1791 665 : x = flag&cmp_IND? gen_indexsort(x, E, CMP): gen_sort(x, E, CMP);
1792 686 : if (flag & cmp_REV)
1793 : { /* reverse order */
1794 35 : GEN y = x;
1795 35 : if (typ(x)==t_LIST) { y = list_data(x); if (!y) return x; }
1796 28 : vecreverse_inplace(y);
1797 : }
1798 679 : return x;
1799 : }
1800 :
1801 : GEN
1802 204793 : indexsort(GEN x) { return gen_indexsort(x, (void*)&gcmp, cmp_nodata); }
1803 : GEN
1804 0 : indexlexsort(GEN x) { return gen_indexsort(x, (void*)&lexcmp, cmp_nodata); }
1805 : GEN
1806 42 : indexvecsort(GEN x, GEN k)
1807 : {
1808 42 : if (typ(k) != t_VECSMALL) pari_err_TYPE("vecsort",k);
1809 42 : return gen_indexsort(x, (void*)k, &veccmp);
1810 : }
1811 :
1812 : GEN
1813 1836219 : sort(GEN x) { return gen_sort(x, (void*)gcmp, cmp_nodata); }
1814 : GEN
1815 127979 : lexsort(GEN x) { return gen_sort(x, (void*)lexcmp, cmp_nodata); }
1816 : GEN
1817 2954 : vecsort(GEN x, GEN k)
1818 : {
1819 2954 : if (typ(k) != t_VECSMALL) pari_err_TYPE("vecsort",k);
1820 2954 : return gen_sort(x, (void*)k, &veccmp);
1821 : }
1822 : /* adapted from gen_search; don't export: keys of T[i] should be precomputed */
1823 : static long
1824 7 : key_search(GEN T, GEN x, GEN code)
1825 : {
1826 7 : long u = lg(T)-1, i, l, s;
1827 :
1828 7 : if (!u) return 0;
1829 7 : l = 1; x = closure_callgen1(code, x);
1830 : do
1831 : {
1832 14 : i = (l+u)>>1; s = lexcmp(x, closure_callgen1(code, gel(T,i)));
1833 14 : if (!s) return i;
1834 7 : if (s<0) u=i-1; else l=i+1;
1835 7 : } while (u>=l);
1836 0 : return 0;
1837 : }
1838 : long
1839 126126 : vecsearch(GEN v, GEN x, GEN k)
1840 : {
1841 126126 : pari_sp av = avma;
1842 : long r;
1843 : void *E;
1844 126126 : int (*CMP)(void*,GEN,GEN) = sort_function(&E, v, k);
1845 126119 : switch(typ(v))
1846 : {
1847 21 : case t_VECSMALL: x = (GEN)itos(x); break;
1848 126077 : case t_VEC: case t_COL: break;
1849 21 : case t_LIST:
1850 21 : if (list_typ(v)==t_LIST_RAW)
1851 : {
1852 21 : v = list_data(v); if (!v) v = cgetg(1, t_VEC);
1853 21 : break;
1854 : }
1855 : /* fall through */
1856 : default:
1857 0 : pari_err_TYPE("vecsearch", v);
1858 : }
1859 126119 : r = CMP? gen_search(v, x, E, CMP): key_search(v, x, k);
1860 126119 : return gc_long(av, r < 0? 0: r);
1861 : }
1862 :
1863 : GEN
1864 3759 : ZV_indexsort(GEN L) { return gen_indexsort(L, (void*)&cmpii, &cmp_nodata); }
1865 : GEN
1866 63 : ZV_sort(GEN L) { return gen_sort(L, (void*)&cmpii, &cmp_nodata); }
1867 : GEN
1868 61299 : ZV_sort_uniq(GEN L) { return gen_sort_uniq(L, (void*)&cmpii, &cmp_nodata); }
1869 : void
1870 1221322 : ZV_sort_inplace(GEN L) { gen_sort_inplace(L, (void*)&cmpii, &cmp_nodata,NULL); }
1871 : GEN
1872 37652 : ZV_sort_uniq_shallow(GEN L)
1873 : {
1874 37652 : GEN v = gen_indexsort_uniq(L, (void*)&cmpii, &cmp_nodata);
1875 37652 : return vecpermute(L, v);
1876 : }
1877 : GEN
1878 1288 : ZV_sort_shallow(GEN L)
1879 : {
1880 1288 : GEN v = gen_indexsort(L, (void*)&cmpii, &cmp_nodata);
1881 1288 : return vecpermute(L, v);
1882 : }
1883 :
1884 : GEN
1885 1127 : vec_equiv(GEN F)
1886 : {
1887 1127 : pari_sp av = avma;
1888 1127 : long j, k, L = lg(F);
1889 1127 : GEN w = cgetg(L, t_VEC);
1890 1127 : GEN perm = gen_indexsort(F, (void*)&cmp_universal, cmp_nodata);
1891 3213 : for (j = k = 1; j < L;)
1892 : {
1893 2086 : GEN v = cgetg(L, t_VECSMALL);
1894 2086 : long l = 1, o = perm[j];
1895 2086 : v[l++] = o;
1896 4921 : for (j++; j < L; v[l++] = perm[j++])
1897 3794 : if (!gequal(gel(F,o), gel(F, perm[j]))) break;
1898 2086 : setlg(v, l); gel(w, k++) = v;
1899 : }
1900 1127 : setlg(w, k); return gerepilecopy(av,w);
1901 : }
1902 :
1903 : GEN
1904 21280 : vec_reduce(GEN v, GEN *pE)
1905 : {
1906 21280 : GEN E, F, P = gen_indexsort(v, (void*)cmp_universal, cmp_nodata);
1907 : long i, m, l;
1908 21280 : F = cgetg_copy(v, &l);
1909 21280 : *pE = E = cgetg(l, t_VECSMALL);
1910 54243 : for (i = m = 1; i < l;)
1911 : {
1912 32963 : GEN u = gel(v, P[i]);
1913 : long k;
1914 38311 : for(k = i + 1; k < l; k++)
1915 17038 : if (cmp_universal(gel(v, P[k]), u)) break;
1916 32963 : E[m] = k - i; gel(F, m) = u; i = k; m++;
1917 : }
1918 21280 : setlg(F, m);
1919 21280 : setlg(E, m); return F;
1920 : }
1921 :
1922 : /********************************************************************/
1923 : /** SEARCH IN SORTED VECTOR **/
1924 : /********************************************************************/
1925 : /* index of x in table T, 0 otherwise */
1926 : long
1927 1126999 : tablesearch(GEN T, GEN x, int (*cmp)(GEN,GEN))
1928 : {
1929 1126999 : long l = 1, u = lg(T)-1, i, s;
1930 :
1931 8170712 : while (u>=l)
1932 : {
1933 8125082 : i = (l+u)>>1; s = cmp(x, gel(T,i));
1934 8125083 : if (!s) return i;
1935 7043713 : if (s<0) u=i-1; else l=i+1;
1936 : }
1937 45630 : return 0;
1938 : }
1939 :
1940 : /* looks if x belongs to the set T and returns the index if yes, 0 if no */
1941 : long
1942 23575090 : gen_search(GEN T, GEN x, void *data, int (*cmp)(void*,GEN,GEN))
1943 : {
1944 23575090 : long u = lg(T)-1, i, l, s;
1945 :
1946 23575090 : if (!u) return -1;
1947 23575062 : l = 1;
1948 : do
1949 : {
1950 110214580 : i = (l+u) >> 1; s = cmp(data, x, gel(T,i));
1951 110214580 : if (!s) return i;
1952 86693033 : if (s < 0) u = i-1; else l = i+1;
1953 86693033 : } while (u >= l);
1954 53515 : return -((s < 0)? i: i+1);
1955 : }
1956 :
1957 : long
1958 1073741 : ZV_search(GEN x, GEN y) { return tablesearch(x, y, cmpii); }
1959 :
1960 : long
1961 12896802 : zv_search(GEN T, long x)
1962 : {
1963 12896802 : long l = 1, u = lg(T)-1;
1964 52631252 : while (u>=l)
1965 : {
1966 42141002 : long i = (l+u)>>1;
1967 42141002 : if (x < T[i]) u = i-1;
1968 26072010 : else if (x > T[i]) l = i+1;
1969 2406552 : else return i;
1970 : }
1971 10490250 : return 0;
1972 : }
1973 :
1974 : /********************************************************************/
1975 : /** COMPARISON FUNCTIONS **/
1976 : /********************************************************************/
1977 : int
1978 645116265 : cmp_nodata(void *data, GEN x, GEN y)
1979 : {
1980 645116265 : int (*cmp)(GEN,GEN)=(int (*)(GEN,GEN)) data;
1981 645116265 : return cmp(x,y);
1982 : }
1983 :
1984 : /* assume x and y come from the same idealprimedec call (uniformizer unique) */
1985 : int
1986 3369147 : cmp_prime_over_p(GEN x, GEN y)
1987 : {
1988 3369147 : long k = pr_get_f(x) - pr_get_f(y); /* diff. between residue degree */
1989 222189 : return k? ((k > 0)? 1: -1)
1990 3591327 : : ZV_cmp(pr_get_gen(x), pr_get_gen(y));
1991 : }
1992 :
1993 : int
1994 673246 : cmp_prime_ideal(GEN x, GEN y)
1995 : {
1996 673246 : int k = cmpii(pr_get_p(x), pr_get_p(y));
1997 673245 : return k? k: cmp_prime_over_p(x,y);
1998 : }
1999 :
2000 : /* assume x and y are t_POL in the same variable whose coeffs can be
2001 : * compared (used to sort polynomial factorizations) */
2002 : int
2003 5522862 : gen_cmp_RgX(void *data, GEN x, GEN y)
2004 : {
2005 5522862 : int (*coeff_cmp)(GEN,GEN)=(int(*)(GEN,GEN))data;
2006 5522862 : long i, lx = lg(x), ly = lg(y);
2007 : int fl;
2008 5522862 : if (lx > ly) return 1;
2009 5485713 : if (lx < ly) return -1;
2010 12293196 : for (i=lx-1; i>1; i--)
2011 11519461 : if ((fl = coeff_cmp(gel(x,i), gel(y,i)))) return fl;
2012 773735 : return 0;
2013 : }
2014 :
2015 : static int
2016 3609 : cmp_RgX_Rg(GEN x, GEN y)
2017 : {
2018 3609 : long lx = lgpol(x), ly;
2019 3609 : if (lx > 1) return 1;
2020 0 : ly = gequal0(y) ? 0:1;
2021 0 : if (lx > ly) return 1;
2022 0 : if (lx < ly) return -1;
2023 0 : if (lx==0) return 0;
2024 0 : return gcmp(gel(x,2), y);
2025 : }
2026 : int
2027 111187 : cmp_RgX(GEN x, GEN y)
2028 : {
2029 111187 : if (typ(x) == t_POLMOD) x = gel(x,2);
2030 111187 : if (typ(y) == t_POLMOD) y = gel(y,2);
2031 111187 : if (typ(x) == t_POL) {
2032 55992 : if (typ(y) != t_POL) return cmp_RgX_Rg(x, y);
2033 : } else {
2034 55195 : if (typ(y) != t_POL) return gcmp(x,y);
2035 3364 : return - cmp_RgX_Rg(y,x);
2036 : }
2037 55747 : return gen_cmp_RgX((void*)&gcmp,x,y);
2038 : }
2039 :
2040 : int
2041 323638 : cmp_Flx(GEN x, GEN y)
2042 : {
2043 323638 : long i, lx = lg(x), ly = lg(y);
2044 323638 : if (lx > ly) return 1;
2045 307623 : if (lx < ly) return -1;
2046 548851 : for (i=lx-1; i>1; i--)
2047 463907 : if (uel(x,i) != uel(y,i)) return uel(x,i)<uel(y,i)? -1: 1;
2048 84944 : return 0;
2049 : }
2050 : /********************************************************************/
2051 : /** MERGE & SORT FACTORIZATIONS **/
2052 : /********************************************************************/
2053 : /* merge fx, fy two factorizations, whose 1st column is sorted in strictly
2054 : * increasing order wrt cmp */
2055 : GEN
2056 692242 : merge_factor(GEN fx, GEN fy, void *data, int (*cmp)(void *,GEN,GEN))
2057 : {
2058 692242 : GEN x = gel(fx,1), e = gel(fx,2), M, E;
2059 692242 : GEN y = gel(fy,1), f = gel(fy,2);
2060 692242 : long ix, iy, m, lx = lg(x), ly = lg(y), l = lx+ly-1;
2061 :
2062 692242 : M = cgetg(l, t_COL);
2063 692242 : E = cgetg(l, t_COL);
2064 :
2065 692242 : m = ix = iy = 1;
2066 10044878 : while (ix<lx && iy<ly)
2067 : {
2068 9352636 : int s = cmp(data, gel(x,ix), gel(y,iy));
2069 9352636 : if (s < 0)
2070 8717968 : { gel(M,m) = gel(x,ix); gel(E,m) = gel(e,ix); ix++; }
2071 634668 : else if (s == 0)
2072 : {
2073 95004 : GEN z = gel(x,ix), g = addii(gel(e,ix), gel(f,iy));
2074 95004 : iy++; ix++; if (!signe(g)) continue;
2075 11067 : gel(M,m) = z; gel(E,m) = g;
2076 : }
2077 : else
2078 539664 : { gel(M,m) = gel(y,iy); gel(E,m) = gel(f,iy); iy++; }
2079 9268699 : m++;
2080 : }
2081 4860282 : while (ix<lx) { gel(M,m) = gel(x,ix); gel(E,m) = gel(e,ix); ix++; m++; }
2082 937907 : while (iy<ly) { gel(M,m) = gel(y,iy); gel(E,m) = gel(f,iy); iy++; m++; }
2083 692242 : setlg(M, m);
2084 692242 : setlg(E, m); return mkmat2(M, E);
2085 : }
2086 :
2087 : GEN
2088 30952 : ZM_merge_factor(GEN A, GEN B)
2089 : {
2090 30952 : return merge_factor(A, B, (void*)&cmpii, cmp_nodata);
2091 : }
2092 :
2093 : /* merge two sorted vectors, removing duplicates. Shallow */
2094 : GEN
2095 451983 : merge_sort_uniq(GEN x, GEN y, void *data, int (*cmp)(void *,GEN,GEN))
2096 : {
2097 451983 : long i, j, k, lx = lg(x), ly = lg(y);
2098 451983 : GEN z = cgetg(lx + ly - 1, typ(x));
2099 451980 : i = j = k = 1;
2100 592061 : while (i<lx && j<ly)
2101 : {
2102 140082 : int s = cmp(data, gel(x,i), gel(y,j));
2103 140081 : if (s < 0)
2104 118975 : gel(z,k++) = gel(x,i++);
2105 21106 : else if (s > 0)
2106 21085 : gel(z,k++) = gel(y,j++);
2107 : else
2108 21 : { gel(z,k++) = gel(x,i++); j++; }
2109 : }
2110 803204 : while (i<lx) gel(z,k++) = gel(x,i++);
2111 579755 : while (j<ly) gel(z,k++) = gel(y,j++);
2112 451979 : setlg(z, k); return z;
2113 : }
2114 : /* in case of equal keys in x,y, take the key from x */
2115 : static GEN
2116 34741 : ZV_union_shallow_t(GEN x, GEN y, long t)
2117 : {
2118 34741 : long i, j, k, lx = lg(x), ly = lg(y);
2119 34741 : GEN z = cgetg(lx + ly - 1, t);
2120 34741 : i = j = k = 1;
2121 78225 : while (i<lx && j<ly)
2122 : {
2123 43484 : int s = cmpii(gel(x,i), gel(y,j));
2124 43484 : if (s < 0)
2125 23632 : gel(z,k++) = gel(x,i++);
2126 19852 : else if (s > 0)
2127 10353 : gel(z,k++) = gel(y,j++);
2128 : else
2129 9499 : { gel(z,k++) = gel(x,i++); j++; }
2130 : }
2131 41867 : while (i < lx) gel(z,k++) = gel(x,i++);
2132 70413 : while (j < ly) gel(z,k++) = gel(y,j++);
2133 34741 : setlg(z, k); return z;
2134 : }
2135 : GEN
2136 34559 : ZV_union_shallow(GEN x, GEN y)
2137 34559 : { return ZV_union_shallow_t(x, y, t_VEC); }
2138 : GEN
2139 182 : ZC_union_shallow(GEN x, GEN y)
2140 182 : { return ZV_union_shallow_t(x, y, t_COL); }
2141 :
2142 : /* sort generic factorization, in place */
2143 : GEN
2144 9886882 : sort_factor(GEN y, void *data, int (*cmp)(void *,GEN,GEN))
2145 : {
2146 : GEN a, b, A, B, w;
2147 : pari_sp av;
2148 : long n, i;
2149 :
2150 9886882 : a = gel(y,1); n = lg(a); if (n == 1) return y;
2151 9864246 : b = gel(y,2); av = avma;
2152 9864246 : A = new_chunk(n);
2153 9864106 : B = new_chunk(n);
2154 9863876 : w = gen_sortspec(a, n-1, data, cmp);
2155 29599456 : for (i=1; i<n; i++) { long k=w[i]; gel(A,i) = gel(a,k); gel(B,i) = gel(b,k); }
2156 29599725 : for (i=1; i<n; i++) { gel(a,i) = gel(A,i); gel(b,i) = gel(B,i); }
2157 9864110 : set_avma(av); return y;
2158 : }
2159 : /* sort polynomial factorization, in place */
2160 : GEN
2161 1988016 : sort_factor_pol(GEN y,int (*cmp)(GEN,GEN))
2162 : {
2163 1988016 : (void)sort_factor(y,(void*)cmp, &gen_cmp_RgX);
2164 1988024 : return y;
2165 : }
2166 :
2167 : /***********************************************************************/
2168 : /* */
2169 : /* SET OPERATIONS */
2170 : /* */
2171 : /***********************************************************************/
2172 : GEN
2173 227115 : gtoset(GEN x)
2174 : {
2175 : long lx;
2176 227115 : if (!x) return cgetg(1, t_VEC);
2177 227115 : switch(typ(x))
2178 : {
2179 227087 : case t_VEC:
2180 227087 : case t_COL: lx = lg(x); break;
2181 14 : case t_LIST:
2182 14 : if (list_typ(x)==t_LIST_MAP) return mapdomain(x);
2183 14 : x = list_data(x); lx = x? lg(x): 1; break;
2184 7 : case t_VECSMALL: lx = lg(x); x = zv_to_ZV(x); break;
2185 7 : default: return mkveccopy(x);
2186 : }
2187 227108 : if (lx==1) return cgetg(1,t_VEC);
2188 226933 : x = gen_sort_uniq(x, (void*)&cmp_universal, cmp_nodata);
2189 226933 : settyp(x, t_VEC); /* it may be t_COL */
2190 226933 : return x;
2191 : }
2192 :
2193 : long
2194 14 : setisset(GEN x)
2195 : {
2196 14 : long i, lx = lg(x);
2197 :
2198 14 : if (typ(x) != t_VEC) return 0;
2199 14 : if (lx == 1) return 1;
2200 70 : for (i=1; i<lx-1; i++)
2201 63 : if (cmp_universal(gel(x,i+1), gel(x,i)) <= 0) return 0;
2202 7 : return 1;
2203 : }
2204 :
2205 : long
2206 83622 : setsearch(GEN T, GEN y, long flag)
2207 : {
2208 : long i, lx;
2209 83622 : switch(typ(T))
2210 : {
2211 83608 : case t_VEC: lx = lg(T); break;
2212 7 : case t_LIST:
2213 7 : if (list_typ(T) != t_LIST_RAW) pari_err_TYPE("setsearch",T);
2214 7 : T = list_data(T); lx = T? lg(T): 1; break;
2215 7 : default: pari_err_TYPE("setsearch",T);
2216 : return 0; /*LCOV_EXCL_LINE*/
2217 : }
2218 83615 : if (lx==1) return flag? 1: 0;
2219 83615 : i = gen_search(T,y,(void*)cmp_universal,cmp_nodata);
2220 83615 : if (i > 0) return flag? 0: i;
2221 56 : return flag ? -i: 0;
2222 : }
2223 :
2224 : GEN
2225 7 : setunion_i(GEN x, GEN y)
2226 7 : { return merge_sort_uniq(x,y, (void*)cmp_universal, cmp_nodata); }
2227 :
2228 : GEN
2229 7 : setunion(GEN x, GEN y)
2230 : {
2231 7 : pari_sp av = avma;
2232 7 : if (typ(x) != t_VEC) pari_err_TYPE("setunion",x);
2233 7 : if (typ(y) != t_VEC) pari_err_TYPE("setunion",y);
2234 7 : return gerepilecopy(av, setunion_i(x, y));
2235 : }
2236 :
2237 : GEN
2238 14 : setdelta(GEN x, GEN y)
2239 : {
2240 14 : long ix = 1, iy = 1, iz = 1, lx = lg(x), ly = lg(y);
2241 14 : pari_sp av = avma;
2242 14 : GEN z = cgetg(lx + ly - 1,t_VEC);
2243 14 : if (typ(x) != t_VEC) pari_err_TYPE("setdelta",x);
2244 14 : if (typ(y) != t_VEC) pari_err_TYPE("setdelta",y);
2245 84 : while (ix < lx && iy < ly)
2246 : {
2247 70 : int c = cmp_universal(gel(x,ix), gel(y,iy));
2248 70 : if (c < 0) gel(z, iz++) = gel(x,ix++);
2249 42 : else if (c > 0) gel(z, iz++) = gel(y,iy++);
2250 28 : else { ix++; iy++; }
2251 : }
2252 21 : while (ix<lx) gel(z,iz++) = gel(x,ix++);
2253 14 : while (iy<ly) gel(z,iz++) = gel(y,iy++);
2254 14 : setlg(z,iz); return gerepilecopy(av,z);
2255 : }
2256 :
2257 : GEN
2258 7 : setintersect(GEN x, GEN y)
2259 : {
2260 7 : long ix = 1, iy = 1, iz = 1, lx = lg(x), ly = lg(y);
2261 7 : pari_sp av = avma;
2262 7 : GEN z = cgetg(lx,t_VEC);
2263 7 : if (typ(x) != t_VEC) pari_err_TYPE("setintersect",x);
2264 7 : if (typ(y) != t_VEC) pari_err_TYPE("setintersect",y);
2265 70 : while (ix < lx && iy < ly)
2266 : {
2267 63 : int c = cmp_universal(gel(x,ix), gel(y,iy));
2268 63 : if (c < 0) ix++;
2269 35 : else if (c > 0) iy++;
2270 21 : else { gel(z, iz++) = gel(x,ix); ix++; iy++; }
2271 : }
2272 7 : setlg(z,iz); return gerepilecopy(av,z);
2273 : }
2274 :
2275 : GEN
2276 259 : gen_setminus(GEN A, GEN B, int (*cmp)(GEN,GEN))
2277 : {
2278 259 : pari_sp ltop = avma;
2279 259 : long i = 1, j = 1, k = 1, lx = lg(A), ly = lg(B);
2280 259 : GEN diff = cgetg(lx,t_VEC);
2281 5481 : while (i < lx && j < ly)
2282 5222 : switch ( cmp(gel(A,i),gel(B,j)) )
2283 : {
2284 938 : case -1: gel(diff,k++) = gel(A,i++); break;
2285 2044 : case 1: j++; break;
2286 2240 : case 0: i++; break;
2287 : }
2288 308 : while (i < lx) gel(diff,k++) = gel(A,i++);
2289 259 : setlg(diff,k);
2290 259 : return gerepilecopy(ltop,diff);
2291 : }
2292 :
2293 : GEN
2294 259 : setminus(GEN x, GEN y)
2295 : {
2296 259 : if (typ(x) != t_VEC) pari_err_TYPE("setminus",x);
2297 259 : if (typ(y) != t_VEC) pari_err_TYPE("setminus",y);
2298 259 : return gen_setminus(x,y,cmp_universal);
2299 : }
2300 :
2301 : GEN
2302 21 : setbinop(GEN f, GEN x, GEN y)
2303 : {
2304 21 : pari_sp av = avma;
2305 21 : long i, j, lx, ly, k = 1;
2306 : GEN z;
2307 21 : if (typ(f) != t_CLOSURE || closure_arity(f) != 2 || closure_is_variadic(f))
2308 7 : pari_err_TYPE("setbinop [function needs exactly 2 arguments]",f);
2309 14 : lx = lg(x);
2310 14 : if (typ(x) != t_VEC) pari_err_TYPE("setbinop", x);
2311 14 : if (y == NULL) { /* assume x = y and f symmetric */
2312 7 : z = cgetg((((lx-1)*lx) >> 1) + 1, t_VEC);
2313 28 : for (i = 1; i < lx; i++)
2314 63 : for (j = i; j < lx; j++)
2315 42 : gel(z, k++) = closure_callgen2(f, gel(x,i),gel(x,j));
2316 : } else {
2317 7 : ly = lg(y);
2318 7 : if (typ(y) != t_VEC) pari_err_TYPE("setbinop", y);
2319 7 : z = cgetg((lx-1)*(ly-1) + 1, t_VEC);
2320 28 : for (i = 1; i < lx; i++)
2321 84 : for (j = 1; j < ly; j++)
2322 63 : gel(z, k++) = closure_callgen2(f, gel(x,i),gel(y,j));
2323 : }
2324 14 : return gerepileupto(av, gtoset(z));
2325 : }
|