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 : /*******************************************************************/
16 : /* */
17 : /* ROOTS OF COMPLEX POLYNOMIALS */
18 : /* (original code contributed by Xavier Gourdon, INRIA RR 1852) */
19 : /* */
20 : /*******************************************************************/
21 : #include "pari.h"
22 : #include "paripriv.h"
23 :
24 : #define DEBUGLEVEL DEBUGLEVEL_polroots
25 :
26 : static const double pariINFINITY = 1./0.;
27 :
28 : static long
29 1224421 : isvalidcoeff(GEN x)
30 : {
31 1224421 : switch (typ(x))
32 : {
33 1201207 : case t_INT: case t_REAL: case t_FRAC: return 1;
34 23200 : case t_COMPLEX: return isvalidcoeff(gel(x,1)) && isvalidcoeff(gel(x,2));
35 : }
36 14 : return 0;
37 : }
38 :
39 : static void
40 266505 : checkvalidpol(GEN p, const char *f)
41 : {
42 266505 : long i,n = lg(p);
43 1444505 : for (i=2; i<n; i++)
44 1178008 : if (!isvalidcoeff(gel(p,i))) pari_err_TYPE(f, gel(p,i));
45 266497 : }
46 :
47 : /********************************************************************/
48 : /** **/
49 : /** FAST ARITHMETIC over Z[i] **/
50 : /** **/
51 : /********************************************************************/
52 :
53 : static GEN
54 16459662 : ZX_to_ZiX(GEN Pr, GEN Pi)
55 : {
56 16459662 : long i, lr = lg(Pr), li = lg(Pi), l = maxss(lr, li), m = minss(lr, li);
57 16462043 : GEN P = cgetg(l, t_POL);
58 16467620 : P[1] = Pr[1];
59 66809921 : for(i = 2; i < m; i++)
60 50345059 : gel(P,i) = signe(gel(Pi,i)) ? mkcomplex(gel(Pr,i), gel(Pi,i))
61 50345059 : : gel(Pr,i);
62 22366882 : for( ; i < lr; i++)
63 5902020 : gel(P,i) = gel(Pr, i);
64 16499848 : for( ; i < li; i++)
65 34986 : gel(P,i) = mkcomplex(gen_0, gel(Pi, i));
66 16464862 : return normalizepol_lg(P, l);
67 : }
68 :
69 : static GEN
70 98790570 : ZiX_sqr(GEN P)
71 : {
72 98790570 : pari_sp av = avma;
73 : GEN Pr2, Pi2, Qr, Qi;
74 98790570 : GEN Pr = real_i(P), Pi = imag_i(P);
75 98785508 : if (signe(Pi)==0) return gerepileupto(av, ZX_sqr(Pr));
76 16524466 : if (signe(Pr)==0) return gerepileupto(av, ZX_neg(ZX_sqr(Pi)));
77 16467792 : Pr2 = ZX_sqr(Pr); Pi2 = ZX_sqr(Pi);
78 16460192 : Qr = ZX_sub(Pr2, Pi2);
79 16461936 : if (degpol(Pr)==degpol(Pi))
80 10706485 : Qi = ZX_sub(ZX_sqr(ZX_add(Pr, Pi)), ZX_add(Pr2, Pi2));
81 : else
82 5758578 : Qi = ZX_shifti(ZX_mul(Pr, Pi), 1);
83 16463293 : return gerepilecopy(av, ZX_to_ZiX(Qr, Qi));
84 : }
85 :
86 : static GEN
87 49388695 : graeffe(GEN p)
88 : {
89 49388695 : pari_sp av = avma;
90 : GEN p0, p1, s0, s1;
91 49388695 : long n = degpol(p);
92 :
93 49393420 : if (!n) return RgX_copy(p);
94 49393420 : RgX_even_odd(p, &p0, &p1);
95 : /* p = p0(x^2) + x p1(x^2) */
96 49401810 : s0 = ZiX_sqr(p0);
97 49404592 : s1 = ZiX_sqr(p1);
98 49403893 : return gerepileupto(av, RgX_sub(s0, RgX_shift_shallow(s1,1)));
99 : }
100 :
101 : GEN
102 5369 : ZX_graeffe(GEN p)
103 : {
104 5369 : pari_sp av = avma;
105 : GEN p0, p1, s0, s1;
106 5369 : long n = degpol(p);
107 :
108 5369 : if (!n) return ZX_copy(p);
109 5369 : RgX_even_odd(p, &p0, &p1);
110 : /* p = p0(x^2) + x p1(x^2) */
111 5369 : s0 = ZX_sqr(p0);
112 5369 : s1 = ZX_sqr(p1);
113 5369 : return gerepileupto(av, ZX_sub(s0, RgX_shift_shallow(s1,1)));
114 : }
115 : GEN
116 14 : polgraeffe(GEN p)
117 : {
118 14 : pari_sp av = avma;
119 : GEN p0, p1, s0, s1;
120 14 : long n = degpol(p);
121 :
122 14 : if (typ(p) != t_POL) pari_err_TYPE("polgraeffe",p);
123 14 : n = degpol(p);
124 14 : if (!n) return gcopy(p);
125 14 : RgX_even_odd(p, &p0, &p1);
126 : /* p = p0(x^2) + x p1(x^2) */
127 14 : s0 = RgX_sqr(p0);
128 14 : s1 = RgX_sqr(p1);
129 14 : return gerepileupto(av, RgX_sub(s0, RgX_shift_shallow(s1,1)));
130 : }
131 :
132 : /********************************************************************/
133 : /** **/
134 : /** MODULUS OF ROOTS **/
135 : /** **/
136 : /********************************************************************/
137 :
138 : /* Quick approximation to log2(|x|); first define y s.t. |y-x| < 2^-32 then
139 : * return y rounded to 2 ulp. In particular, if result < 2^21, absolute error
140 : * is bounded by 2^-31. If result > 2^21, it is correct to 2 ulp */
141 : static double
142 223849492 : mydbllog2i(GEN x)
143 : {
144 : #ifdef LONG_IS_64BIT
145 192552140 : const double W = 1/(4294967296. * 4294967296.); /* 2^-64 */
146 : #else
147 31297352 : const double W = 1/4294967296.; /*2^-32*/
148 : #endif
149 : GEN m;
150 223849492 : long lx = lgefint(x);
151 : double l;
152 223849492 : if (lx == 2) return -pariINFINITY;
153 223033717 : m = int_MSW(x);
154 223033717 : l = (double)(ulong)*m;
155 223033717 : if (lx == 3) return log2(l);
156 69920870 : l += ((double)(ulong)*int_precW(m)) * W;
157 : /* at least m = min(53,BIL) bits are correct in the mantissa, thus log2
158 : * is correct with error < log(1 + 2^-m) ~ 2^-m. Adding the correct
159 : * exponent BIL(lx-3) causes 1ulp further round-off error */
160 69920870 : return log2(l) + (double)(BITS_IN_LONG*(lx-3));
161 : }
162 :
163 : /* return log(|x|) or -pariINFINITY */
164 : static double
165 9526568 : mydbllogr(GEN x) {
166 9526568 : if (!signe(x)) return -pariINFINITY;
167 9526568 : return M_LN2*dbllog2r(x);
168 : }
169 :
170 : /* return log2(|x|) or -pariINFINITY */
171 : static double
172 55810054 : mydbllog2r(GEN x) {
173 55810054 : if (!signe(x)) return -pariINFINITY;
174 55374174 : return dbllog2r(x);
175 : }
176 : double
177 299691593 : dbllog2(GEN z)
178 : {
179 : double x, y;
180 299691593 : switch(typ(z))
181 : {
182 223742873 : case t_INT: return mydbllog2i(z);
183 22364 : case t_FRAC: return mydbllog2i(gel(z,1))-mydbllog2i(gel(z,2));
184 48879570 : case t_REAL: return mydbllog2r(z);
185 27046786 : default: /*t_COMPLEX*/
186 27046786 : x = dbllog2(gel(z,1));
187 27130121 : y = dbllog2(gel(z,2));
188 27129992 : if (x == -pariINFINITY) return y;
189 26886827 : if (y == -pariINFINITY) return x;
190 26685198 : if (fabs(x-y) > 10) return maxdd(x,y);
191 26069599 : return x + 0.5*log2(1 + exp2(2*(y-x)));
192 : }
193 : }
194 : static GEN /* beware overflow */
195 6532166 : dblexp(double x) { return fabs(x) < 100.? dbltor(exp(x)): mpexp(dbltor(x)); }
196 :
197 : /* find s such that A_h <= 2^s <= 2 A_i for one h and all i < n = deg(p),
198 : * with A_i := (binom(n,i) lc(p) / p_i) ^ 1/(n-i), and p = sum p_i X^i */
199 : static long
200 40862707 : findpower(GEN p)
201 : {
202 40862707 : double x, L, mins = pariINFINITY;
203 40862707 : long n = degpol(p),i;
204 :
205 40862010 : L = dbllog2(gel(p,n+2)); /* log2(lc * binom(n,i)) */
206 164255316 : for (i=n-1; i>=0; i--)
207 : {
208 123390954 : L += log2((double)(i+1) / (double)(n-i));
209 123390954 : x = dbllog2(gel(p,i+2));
210 123388767 : if (x != -pariINFINITY)
211 : {
212 122599899 : double s = (L - x) / (double)(n-i);
213 122599899 : if (s < mins) mins = s;
214 : }
215 : }
216 40864362 : i = (long)ceil(mins);
217 40864362 : if (i - mins > 1 - 1e-12) i--;
218 40864362 : return i;
219 : }
220 :
221 : /* returns the exponent for logmodulus(), from the Newton diagram */
222 : static long
223 5489266 : newton_polygon(GEN p, long k)
224 : {
225 5489266 : pari_sp av = avma;
226 5489266 : long n = degpol(p), i, j, h, l, *vertex = (long*)new_chunk(n+1);
227 5489249 : double *L = (double*)stack_malloc_align((n+1)*sizeof(double), sizeof(double));
228 :
229 : /* vertex[i] = 1 if i a vertex of convex hull, 0 otherwise */
230 26367859 : for (i=0; i<=n; i++) { L[i] = dbllog2(gel(p,2+i)); vertex[i] = 0; }
231 5489265 : vertex[0] = 1; /* sentinel */
232 19540625 : for (i=0; i < n; i=h)
233 : {
234 : double slope;
235 14051360 : h = i+1;
236 14056219 : while (L[i] == -pariINFINITY) { vertex[h] = 1; i = h; h = i+1; }
237 14051360 : slope = L[h] - L[i];
238 37849128 : for (j = i+2; j<=n; j++) if (L[j] != -pariINFINITY)
239 : {
240 23791572 : double pij = (L[j] - L[i])/(double)(j - i);
241 23791572 : if (slope < pij) { slope = pij; h = j; }
242 : }
243 14051360 : vertex[h] = 1;
244 : }
245 6202870 : h = k; while (!vertex[h]) h++;
246 5677121 : l = k-1; while (!vertex[l]) l--;
247 5489265 : set_avma(av);
248 5489307 : return (long)floor((L[h]-L[l])/(double)(h-l) + 0.5);
249 : }
250 :
251 : /* change z into z*2^e, where z is real or complex of real */
252 : static void
253 35405508 : myshiftrc(GEN z, long e)
254 : {
255 35405508 : if (typ(z)==t_COMPLEX)
256 : {
257 6402514 : if (signe(gel(z,1))) shiftr_inplace(gel(z,1), e);
258 6402533 : if (signe(gel(z,2))) shiftr_inplace(gel(z,2), e);
259 : }
260 : else
261 29002994 : if (signe(z)) shiftr_inplace(z, e);
262 35405623 : }
263 :
264 : /* return z*2^e, where z is integer or complex of integer (destroy z) */
265 : static GEN
266 135461363 : myshiftic(GEN z, long e)
267 : {
268 135461363 : if (typ(z)==t_COMPLEX)
269 : {
270 17910309 : gel(z,1) = signe(gel(z,1))? mpshift(gel(z,1),e): gen_0;
271 17908730 : gel(z,2) = mpshift(gel(z,2),e);
272 17906781 : return z;
273 : }
274 117551054 : return signe(z)? mpshift(z,e): gen_0;
275 : }
276 :
277 : static GEN
278 7003067 : RgX_gtofp_bit(GEN q, long bit) { return RgX_gtofp(q, nbits2prec(bit)); }
279 :
280 : static GEN
281 214127278 : mygprecrc(GEN x, long prec, long e)
282 : {
283 : GEN y;
284 214127278 : switch(typ(x))
285 : {
286 159823603 : case t_REAL:
287 159823603 : if (!signe(x)) return real_0_bit(e);
288 156609610 : return realprec(x) == prec? x: rtor(x, prec);
289 36994371 : case t_COMPLEX:
290 36994371 : y = cgetg(3,t_COMPLEX);
291 36993583 : gel(y,1) = mygprecrc(gel(x,1),prec,e);
292 36993787 : gel(y,2) = mygprecrc(gel(x,2),prec,e);
293 36994132 : return y;
294 17309304 : default: return x;
295 : }
296 : }
297 :
298 : /* gprec behaves badly with the zero for polynomials.
299 : The second parameter in mygprec is the precision in base 2 */
300 : static GEN
301 63485526 : mygprec(GEN x, long bit)
302 : {
303 : long e, prec;
304 63485526 : if (bit < 0) bit = 0; /* should rarely happen */
305 63485526 : e = gexpo(x) - bit;
306 63487059 : prec = nbits2prec(bit);
307 63492779 : switch(typ(x))
308 : {
309 167434876 : case t_POL: pari_APPLY_pol_normalized(mygprecrc(gel(x,i),prec,e));
310 17169264 : default: return mygprecrc(x,prec,e);
311 : }
312 : }
313 :
314 : /* normalize a polynomial x, that is change it with coefficients in Z[i],
315 : after making product by 2^shift */
316 : static GEN
317 17543066 : pol_to_gaussint(GEN x, long shift)
318 100062734 : { pari_APPLY_pol_normalized(gtrunc2n(gel(x,i), shift)); }
319 :
320 : /* returns a polynomial q in Z[i][x] keeping bit bits of p */
321 : static GEN
322 13102292 : eval_rel_pol(GEN p, long bit)
323 : {
324 : long i;
325 73802349 : for (i = 2; i < lg(p); i++)
326 60700091 : if (gequal0(gel(p,i))) gel(p,i) = gen_0; /* bad behavior of gexpo */
327 13102258 : return pol_to_gaussint(p, bit-gexpo(p)+1);
328 : }
329 :
330 : /* returns p(R*x)/R^n (in R or R[i]), R = exp(lrho), bit bits of precision */
331 : static GEN
332 1605621 : homothetie(GEN p, double lrho, long bit)
333 : {
334 : GEN q, r, t, iR;
335 1605621 : long n = degpol(p), i;
336 :
337 1605617 : iR = mygprec(dblexp(-lrho),bit);
338 1605590 : q = mygprec(p, bit);
339 1605619 : r = cgetg(n+3,t_POL); r[1] = p[1];
340 1605615 : t = iR; r[n+2] = q[n+2];
341 6774153 : for (i=n-1; i>0; i--)
342 : {
343 5168568 : gel(r,i+2) = gmul(t, gel(q,i+2));
344 5168482 : t = mulrr(t, iR);
345 : }
346 1605585 : gel(r,2) = gmul(t, gel(q,2)); return r;
347 : }
348 :
349 : /* change q in 2^(n*e) p(x*2^(-e)), n=deg(q) [ ~as above with R = 2^-e ]*/
350 : static void
351 9929751 : homothetie2n(GEN p, long e)
352 : {
353 9929751 : if (e)
354 : {
355 7952872 : long i,n = lg(p)-1;
356 43358400 : for (i=2; i<=n; i++) myshiftrc(gel(p,i), (n-i)*e);
357 : }
358 9929869 : }
359 :
360 : /* return 2^f * 2^(n*e) p(x*2^(-e)), n=deg(q) */
361 : static void
362 36413272 : homothetie_gauss(GEN p, long e, long f)
363 : {
364 36413272 : if (e || f)
365 : {
366 32528382 : long i, n = lg(p)-1;
367 167935576 : for (i=2; i<=n; i++) gel(p,i) = myshiftic(gel(p,i), f+(n-i)*e);
368 : }
369 36359077 : }
370 :
371 : /* Lower bound on the modulus of the largest root z_0
372 : * k is set to an upper bound for #{z roots, |z-z_0| < eps} */
373 : static double
374 40862406 : lower_bound(GEN p, long *k, double eps)
375 : {
376 40862406 : long n = degpol(p), i, j;
377 40862907 : pari_sp ltop = avma;
378 : GEN a, s, S, ilc;
379 : double r, R, rho;
380 :
381 40862907 : if (n < 4) { *k = n; return 0.; }
382 8207167 : S = cgetg(5,t_VEC);
383 8208006 : a = cgetg(5,t_VEC); ilc = gdiv(real_1(DEFAULTPREC), gel(p,n+2));
384 41016165 : for (i=1; i<=4; i++) gel(a,i) = gmul(ilc,gel(p,n+2-i));
385 : /* i = 1 split out from next loop for efficiency and initialization */
386 8206668 : s = gel(a,1);
387 8206668 : gel(S,1) = gneg(s); /* Newton sum S_i */
388 8207363 : rho = r = gtodouble(gabs(s,3));
389 8207964 : R = r / n;
390 32824988 : for (i=2; i<=4; i++)
391 : {
392 24617511 : s = gmulsg(i,gel(a,i));
393 73769535 : for (j=1; j<i; j++) s = gadd(s, gmul(gel(S,j),gel(a,i-j)));
394 24600854 : gel(S,i) = gneg(s); /* Newton sum S_i */
395 24611970 : r = gtodouble(gabs(s,3));
396 24617024 : if (r > 0.)
397 : {
398 24571255 : r = exp(log(r/n) / (double)i);
399 24571255 : if (r > R) R = r;
400 : }
401 : }
402 8207477 : if (R > 0. && eps < 1.2)
403 8203740 : *k = (long)floor((rho/R + n) / (1 + exp(-eps)*cos(eps)));
404 : else
405 3737 : *k = n;
406 8207477 : return gc_double(ltop, R);
407 : }
408 :
409 : /* return R such that exp(R - tau) <= rho_n(P) <= exp(R + tau)
410 : * P(0) != 0 and P non constant */
411 : static double
412 4440807 : logmax_modulus(GEN p, double tau)
413 : {
414 : GEN r, q, aux, gunr;
415 4440807 : pari_sp av, ltop = avma;
416 4440807 : long i,k,n=degpol(p),nn,bit,M,e;
417 4440804 : double rho,eps, tau2 = (tau > 3.0)? 0.5: tau/6.;
418 :
419 4440804 : r = cgeti(BIGDEFAULTPREC);
420 4440786 : av = avma;
421 :
422 4440786 : eps = - 1/log(1.5*tau2); /* > 0 */
423 4440786 : bit = (long) ((double) n*log2(1./tau2)+3*log2((double) n))+1;
424 4440786 : gunr = real_1_bit(bit+2*n);
425 4440768 : aux = gdiv(gunr, gel(p,2+n));
426 4440762 : q = RgX_Rg_mul(p, aux); gel(q,2+n) = gunr;
427 4440610 : e = findpower(q);
428 4440750 : homothetie2n(q,e);
429 4440794 : affsi(e, r);
430 4440796 : q = pol_to_gaussint(q, bit);
431 4440349 : M = (long) (log2( log(4.*n) / (2*tau2) )) + 2;
432 4440349 : nn = n;
433 4440349 : for (i=0,e=0;;)
434 : { /* nn = deg(q) */
435 40865408 : rho = lower_bound(q, &k, eps);
436 40862679 : if (rho > exp2(-(double)e)) e = (long)-floor(log2(rho));
437 40862679 : affii(shifti(addis(r,e), 1), r);
438 40818456 : if (++i == M) break;
439 :
440 36378435 : bit = (long) ((double)k * log2(1./tau2) +
441 36378435 : (double)(nn-k)*log2(1./eps) + 3*log2((double)nn)) + 1;
442 36378435 : homothetie_gauss(q, e, bit-(long)floor(dbllog2(gel(q,2+nn))+0.5));
443 36394338 : nn -= RgX_valrem(q, &q);
444 36401297 : q = gerepileupto(av, graeffe(q));
445 36423444 : tau2 *= 1.5; if (tau2 > 0.9) tau2 = 0.5;
446 36423444 : eps = -1/log(tau2); /* > 0 */
447 36423444 : e = findpower(q);
448 : }
449 4440021 : if (!signe(r)) return gc_double(ltop,0.);
450 4019297 : r = itor(r, DEFAULTPREC); shiftr_inplace(r, -M);
451 4019765 : return gc_double(ltop, -rtodbl(r) * M_LN2); /* -log(2) sum e_i 2^-i */
452 : }
453 :
454 : static GEN
455 35451 : RgX_normalize1(GEN x)
456 : {
457 35451 : long i, n = lg(x)-1;
458 : GEN y;
459 35465 : for (i = n; i > 1; i--)
460 35458 : if (!gequal0( gel(x,i) )) break;
461 35451 : if (i == n) return x;
462 14 : pari_warn(warner,"normalizing a polynomial with 0 leading term");
463 14 : if (i == 1) pari_err_ROOTS0("roots");
464 14 : y = cgetg(i+1, t_POL); y[1] = x[1];
465 42 : for (; i > 1; i--) gel(y,i) = gel(x,i);
466 14 : return y;
467 : }
468 :
469 : static GEN
470 27195 : polrootsbound_i(GEN P, double TAU)
471 : {
472 27195 : pari_sp av = avma;
473 : double d;
474 27195 : (void)RgX_valrem_inexact(P,&P);
475 27195 : P = RgX_normalize1(P);
476 27195 : switch(degpol(P))
477 : {
478 7 : case -1: pari_err_ROOTS0("roots");
479 126 : case 0: set_avma(av); return gen_0;
480 : }
481 27062 : d = logmax_modulus(P, TAU) + TAU;
482 : /* not dblexp: result differs on ARM emulator */
483 27062 : return gerepileuptoleaf(av, mpexp(dbltor(d)));
484 : }
485 : GEN
486 27202 : polrootsbound(GEN P, GEN tau)
487 : {
488 27202 : if (typ(P) != t_POL) pari_err_TYPE("polrootsbound",P);
489 27195 : checkvalidpol(P, "polrootsbound");
490 27195 : return polrootsbound_i(P, tau? gtodouble(tau): 0.01);
491 : }
492 :
493 : /* log of modulus of the smallest root of p, with relative error tau */
494 : static double
495 1611321 : logmin_modulus(GEN p, double tau)
496 : {
497 1611321 : pari_sp av = avma;
498 1611321 : if (gequal0(gel(p,2))) return -pariINFINITY;
499 1611317 : return gc_double(av, - logmax_modulus(RgX_recip_i(p),tau));
500 : }
501 :
502 : /* return the log of the k-th modulus (ascending order) of p, rel. error tau*/
503 : static double
504 604811 : logmodulus(GEN p, long k, double tau)
505 : {
506 : GEN q;
507 604811 : long i, kk = k, imax, n = degpol(p), nn, bit, e;
508 604811 : pari_sp av, ltop=avma;
509 604811 : double r, tau2 = tau/6;
510 :
511 604811 : bit = (long)(n * (2. + log2(3.*n/tau2)));
512 604811 : av = avma;
513 604811 : q = gprec_w(p, nbits2prec(bit));
514 604812 : q = RgX_gtofp_bit(q, bit);
515 604811 : e = newton_polygon(q,k);
516 604811 : r = (double)e;
517 604811 : homothetie2n(q,e);
518 604817 : imax = (long)(log2(3./tau) + log2(log(4.*n)))+1;
519 5489272 : for (i=1; i<imax; i++)
520 : {
521 4884466 : q = eval_rel_pol(q,bit);
522 4884152 : kk -= RgX_valrem(q, &q);
523 4884304 : nn = degpol(q);
524 :
525 4884302 : q = gerepileupto(av, graeffe(q));
526 4884458 : e = newton_polygon(q,kk);
527 4884510 : r += e / exp2((double)i);
528 4884510 : q = RgX_gtofp_bit(q, bit);
529 4884400 : homothetie2n(q,e);
530 :
531 4884455 : tau2 *= 1.5; if (tau2 > 1.) tau2 = 1.;
532 4884455 : bit = 1 + (long)(nn*(2. + log2(3.*nn/tau2)));
533 : }
534 604806 : return gc_double(ltop, -r * M_LN2);
535 : }
536 :
537 : /* return the log of the k-th modulus r_k of p, rel. error tau, knowing that
538 : * rmin < r_k < rmax. This information helps because we may reduce precision
539 : * quicker */
540 : static double
541 604806 : logpre_modulus(GEN p, long k, double tau, double lrmin, double lrmax)
542 : {
543 : GEN q;
544 604806 : long n = degpol(p), i, imax, imax2, bit;
545 604805 : pari_sp ltop = avma, av;
546 604805 : double lrho, aux, tau2 = tau/6.;
547 :
548 604805 : aux = (lrmax - lrmin) / 2. + 4*tau2;
549 604805 : imax = (long) log2(log((double)n)/ aux);
550 604805 : if (imax <= 0) return logmodulus(p,k,tau);
551 :
552 595684 : lrho = (lrmin + lrmax) / 2;
553 595684 : av = avma;
554 595684 : bit = (long)(n*(2. + aux / M_LN2 - log2(tau2)));
555 595684 : q = homothetie(p, lrho, bit);
556 595679 : imax2 = (long)(log2(3./tau * log(4.*n))) + 1;
557 595679 : if (imax > imax2) imax = imax2;
558 :
559 1580219 : for (i=0; i<imax; i++)
560 : {
561 984529 : q = eval_rel_pol(q,bit);
562 984521 : q = gerepileupto(av, graeffe(q));
563 984536 : aux = 2*aux + 2*tau2;
564 984536 : tau2 *= 1.5;
565 984536 : bit = (long)(n*(2. + aux / M_LN2 - log2(1-exp(-tau2))));
566 984536 : q = RgX_gtofp_bit(q, bit);
567 : }
568 595690 : aux = exp2((double)imax);
569 595690 : return gc_double(ltop, lrho + logmodulus(q,k, aux*tau/3.) / aux);
570 : }
571 :
572 : static double
573 900859 : ind_maxlog2(GEN q)
574 : {
575 900859 : long i, k = -1;
576 900859 : double L = - pariINFINITY;
577 2218523 : for (i=0; i<=degpol(q); i++)
578 : {
579 1317654 : double d = dbllog2(gel(q,2+i));
580 1317664 : if (d > L) { L = d; k = i; }
581 : }
582 900862 : return k;
583 : }
584 :
585 : /* Returns k such that r_k e^(-tau) < R < r_{k+1} e^tau.
586 : * Assume that l <= k <= n-l */
587 : static long
588 1009939 : dual_modulus(GEN p, double lrho, double tau, long l)
589 : {
590 1009939 : long i, imax, delta_k = 0, n = degpol(p), nn, v2, v, bit, ll = l;
591 1009938 : double tau2 = tau * 7./8.;
592 1009938 : pari_sp av = avma;
593 : GEN q;
594 :
595 1009938 : bit = 6*n - 5*l + (long)(n*(-log2(tau2) + tau2 * 8./7.));
596 1009938 : q = homothetie(p, lrho, bit);
597 1009933 : imax = (long)(log(log(2.*n)/tau2)/log(7./4.)+1);
598 :
599 8134479 : for (i=0; i<imax; i++)
600 : {
601 7233620 : q = eval_rel_pol(q,bit); v2 = n - degpol(q);
602 7232399 : v = RgX_valrem(q, &q);
603 7233050 : ll -= maxss(v, v2); if (ll < 0) ll = 0;
604 :
605 7233187 : nn = degpol(q); delta_k += v;
606 7233186 : if (!nn) return delta_k;
607 :
608 7124107 : q = gerepileupto(av, graeffe(q));
609 7124546 : tau2 *= 7./4.;
610 7124546 : bit = 6*nn - 5*ll + (long)(nn*(-log2(tau2) + tau2 * 8./7.));
611 : }
612 900859 : return gc_long(av, delta_k + (long)ind_maxlog2(q));
613 : }
614 :
615 : /********************************************************************/
616 : /** **/
617 : /** FACTORS THROUGH CIRCLE INTEGRATION **/
618 : /** **/
619 : /********************************************************************/
620 : /* l power of 2, W[step*j] = w_j; set f[j] = p(w_j)
621 : * if inv, w_j = exp(2IPi*j/l), else exp(-2IPi*j/l) */
622 :
623 : static void
624 7462 : fft2(GEN W, GEN p, GEN f, long step, long l)
625 : {
626 : pari_sp av;
627 : long i, s1, l1, step2;
628 :
629 7462 : if (l == 2)
630 : {
631 3766 : gel(f,0) = gadd(gel(p,0), gel(p,step));
632 3766 : gel(f,1) = gsub(gel(p,0), gel(p,step)); return;
633 : }
634 3696 : av = avma;
635 3696 : l1 = l>>1; step2 = step<<1;
636 3696 : fft2(W,p, f, step2,l1);
637 3696 : fft2(W,p+step, f+l1,step2,l1);
638 32760 : for (i = s1 = 0; i < l1; i++, s1 += step)
639 : {
640 29064 : GEN f0 = gel(f,i);
641 29064 : GEN f1 = gmul(gel(W,s1), gel(f,i+l1));
642 29064 : gel(f,i) = gadd(f0, f1);
643 29064 : gel(f,i+l1) = gsub(f0, f1);
644 : }
645 3696 : gerepilecoeffs(av, f, l);
646 : }
647 :
648 : static void
649 14114560 : fft(GEN W, GEN p, GEN f, long step, long l, long inv)
650 : {
651 : pari_sp av;
652 : long i, s1, l1, l2, l3, step4;
653 : GEN f1, f2, f3, f02;
654 :
655 14114560 : if (l == 2)
656 : {
657 6612808 : gel(f,0) = gadd(gel(p,0), gel(p,step));
658 6612458 : gel(f,1) = gsub(gel(p,0), gel(p,step)); return;
659 : }
660 7501752 : av = avma;
661 7501752 : if (l == 4)
662 : {
663 : pari_sp av2;
664 5308392 : f1 = gadd(gel(p,0), gel(p,step<<1));
665 5307848 : f2 = gsub(gel(p,0), gel(p,step<<1));
666 5307908 : f3 = gadd(gel(p,step), gel(p,3*step));
667 5307870 : f02 = gsub(gel(p,step), gel(p,3*step));
668 5307845 : f02 = inv? mulcxI(f02): mulcxmI(f02);
669 5308264 : av2 = avma;
670 5308264 : gel(f,0) = gadd(f1, f3);
671 5307674 : gel(f,1) = gadd(f2, f02);
672 5307860 : gel(f,2) = gsub(f1, f3);
673 5307733 : gel(f,3) = gsub(f2, f02);
674 5307971 : gerepileallsp(av,av2,4,&gel(f,0),&gel(f,1),&gel(f,2),&gel(f,3));
675 5308536 : return;
676 : }
677 2193360 : l1 = l>>2; l2 = 2*l1; l3 = l1+l2; step4 = step<<2;
678 2193360 : fft(W,p, f, step4,l1,inv);
679 2193813 : fft(W,p+step, f+l1,step4,l1,inv);
680 2193814 : fft(W,p+(step<<1),f+l2,step4,l1,inv);
681 2193826 : fft(W,p+3*step, f+l3,step4,l1,inv);
682 8171763 : for (i = s1 = 0; i < l1; i++, s1 += step)
683 : {
684 5978009 : long s2 = s1 << 1, s3 = s1 + s2;
685 : GEN g02, g13, f13;
686 5978009 : f1 = gmul(gel(W,s1), gel(f,i+l1));
687 5978178 : f2 = gmul(gel(W,s2), gel(f,i+l2));
688 5978052 : f3 = gmul(gel(W,s3), gel(f,i+l3));
689 :
690 5978110 : f02 = gadd(gel(f,i),f2);
691 5977624 : g02 = gsub(gel(f,i),f2);
692 5977700 : f13 = gadd(f1,f3);
693 5977619 : g13 = gsub(f1,f3); g13 = inv? mulcxI(g13): mulcxmI(g13);
694 :
695 5978044 : gel(f,i) = gadd(f02, f13);
696 5977670 : gel(f,i+l1) = gadd(g02, g13);
697 5977694 : gel(f,i+l2) = gsub(f02, f13);
698 5977693 : gel(f,i+l3) = gsub(g02, g13);
699 : }
700 2193754 : gerepilecoeffs(av, f, l);
701 : }
702 :
703 : #define code(t1,t2) ((t1 << 6) | t2)
704 :
705 : static GEN
706 98 : FFT_i(GEN W, GEN x)
707 : {
708 98 : long i, l = lg(W), n = lg(x), tx = typ(x), tw, pa;
709 : GEN y, z, p, pol;
710 98 : if (l==1 || ((l-1) & (l-2))) pari_err_DIM("fft");
711 84 : tw = RgV_type(W, &p, &pol, &pa);
712 84 : if (tx == t_POL) { x++; n--; }
713 49 : else if (!is_vec_t(tx)) pari_err_TYPE("fft",x);
714 84 : if (n > l) pari_err_DIM("fft");
715 84 : if (n < l) {
716 0 : z = cgetg(l, t_VECSMALL); /* cf stackdummy */
717 0 : for (i = 1; i < n; i++) gel(z,i) = gel(x,i);
718 0 : for ( ; i < l; i++) gel(z,i) = gen_0;
719 : }
720 84 : else z = x;
721 84 : if (l == 2) return mkveccopy(gel(z,1));
722 70 : y = cgetg(l, t_VEC);
723 70 : if (tw==code(t_COMPLEX,t_INT) || tw==code(t_COMPLEX,t_REAL))
724 0 : {
725 0 : long inv = (l >= 5 && signe(imag_i(gel(W,1+(l>>2))))==1) ? 1 : 0;
726 0 : fft(W+1, z+1, y+1, 1, l-1, inv);
727 : } else
728 70 : fft2(W+1, z+1, y+1, 1, l-1);
729 70 : return y;
730 : }
731 :
732 : #undef code
733 :
734 : GEN
735 56 : FFT(GEN W, GEN x)
736 : {
737 56 : if (!is_vec_t(typ(W))) pari_err_TYPE("fft",W);
738 56 : return FFT_i(W, x);
739 : }
740 :
741 : GEN
742 56 : FFTinv(GEN W, GEN x)
743 : {
744 56 : long l = lg(W), i;
745 : GEN w;
746 56 : if (!is_vec_t(typ(W))) pari_err_TYPE("fft",W);
747 56 : if (l==1 || ((l-1) & (l-2))) pari_err_DIM("fft");
748 42 : w = cgetg(l, t_VECSMALL); /* cf stackdummy */
749 42 : gel(w,1) = gel(W,1); /* w = gconj(W), faster */
750 3773 : for (i = 2; i < l; i++) gel(w, i) = gel(W, l-i+1);
751 42 : return FFT_i(w, x);
752 : }
753 :
754 : /* returns 1 if p has only real coefficients, 0 else */
755 : static int
756 959491 : isreal(GEN p)
757 : {
758 : long i;
759 4848273 : for (i = lg(p)-1; i > 1; i--)
760 4049400 : if (typ(gel(p,i)) == t_COMPLEX) return 0;
761 798873 : return 1;
762 : }
763 :
764 : /* x non complex */
765 : static GEN
766 776173 : abs_update_r(GEN x, double *mu) {
767 776173 : GEN y = gtofp(x, DEFAULTPREC);
768 776175 : double ly = mydbllogr(y); if (ly < *mu) *mu = ly;
769 776175 : setabssign(y); return y;
770 : }
771 : /* return |x|, low accuracy. Set *mu = min(log(y), *mu) */
772 : static GEN
773 7988759 : abs_update(GEN x, double *mu) {
774 : GEN y, xr, yr;
775 : double ly;
776 7988759 : if (typ(x) != t_COMPLEX) return abs_update_r(x, mu);
777 7225116 : xr = gel(x,1);
778 7225116 : yr = gel(x,2);
779 7225116 : if (gequal0(xr)) return abs_update_r(yr,mu);
780 7223183 : if (gequal0(yr)) return abs_update_r(xr,mu);
781 : /* have to treat 0 specially: 0E-10 + 1e-20 = 0E-10 */
782 7212793 : xr = gtofp(xr, DEFAULTPREC);
783 7213615 : yr = gtofp(yr, DEFAULTPREC);
784 7213679 : y = sqrtr(addrr(sqrr(xr), sqrr(yr)));
785 7213518 : ly = mydbllogr(y); if (ly < *mu) *mu = ly;
786 7213753 : return y;
787 : }
788 :
789 : static void
790 993578 : initdft(GEN *Omega, GEN *prim, long N, long Lmax, long bit)
791 : {
792 993578 : long prec = nbits2prec(bit);
793 993576 : *Omega = grootsof1(Lmax, prec) + 1;
794 993578 : *prim = rootsof1u_cx(N, prec);
795 993577 : }
796 :
797 : static void
798 492669 : parameters(GEN p, long *LMAX, double *mu, double *gamma,
799 : int polreal, double param, double param2)
800 : {
801 : GEN q, pc, Omega, A, RU, prim, g, TWO;
802 492669 : long n = degpol(p), bit, NN, K, i, j, Lmax;
803 492669 : pari_sp av2, av = avma;
804 :
805 492669 : bit = gexpo(p) + (long)param2+8;
806 681640 : Lmax = 4; while (Lmax <= n) Lmax <<= 1;
807 492670 : NN = (long)(param*3.14)+1; if (NN < Lmax) NN = Lmax;
808 492670 : K = NN/Lmax; if (K & 1) K++;
809 492670 : NN = Lmax*K;
810 492670 : if (polreal) K = K/2+1;
811 :
812 492670 : initdft(&Omega, &prim, NN, Lmax, bit);
813 492672 : q = mygprec(p,bit) + 2;
814 492671 : A = cgetg(Lmax+1,t_VEC); A++;
815 492671 : pc= cgetg(Lmax+1,t_VEC); pc++;
816 2948568 : for (i=0; i <= n; i++) gel(pc,i)= gel(q,i);
817 964985 : for ( ; i<Lmax; i++) gel(pc,i) = gen_0;
818 :
819 492671 : *mu = pariINFINITY;
820 492671 : g = real_0_bit(-bit);
821 492671 : TWO = real2n(1, DEFAULTPREC);
822 492674 : av2 = avma;
823 492674 : RU = gen_1;
824 1733912 : for (i=0; i<K; i++)
825 : {
826 1241242 : if (i) {
827 748577 : GEN z = RU;
828 3437330 : for (j=1; j<n; j++)
829 : {
830 2688750 : gel(pc,j) = gmul(gel(q,j),z);
831 2688742 : z = gmul(z,RU); /* RU = prim^i, z=prim^(ij) */
832 : }
833 748580 : gel(pc,n) = gmul(gel(q,n),z);
834 : }
835 :
836 1241243 : fft(Omega,pc,A,1,Lmax,1);
837 1241257 : if (polreal && i>0 && i<K-1)
838 1137500 : for (j=0; j<Lmax; j++) g = addrr(g, divrr(TWO, abs_update(gel(A,j),mu)));
839 : else
840 8092434 : for (j=0; j<Lmax; j++) g = addrr(g, invr(abs_update(gel(A,j),mu)));
841 1240976 : RU = gmul(RU, prim);
842 1241239 : if (gc_needed(av,1))
843 : {
844 0 : if(DEBUGMEM>1) pari_warn(warnmem,"parameters");
845 0 : gerepileall(av2,2, &g,&RU);
846 : }
847 : }
848 492670 : *gamma = mydbllog2r(divru(g,NN));
849 492666 : *LMAX = Lmax; set_avma(av);
850 492664 : }
851 :
852 : /* NN is a multiple of Lmax */
853 : static void
854 500907 : dft(GEN p, long k, long NN, long Lmax, long bit, GEN F, GEN H, long polreal)
855 : {
856 : GEN Omega, q, qd, pc, pd, A, B, C, RU, aux, U, W, prim, prim2;
857 500907 : long n = degpol(p), i, j, K;
858 : pari_sp ltop;
859 :
860 500907 : initdft(&Omega, &prim, NN, Lmax, bit);
861 500906 : RU = cgetg(n+2,t_VEC) + 1;
862 :
863 500907 : K = NN/Lmax; if (polreal) K = K/2+1;
864 500907 : q = mygprec(p,bit);
865 500907 : qd = RgX_deriv(q);
866 :
867 500905 : A = cgetg(Lmax+1,t_VEC); A++;
868 500905 : B = cgetg(Lmax+1,t_VEC); B++;
869 500905 : C = cgetg(Lmax+1,t_VEC); C++;
870 500905 : pc = cgetg(Lmax+1,t_VEC); pc++;
871 500907 : pd = cgetg(Lmax+1,t_VEC); pd++;
872 1015264 : gel(pc,0) = gel(q,2); for (i=n+1; i<Lmax; i++) gel(pc,i) = gen_0;
873 1516171 : gel(pd,0) = gel(qd,2); for (i=n; i<Lmax; i++) gel(pd,i) = gen_0;
874 :
875 500907 : ltop = avma;
876 500907 : W = cgetg(k+1,t_VEC);
877 500907 : U = cgetg(k+1,t_VEC);
878 1197910 : for (i=1; i<=k; i++) gel(W,i) = gel(U,i) = gen_0;
879 :
880 500907 : gel(RU,0) = gen_1;
881 500907 : prim2 = gen_1;
882 1525687 : for (i=0; i<K; i++)
883 : {
884 1024778 : gel(RU,1) = prim2;
885 4432365 : for (j=1; j<n; j++) gel(RU,j+1) = gmul(gel(RU,j),prim2);
886 : /* RU[j] = prim^(ij)= prim2^j */
887 :
888 4432314 : for (j=1; j<n; j++) gel(pd,j) = gmul(gel(qd,j+2),gel(RU,j));
889 1024732 : fft(Omega,pd,A,1,Lmax,1);
890 5457038 : for (j=1; j<=n; j++) gel(pc,j) = gmul(gel(q,j+2),gel(RU,j));
891 1024721 : fft(Omega,pc,B,1,Lmax,1);
892 7642062 : for (j=0; j<Lmax; j++) gel(C,j) = ginv(gel(B,j));
893 7642022 : for (j=0; j<Lmax; j++) gel(B,j) = gmul(gel(A,j),gel(C,j));
894 1024648 : fft(Omega,B,A,1,Lmax,1);
895 1024777 : fft(Omega,C,B,1,Lmax,1);
896 :
897 1024781 : if (polreal) /* p has real coefficients */
898 : {
899 794954 : if (i>0 && i<K-1)
900 : {
901 102716 : for (j=1; j<=k; j++)
902 : {
903 86136 : gel(W,j) = gadd(gel(W,j), gshift(mulreal(gel(A,j+1),gel(RU,j+1)),1));
904 86136 : gel(U,j) = gadd(gel(U,j), gshift(mulreal(gel(B,j),gel(RU,j)),1));
905 : }
906 : }
907 : else
908 : {
909 1824912 : for (j=1; j<=k; j++)
910 : {
911 1046556 : gel(W,j) = gadd(gel(W,j), mulreal(gel(A,j+1),gel(RU,j+1)));
912 1046534 : gel(U,j) = gadd(gel(U,j), mulreal(gel(B,j),gel(RU,j)));
913 : }
914 : }
915 : }
916 : else
917 : {
918 601922 : for (j=1; j<=k; j++)
919 : {
920 372099 : gel(W,j) = gadd(gel(W,j), gmul(gel(A,j+1),gel(RU,j+1)));
921 372095 : gel(U,j) = gadd(gel(U,j), gmul(gel(B,j),gel(RU,j)));
922 : }
923 : }
924 1024759 : prim2 = gmul(prim2,prim);
925 1024768 : gerepileall(ltop,3, &W,&U,&prim2);
926 : }
927 :
928 1197905 : for (i=1; i<=k; i++)
929 : {
930 696999 : aux=gel(W,i);
931 1093591 : for (j=1; j<i; j++) aux = gadd(aux, gmul(gel(W,i-j),gel(F,k+2-j)));
932 697001 : gel(F,k+2-i) = gdivgs(aux,-i*NN);
933 : }
934 1197900 : for (i=0; i<k; i++)
935 : {
936 696996 : aux=gel(U,k-i);
937 1093583 : for (j=1+i; j<k; j++) aux = gadd(aux,gmul(gel(F,2+j),gel(U,j-i)));
938 696992 : gel(H,i+2) = gdivgu(aux,NN);
939 : }
940 500904 : }
941 :
942 : #define NEWTON_MAX 10
943 : static GEN
944 2452199 : refine_H(GEN F, GEN G, GEN HH, long bit, long Sbit)
945 : {
946 2452199 : GEN H = HH, D, aux;
947 2452199 : pari_sp ltop = avma;
948 : long error, i, bit1, bit2;
949 :
950 2452199 : D = Rg_RgX_sub(gen_1, RgX_rem(RgX_mul(H,G),F)); error = gexpo(D);
951 2452170 : bit2 = bit + Sbit;
952 4481886 : for (i=0; error>-bit && i<NEWTON_MAX && error<=0; i++)
953 : {
954 2029730 : if (gc_needed(ltop,1))
955 : {
956 0 : if(DEBUGMEM>1) pari_warn(warnmem,"refine_H");
957 0 : gerepileall(ltop,2, &D,&H);
958 : }
959 2029730 : bit1 = -error + Sbit;
960 2029730 : aux = RgX_mul(mygprec(H,bit1), mygprec(D,bit1));
961 2029730 : aux = RgX_rem(mygprec(aux,bit1), mygprec(F,bit1));
962 :
963 2029745 : bit1 = -error*2 + Sbit; if (bit1 > bit2) bit1 = bit2;
964 2029745 : H = RgX_add(mygprec(H,bit1), aux);
965 2029691 : D = Rg_RgX_sub(gen_1, RgX_rem(RgX_mul(H,G),F));
966 2029716 : error = gexpo(D); if (error < -bit1) error = -bit1;
967 : }
968 2452156 : if (error > -bit/2) return NULL; /* FAIL */
969 2451824 : return gerepilecopy(ltop,H);
970 : }
971 :
972 : /* return 0 if fails, 1 else */
973 : static long
974 500901 : refine_F(GEN p, GEN *F, GEN *G, GEN H, long bit, double gamma)
975 : {
976 500901 : GEN f0, FF, GG, r, HH = H;
977 500901 : long error, i, bit1 = 0, bit2, Sbit, Sbit2, enh, normF, normG, n = degpol(p);
978 500901 : pari_sp av = avma;
979 :
980 500901 : FF = *F; GG = RgX_divrem(p, FF, &r);
981 500909 : error = gexpo(r); if (error <= -bit) error = 1-bit;
982 500909 : normF = gexpo(FF);
983 500908 : normG = gexpo(GG);
984 500909 : enh = gexpo(H); if (enh < 0) enh = 0;
985 500909 : Sbit = normF + 2*normG + enh + (long)(4.*log2((double)n)+gamma) + 1;
986 500909 : Sbit2 = enh + 2*(normF+normG) + (long)(2.*gamma+5.*log2((double)n)) + 1;
987 500909 : bit2 = bit + Sbit;
988 2952745 : for (i=0; error>-bit && i<NEWTON_MAX && error<=0; i++)
989 : {
990 2452184 : if (bit1 == bit2 && i >= 2) { Sbit += n; Sbit2 += n; bit2 += n; }
991 2452184 : if (gc_needed(av,1))
992 : {
993 0 : if(DEBUGMEM>1) pari_warn(warnmem,"refine_F");
994 0 : gerepileall(av,4, &FF,&GG,&r,&HH);
995 : }
996 :
997 2452184 : bit1 = -error + Sbit2;
998 2452184 : HH = refine_H(mygprec(FF,bit1), mygprec(GG,bit1), mygprec(HH,bit1),
999 : 1-error, Sbit2);
1000 2452201 : if (!HH) return 0; /* FAIL */
1001 :
1002 2451869 : bit1 = -error + Sbit;
1003 2451869 : r = RgX_mul(mygprec(HH,bit1), mygprec(r,bit1));
1004 2451832 : f0 = RgX_rem(mygprec(r,bit1), mygprec(FF,bit1));
1005 :
1006 2451847 : bit1 = -2*error + Sbit; if (bit1 > bit2) bit1 = bit2;
1007 2451847 : FF = gadd(mygprec(FF,bit1),f0);
1008 :
1009 2451819 : bit1 = -3*error + Sbit; if (bit1 > bit2) bit1 = bit2;
1010 2451819 : GG = RgX_divrem(mygprec(p,bit1), mygprec(FF,bit1), &r);
1011 2451830 : error = gexpo(r); if (error < -bit1) error = -bit1;
1012 : }
1013 500561 : if (error>-bit) return 0; /* FAIL */
1014 492657 : *F = FF; *G = GG; return 1;
1015 : }
1016 :
1017 : /* returns F and G from the unit circle U such that |p-FG|<2^(-bit) |cd|,
1018 : where cd is the leading coefficient of p */
1019 : static void
1020 492673 : split_fromU(GEN p, long k, double delta, long bit,
1021 : GEN *F, GEN *G, double param, double param2)
1022 : {
1023 : GEN pp, FF, GG, H;
1024 492673 : long n = degpol(p), NN, bit2, Lmax;
1025 492673 : int polreal = isreal(p);
1026 : pari_sp ltop;
1027 : double mu, gamma;
1028 :
1029 492673 : pp = gdiv(p, gel(p,2+n));
1030 492669 : parameters(pp, &Lmax,&mu,&gamma, polreal,param,param2);
1031 :
1032 492664 : H = cgetg(k+2,t_POL); H[1] = p[1];
1033 492666 : FF = cgetg(k+3,t_POL); FF[1]= p[1];
1034 492668 : gel(FF,k+2) = gen_1;
1035 :
1036 492668 : NN = (long)(0.5/delta); NN |= 1; if (NN < 2) NN = 2;
1037 492668 : NN *= Lmax; ltop = avma;
1038 : for(;;)
1039 : {
1040 500904 : bit2 = (long)(((double)NN*delta-mu)/M_LN2) + gexpo(pp) + 8;
1041 500907 : dft(pp, k, NN, Lmax, bit2, FF, H, polreal);
1042 500901 : if (refine_F(pp,&FF,&GG,H,bit,gamma)) break;
1043 8236 : NN <<= 1; set_avma(ltop);
1044 : }
1045 492668 : *G = gmul(GG,gel(p,2+n)); *F = FF;
1046 492670 : }
1047 :
1048 : static void
1049 492673 : optimize_split(GEN p, long k, double delta, long bit,
1050 : GEN *F, GEN *G, double param, double param2)
1051 : {
1052 492673 : long n = degpol(p);
1053 : GEN FF, GG;
1054 :
1055 492673 : if (k <= n/2)
1056 381888 : split_fromU(p,k,delta,bit,F,G,param,param2);
1057 : else
1058 : {
1059 110785 : split_fromU(RgX_recip_i(p),n-k,delta,bit,&FF,&GG,param,param2);
1060 110785 : *F = RgX_recip_i(GG);
1061 110785 : *G = RgX_recip_i(FF);
1062 : }
1063 492670 : }
1064 :
1065 : /********************************************************************/
1066 : /** **/
1067 : /** SEARCH FOR SEPARATING CIRCLE **/
1068 : /** **/
1069 : /********************************************************************/
1070 :
1071 : /* return p(2^e*x) *2^(-n*e) */
1072 : static void
1073 0 : scalepol2n(GEN p, long e)
1074 : {
1075 0 : long i,n=lg(p)-1;
1076 0 : for (i=2; i<=n; i++) gel(p,i) = gmul2n(gel(p,i),(i-n)*e);
1077 0 : }
1078 :
1079 : /* returns p(x/R)*R^n; assume R is at the correct accuracy */
1080 : static GEN
1081 4278836 : scalepol(GEN p, GEN R, long bit)
1082 4278836 : { return RgX_rescale(mygprec(p, bit), R); }
1083 :
1084 : /* return (conj(a)X-1)^n * p[ (X-a) / (conj(a)X-1) ] */
1085 : static GEN
1086 1400571 : conformal_basecase(GEN p, GEN a)
1087 : {
1088 : GEN z, r, ma, ca;
1089 1400571 : long i, n = degpol(p);
1090 : pari_sp av;
1091 :
1092 1400571 : if (n <= 0) return p;
1093 1400571 : ma = gneg(a); ca = conj_i(a);
1094 1400572 : av = avma;
1095 1400572 : z = deg1pol_shallow(ca, gen_m1, 0);
1096 1400571 : r = scalarpol_shallow(gel(p,2+n), 0);
1097 3631778 : for (i=n-1; ; i--)
1098 : {
1099 3631778 : r = RgX_addmulXn_shallow(r, gmul(ma,r), 1); /* r *= (X - a) */
1100 3631735 : r = gadd(r, gmul(z, gel(p,2+i)));
1101 3631743 : if (i == 0) return gerepileupto(av, r);
1102 2231186 : z = RgX_addmulXn_shallow(gmul(z,ca), gneg(z), 1); /* z *= conj(a)X - 1 */
1103 2231211 : if (gc_needed(av,2))
1104 : {
1105 0 : if(DEBUGMEM>1) pari_warn(warnmem,"conformal_pol (%ld/%ld)",n-i, n);
1106 0 : gerepileall(av,2, &r,&z);
1107 : }
1108 : }
1109 : }
1110 : static GEN
1111 1400691 : conformal_pol(GEN p, GEN a)
1112 : {
1113 1400691 : pari_sp av = avma;
1114 1400691 : long d, nR, n = degpol(p), v;
1115 : GEN Q, R, S, T;
1116 1400690 : if (n < 35) return conformal_basecase(p, a);
1117 119 : d = (n+1) >> 1; v = varn(p);
1118 119 : Q = RgX_shift_shallow(p, -d);
1119 119 : R = RgXn_red_shallow(p, d);
1120 119 : Q = conformal_pol(Q, a);
1121 119 : R = conformal_pol(R, a);
1122 119 : S = gpowgs(deg1pol_shallow(gen_1, gneg(a), v), d);
1123 119 : T = RgX_recip_i(S);
1124 119 : if (typ(a) == t_COMPLEX) T = gconj(T);
1125 119 : if (odd(d)) T = RgX_neg(T);
1126 : /* S = (X - a)^d, T = (conj(a) X - 1)^d */
1127 119 : nR = n - degpol(R) - d; /* >= 0 */
1128 119 : if (nR) T = RgX_mul(T, gpowgs(deg1pol_shallow(gconj(a), gen_m1, v), nR));
1129 119 : return gerepileupto(av, RgX_add(RgX_mul(Q, S), RgX_mul(R, T)));
1130 : }
1131 :
1132 : static const double UNDEF = -100000.;
1133 :
1134 : static double
1135 492663 : logradius(double *radii, GEN p, long k, double aux, double *delta)
1136 : {
1137 492663 : long i, n = degpol(p);
1138 : double lrho, lrmin, lrmax;
1139 492664 : if (k > 1)
1140 : {
1141 281275 : i = k-1; while (i>0 && radii[i] == UNDEF) i--;
1142 206404 : lrmin = logpre_modulus(p,k,aux, radii[i], radii[k]);
1143 : }
1144 : else /* k=1 */
1145 286260 : lrmin = logmin_modulus(p,aux);
1146 492669 : radii[k] = lrmin;
1147 :
1148 492669 : if (k+1<n)
1149 : {
1150 589456 : i = k+2; while (i<=n && radii[i] == UNDEF) i++;
1151 398402 : lrmax = logpre_modulus(p,k+1,aux, radii[k+1], radii[i]);
1152 : }
1153 : else /* k+1=n */
1154 94267 : lrmax = logmax_modulus(p,aux);
1155 492668 : radii[k+1] = lrmax;
1156 :
1157 492668 : lrho = radii[k];
1158 821631 : for (i=k-1; i>=1; i--)
1159 : {
1160 328963 : if (radii[i] == UNDEF || radii[i] > lrho)
1161 241131 : radii[i] = lrho;
1162 : else
1163 87832 : lrho = radii[i];
1164 : }
1165 492668 : lrho = radii[k+1];
1166 1634259 : for (i=k+1; i<=n; i++)
1167 : {
1168 1141591 : if (radii[i] == UNDEF || radii[i] < lrho)
1169 566263 : radii[i] = lrho;
1170 : else
1171 575328 : lrho = radii[i];
1172 : }
1173 492668 : *delta = (lrmax - lrmin) / 2;
1174 492668 : if (*delta > 1.) *delta = 1.;
1175 492668 : return (lrmin + lrmax) / 2;
1176 : }
1177 :
1178 : static void
1179 492668 : update_radius(long n, double *radii, double lrho, double *par, double *par2)
1180 : {
1181 492668 : double t, param = 0., param2 = 0.;
1182 : long i;
1183 2455836 : for (i=1; i<=n; i++)
1184 : {
1185 1963184 : radii[i] -= lrho;
1186 1963184 : t = fabs(rtodbl( invr(subsr(1, dblexp(radii[i]))) ));
1187 1963168 : param += t; if (t > 1.) param2 += log2(t);
1188 : }
1189 492652 : *par = param; *par2 = param2;
1190 492652 : }
1191 :
1192 : /* apply the conformal mapping then split from U */
1193 : static void
1194 466815 : conformal_mapping(double *radii, GEN ctr, GEN p, long k, long bit,
1195 : double aux, GEN *F,GEN *G)
1196 : {
1197 466815 : long bit2, n = degpol(p), i;
1198 466815 : pari_sp ltop = avma, av;
1199 : GEN q, FF, GG, a, R;
1200 : double lrho, delta, param, param2;
1201 : /* n * (2.*log2(2.732)+log2(1.5)) + 1 */
1202 466815 : bit2 = bit + (long)(n*3.4848775) + 1;
1203 466815 : a = sqrtr_abs( utor(3, precdbl(MEDDEFAULTPREC)) );
1204 466818 : a = divrs(a, -6);
1205 466816 : a = gmul(mygprec(a,bit2), mygprec(ctr,bit2)); /* a = -ctr/2sqrt(3) */
1206 :
1207 466819 : av = avma;
1208 466819 : q = conformal_pol(mygprec(p,bit2), a);
1209 2282740 : for (i=1; i<=n; i++)
1210 1815930 : if (radii[i] != UNDEF) /* update array radii */
1211 : {
1212 1537326 : pari_sp av2 = avma;
1213 1537326 : GEN t, r = dblexp(radii[i]), r2 = sqrr(r);
1214 : /* 2(r^2 - 1) / (r^2 - 3(r-1)) */
1215 1537287 : t = divrr(shiftr((subrs(r2,1)),1), subrr(r2, mulur(3,subrs(r,1))));
1216 1537327 : radii[i] = mydbllogr(addsr(1,t)) / 2;
1217 1537301 : set_avma(av2);
1218 : }
1219 466810 : lrho = logradius(radii, q,k,aux/10., &delta);
1220 466815 : update_radius(n, radii, lrho, ¶m, ¶m2);
1221 :
1222 466818 : bit2 += (long)(n * fabs(lrho)/M_LN2 + 1.);
1223 466818 : R = mygprec(dblexp(-lrho), bit2);
1224 466816 : q = scalepol(q,R,bit2);
1225 466820 : gerepileall(av,2, &q,&R);
1226 :
1227 466820 : optimize_split(q,k,delta,bit2,&FF,&GG,param,param2);
1228 466818 : bit2 += n; R = invr(R);
1229 466816 : FF = scalepol(FF,R,bit2);
1230 466819 : GG = scalepol(GG,R,bit2);
1231 :
1232 466818 : a = mygprec(a,bit2);
1233 466819 : FF = conformal_pol(FF,a);
1234 466819 : GG = conformal_pol(GG,a);
1235 :
1236 466818 : a = invr(subsr(1, gnorm(a)));
1237 466817 : FF = RgX_Rg_mul(FF, powru(a,k));
1238 466820 : GG = RgX_Rg_mul(GG, powru(a,n-k));
1239 :
1240 466820 : *F = mygprec(FF,bit+n);
1241 466819 : *G = mygprec(GG,bit+n); gerepileall(ltop,2, F,G);
1242 466820 : }
1243 :
1244 : /* split p, this time without scaling. returns in F and G two polynomials
1245 : * such that |p-FG|< 2^(-bit)|p| */
1246 : static void
1247 492671 : split_2(GEN p, long bit, GEN ctr, double thickness, GEN *F, GEN *G)
1248 : {
1249 : GEN q, FF, GG, R;
1250 : double aux, delta, param, param2;
1251 492671 : long n = degpol(p), i, j, k, bit2;
1252 : double lrmin, lrmax, lrho, *radii;
1253 :
1254 492671 : radii = (double*) stack_malloc_align((n+1) * sizeof(double), sizeof(double));
1255 :
1256 1470558 : for (i=2; i<n; i++) radii[i] = UNDEF;
1257 492670 : aux = thickness/(double)(4 * n);
1258 492670 : lrmin = logmin_modulus(p, aux);
1259 492669 : lrmax = logmax_modulus(p, aux);
1260 492668 : radii[1] = lrmin;
1261 492668 : radii[n] = lrmax;
1262 492668 : i = 1; j = n;
1263 492668 : lrho = (lrmin + lrmax) / 2;
1264 492668 : k = dual_modulus(p, lrho, aux, 1);
1265 492667 : if (5*k < n || (n < 2*k && 5*k < 4*n))
1266 77867 : { lrmax = lrho; j=k+1; radii[j] = lrho; }
1267 : else
1268 414800 : { lrmin = lrho; i=k; radii[i] = lrho; }
1269 1009942 : while (j > i+1)
1270 : {
1271 517273 : if (i+j == n+1)
1272 371788 : lrho = (lrmin + lrmax) / 2;
1273 : else
1274 : {
1275 145485 : double kappa = 2. - log(1. + minss(i,n-j)) / log(1. + minss(j,n-i));
1276 145485 : if (i+j < n+1) lrho = lrmax * kappa + lrmin;
1277 116034 : else lrho = lrmin * kappa + lrmax;
1278 145485 : lrho /= 1+kappa;
1279 : }
1280 517273 : aux = (lrmax - lrmin) / (4*(j-i));
1281 517273 : k = dual_modulus(p, lrho, aux, minss(i,n+1-j));
1282 517275 : if (k-i < j-k-1 || (k-i == j-k-1 && 2*k > n))
1283 386061 : { lrmax = lrho; j=k+1; radii[j] = lrho - aux; }
1284 : else
1285 131214 : { lrmin = lrho; i=k; radii[i] = lrho + aux; }
1286 : }
1287 492669 : aux = lrmax - lrmin;
1288 :
1289 492669 : if (ctr)
1290 : {
1291 466820 : lrho = (lrmax + lrmin) / 2;
1292 2282792 : for (i=1; i<=n; i++)
1293 1815972 : if (radii[i] != UNDEF) radii[i] -= lrho;
1294 :
1295 466820 : bit2 = bit + (long)(n * fabs(lrho)/M_LN2 + 1.);
1296 466820 : R = mygprec(dblexp(-lrho), bit2);
1297 466813 : q = scalepol(p,R,bit2);
1298 466817 : conformal_mapping(radii, ctr, q, k, bit2, aux, &FF, &GG);
1299 : }
1300 : else
1301 : {
1302 25849 : lrho = logradius(radii, p, k, aux/10., &delta);
1303 25853 : update_radius(n, radii, lrho, ¶m, ¶m2);
1304 :
1305 25853 : bit2 = bit + (long)(n * fabs(lrho)/M_LN2 + 1.);
1306 25853 : R = mygprec(dblexp(-lrho), bit2);
1307 25853 : q = scalepol(p,R,bit2);
1308 25853 : optimize_split(q, k, delta, bit2, &FF, &GG, param, param2);
1309 : }
1310 492673 : bit += n;
1311 492673 : bit2 += n; R = invr(mygprec(R,bit2));
1312 492672 : *F = mygprec(scalepol(FF,R,bit2), bit);
1313 492672 : *G = mygprec(scalepol(GG,R,bit2), bit);
1314 492670 : }
1315 :
1316 : /* procedure corresponding to steps 5,6,.. page 44 in RR n. 1852 */
1317 : /* put in F and G two polynomial such that |p-FG|<2^(-bit)|p|
1318 : * where the maximum modulus of the roots of p is <=1.
1319 : * Assume sum of roots is 0. */
1320 : static void
1321 466818 : split_1(GEN p, long bit, GEN *F, GEN *G)
1322 : {
1323 466818 : long i, imax, n = degpol(p), polreal = isreal(p), ep = gexpo(p), bit2 = bit+n;
1324 : GEN ctr, q, qq, FF, GG, v, gr, r, newq;
1325 : double lrmin, lrmax, lthick;
1326 466819 : const double LOG3 = 1.098613;
1327 :
1328 466819 : lrmax = logmax_modulus(p, 0.01);
1329 466815 : gr = mygprec(dblexp(-lrmax), bit2);
1330 466813 : q = scalepol(p,gr,bit2);
1331 :
1332 466819 : bit2 = bit + gexpo(q) - ep + (long)((double)n*2.*log2(3.)+1);
1333 466819 : v = cgetg(5,t_VEC);
1334 466819 : gel(v,1) = gen_2;
1335 466819 : gel(v,2) = gen_m2;
1336 466819 : gel(v,3) = mkcomplex(gen_0, gel(v,1));
1337 466820 : gel(v,4) = mkcomplex(gen_0, gel(v,2));
1338 466820 : q = mygprec(q,bit2); lthick = 0;
1339 466820 : newq = ctr = NULL; /* -Wall */
1340 466820 : imax = polreal? 3: 4;
1341 838934 : for (i=1; i<=imax; i++)
1342 : {
1343 832393 : qq = RgX_translate(q, gel(v,i));
1344 832397 : lrmin = logmin_modulus(qq,0.05);
1345 832377 : if (LOG3 > lrmin + lthick)
1346 : {
1347 820185 : double lquo = logmax_modulus(qq,0.05) - lrmin;
1348 820196 : if (lquo > lthick) { lthick = lquo; newq = qq; ctr = gel(v,i); }
1349 : }
1350 832388 : if (lthick > M_LN2) break;
1351 423045 : if (polreal && i==2 && lthick > LOG3 - M_LN2) break;
1352 : }
1353 466815 : bit2 = bit + gexpo(newq) - ep + (long)(n*LOG3/M_LN2 + 1);
1354 466818 : split_2(newq, bit2, ctr, lthick, &FF, &GG);
1355 466817 : r = gneg(mygprec(ctr,bit2));
1356 466818 : FF = RgX_translate(FF,r);
1357 466819 : GG = RgX_translate(GG,r);
1358 :
1359 466820 : gr = invr(gr); bit2 = bit - ep + gexpo(FF)+gexpo(GG);
1360 466820 : *F = scalepol(FF,gr,bit2);
1361 466820 : *G = scalepol(GG,gr,bit2);
1362 466814 : }
1363 :
1364 : /* put in F and G two polynomials such that |P-FG|<2^(-bit)|P|,
1365 : where the maximum modulus of the roots of p is < 0.5 */
1366 : static int
1367 467138 : split_0_2(GEN p, long bit, GEN *F, GEN *G)
1368 : {
1369 : GEN q, b;
1370 467138 : long n = degpol(p), k, bit2, eq;
1371 467138 : double aux0 = dbllog2(gel(p,n+2)); /* != -oo */
1372 467139 : double aux1 = dbllog2(gel(p,n+1)), aux;
1373 :
1374 467141 : if (aux1 == -pariINFINITY) /* p1 = 0 */
1375 9874 : aux = 0;
1376 : else
1377 : {
1378 457267 : aux = aux1 - aux0; /* log2(p1/p0) */
1379 : /* beware double overflow */
1380 457267 : if (aux >= 0 && (aux > 1e4 || exp2(aux) > 2.5*n)) return 0;
1381 457267 : aux = (aux < -300)? 0.: n*log2(1 + exp2(aux)/(double)n);
1382 : }
1383 467141 : bit2 = bit+1 + (long)(log2((double)n) + aux);
1384 467141 : q = mygprec(p,bit2);
1385 467141 : if (aux1 == -pariINFINITY) b = NULL;
1386 : else
1387 : {
1388 457267 : b = gdivgs(gdiv(gel(q,n+1),gel(q,n+2)),-n);
1389 457268 : q = RgX_translate(q,b);
1390 : }
1391 467144 : gel(q,n+1) = gen_0; eq = gexpo(q);
1392 467143 : k = 0;
1393 467697 : while (k <= n/2 && (- gexpo(gel(q,k+2)) > bit2 + 2*(n-k) + eq
1394 467569 : || gequal0(gel(q,k+2)))) k++;
1395 467142 : if (k > 0)
1396 : {
1397 324 : if (k > n/2) k = n/2;
1398 324 : bit2 += k<<1;
1399 324 : *F = pol_xn(k, 0);
1400 324 : *G = RgX_shift_shallow(q, -k);
1401 : }
1402 : else
1403 : {
1404 466818 : split_1(q,bit2,F,G);
1405 466814 : bit2 = bit + gexpo(*F) + gexpo(*G) - gexpo(p) + (long)aux+1;
1406 466819 : *F = mygprec(*F,bit2);
1407 : }
1408 467141 : *G = mygprec(*G,bit2);
1409 467142 : if (b)
1410 : {
1411 457268 : GEN mb = mygprec(gneg(b), bit2);
1412 457267 : *F = RgX_translate(*F, mb);
1413 457270 : *G = RgX_translate(*G, mb);
1414 : }
1415 467143 : return 1;
1416 : }
1417 :
1418 : /* put in F and G two polynomials such that |P-FG|<2^(-bit)|P|.
1419 : * Assume max_modulus(p) < 2 */
1420 : static void
1421 467138 : split_0_1(GEN p, long bit, GEN *F, GEN *G)
1422 : {
1423 : GEN FF, GG;
1424 : long n, bit2, normp;
1425 :
1426 467138 : if (split_0_2(p,bit,F,G)) return;
1427 :
1428 0 : normp = gexpo(p);
1429 0 : scalepol2n(p,2); /* p := 4^(-n) p(4*x) */
1430 0 : n = degpol(p); bit2 = bit + 2*n + gexpo(p) - normp;
1431 0 : split_1(mygprec(p,bit2), bit2,&FF,&GG);
1432 0 : scalepol2n(FF,-2);
1433 0 : scalepol2n(GG,-2); bit2 = bit + gexpo(FF) + gexpo(GG) - normp;
1434 0 : *F = mygprec(FF,bit2);
1435 0 : *G = mygprec(GG,bit2);
1436 : }
1437 :
1438 : /* put in F and G two polynomials such that |P-FG|<2^(-bit)|P| */
1439 : static void
1440 492994 : split_0(GEN p, long bit, GEN *F, GEN *G)
1441 : {
1442 492994 : const double LOG1_9 = 0.6418539;
1443 492994 : long n = degpol(p), k = 0;
1444 : GEN q;
1445 :
1446 492993 : while (gexpo(gel(p,k+2)) < -bit && k <= n/2) k++;
1447 492993 : if (k > 0)
1448 : {
1449 0 : if (k > n/2) k = n/2;
1450 0 : *F = pol_xn(k, 0);
1451 0 : *G = RgX_shift_shallow(p, -k);
1452 : }
1453 : else
1454 : {
1455 492993 : double lr = logmax_modulus(p, 0.05);
1456 492994 : if (lr < LOG1_9) split_0_1(p, bit, F, G);
1457 : else
1458 : {
1459 435568 : q = RgX_recip_i(p);
1460 435569 : lr = logmax_modulus(q,0.05);
1461 435565 : if (lr < LOG1_9)
1462 : {
1463 409712 : split_0_1(q, bit, F, G);
1464 409717 : *F = RgX_recip_i(*F);
1465 409717 : *G = RgX_recip_i(*G);
1466 : }
1467 : else
1468 25853 : split_2(p,bit,NULL, 1.2837,F,G);
1469 : }
1470 : }
1471 492996 : }
1472 :
1473 : /********************************************************************/
1474 : /** **/
1475 : /** ERROR ESTIMATE FOR THE ROOTS **/
1476 : /** **/
1477 : /********************************************************************/
1478 :
1479 : static GEN
1480 1892874 : root_error(long n, long k, GEN roots_pol, long err, GEN shatzle)
1481 : {
1482 1892874 : GEN rho, d, eps, epsbis, eps2, aux, rap = NULL;
1483 : long i, j;
1484 :
1485 1892874 : d = cgetg(n+1,t_VEC);
1486 12062035 : for (i=1; i<=n; i++)
1487 : {
1488 10169311 : if (i!=k)
1489 : {
1490 8276525 : aux = gsub(gel(roots_pol,i), gel(roots_pol,k));
1491 8275928 : gel(d,i) = gabs(mygprec(aux,31), DEFAULTPREC);
1492 : }
1493 : }
1494 1892724 : rho = gabs(mygprec(gel(roots_pol,k),31), DEFAULTPREC);
1495 1892896 : if (expo(rho) < 0) rho = real_1(DEFAULTPREC);
1496 1892896 : eps = mulrr(rho, shatzle);
1497 1892822 : aux = shiftr(powru(rho,n), err);
1498 :
1499 5745279 : for (j=1; j<=2 || (j<=5 && cmprr(rap, dbltor(1.2)) > 0); j++)
1500 : {
1501 3852537 : GEN prod = NULL; /* 1. */
1502 3852537 : long m = n;
1503 3852537 : epsbis = mulrr(eps, dbltor(1.25));
1504 26372375 : for (i=1; i<=n; i++)
1505 : {
1506 22519700 : if (i != k && cmprr(gel(d,i),epsbis) > 0)
1507 : {
1508 18627799 : GEN dif = subrr(gel(d,i),eps);
1509 18625027 : prod = prod? mulrr(prod, dif): dif;
1510 18626980 : m--;
1511 : }
1512 : }
1513 3852675 : eps2 = prod? divrr(aux, prod): aux;
1514 3852581 : if (m > 1) eps2 = sqrtnr(shiftr(eps2, 2*m-2), m);
1515 3852581 : rap = divrr(eps,eps2); eps = eps2;
1516 : }
1517 1892640 : return eps;
1518 : }
1519 :
1520 : /* round a complex or real number x to an absolute value of 2^(-bit) */
1521 : static GEN
1522 4284656 : mygprec_absolute(GEN x, long bit)
1523 : {
1524 : long e;
1525 : GEN y;
1526 :
1527 4284656 : switch(typ(x))
1528 : {
1529 2942188 : case t_REAL:
1530 2942188 : e = expo(x) + bit;
1531 2942188 : return (e <= 0 || !signe(x))? real_0_bit(-bit): rtor(x, nbits2prec(e));
1532 1213016 : case t_COMPLEX:
1533 1213016 : if (gexpo(gel(x,2)) < -bit) return mygprec_absolute(gel(x,1),bit);
1534 1178877 : y = cgetg(3,t_COMPLEX);
1535 1178880 : gel(y,1) = mygprec_absolute(gel(x,1),bit);
1536 1178882 : gel(y,2) = mygprec_absolute(gel(x,2),bit);
1537 1178891 : return y;
1538 129452 : default: return x;
1539 : }
1540 : }
1541 :
1542 : static long
1543 529222 : a_posteriori_errors(GEN p, GEN roots_pol, long err)
1544 : {
1545 529222 : long i, n = degpol(p), e_max = -(long)EXPOBITS;
1546 : GEN sigma, shatzle;
1547 :
1548 529222 : err += (long)log2((double)n) + 1;
1549 529222 : if (err > -2) return 0;
1550 529222 : sigma = real2n(-err, LOWDEFAULTPREC);
1551 : /* 2 / ((s - 1)^(1/n) - 1) */
1552 529223 : shatzle = divur(2, subrs(sqrtnr(subrs(sigma,1),n), 1));
1553 2422097 : for (i=1; i<=n; i++)
1554 : {
1555 1892871 : pari_sp av = avma;
1556 1892871 : GEN x = root_error(n,i,roots_pol,err,shatzle);
1557 1892634 : long e = gexpo(x);
1558 1892700 : set_avma(av); if (e > e_max) e_max = e;
1559 1892789 : gel(roots_pol,i) = mygprec_absolute(gel(roots_pol,i), -e);
1560 : }
1561 529226 : return e_max;
1562 : }
1563 :
1564 : /********************************************************************/
1565 : /** **/
1566 : /** MAIN **/
1567 : /** **/
1568 : /********************************************************************/
1569 : static GEN
1570 1598461 : append_clone(GEN r, GEN a) { a = gclone(a); vectrunc_append(r, a); return a; }
1571 :
1572 : /* put roots in placeholder roots_pol so that |P - L_1...L_n| < 2^(-bit)|P|
1573 : * returns prod (x-roots_pol[i]) */
1574 : static GEN
1575 1515210 : split_complete(GEN p, long bit, GEN roots_pol)
1576 : {
1577 1515210 : long n = degpol(p);
1578 : pari_sp ltop;
1579 : GEN p1, F, G, a, b, m1, m2;
1580 :
1581 1515210 : if (n == 1)
1582 : {
1583 445975 : a = gneg_i(gdiv(gel(p,2), gel(p,3)));
1584 445978 : (void)append_clone(roots_pol,a); return p;
1585 : }
1586 1069235 : ltop = avma;
1587 1069235 : if (n == 2)
1588 : {
1589 576241 : F = gsub(gsqr(gel(p,3)), gmul2n(gmul(gel(p,2),gel(p,4)), 2));
1590 576238 : F = gsqrt(F, nbits2prec(bit));
1591 576244 : p1 = ginv(gmul2n(gel(p,4),1));
1592 576242 : a = gneg_i(gmul(gadd(F,gel(p,3)), p1));
1593 576246 : b = gmul(gsub(F,gel(p,3)), p1);
1594 576246 : a = append_clone(roots_pol,a);
1595 576245 : b = append_clone(roots_pol,b); set_avma(ltop);
1596 576246 : a = mygprec(a, 3*bit);
1597 576247 : b = mygprec(b, 3*bit);
1598 576246 : return gmul(gel(p,4), mkpoln(3, gen_1, gneg(gadd(a,b)), gmul(a,b)));
1599 : }
1600 492994 : split_0(p,bit,&F,&G);
1601 492996 : m1 = split_complete(F,bit,roots_pol);
1602 492996 : m2 = split_complete(G,bit,roots_pol);
1603 492995 : return gerepileupto(ltop, gmul(m1,m2));
1604 : }
1605 :
1606 : static GEN
1607 6439073 : quicktofp(GEN x)
1608 : {
1609 6439073 : const long prec = DEFAULTPREC;
1610 6439073 : switch(typ(x))
1611 : {
1612 6418627 : case t_INT: return itor(x, prec);
1613 8809 : case t_REAL: return rtor(x, prec);
1614 0 : case t_FRAC: return fractor(x, prec);
1615 11636 : case t_COMPLEX: {
1616 11636 : GEN a = gel(x,1), b = gel(x,2);
1617 : /* avoid problem with 0, e.g. x = 0 + I*1e-100. We don't want |x| = 0. */
1618 11636 : if (isintzero(a)) return cxcompotor(b, prec);
1619 11594 : if (isintzero(b)) return cxcompotor(a, prec);
1620 11594 : a = cxcompotor(a, prec);
1621 11594 : b = cxcompotor(b, prec); return sqrtr(addrr(sqrr(a), sqrr(b)));
1622 : }
1623 1 : default: pari_err_TYPE("quicktofp",x);
1624 : return NULL;/*LCOV_EXCL_LINE*/
1625 : }
1626 :
1627 : }
1628 :
1629 : /* bound log_2 |largest root of p| (Fujiwara's bound) */
1630 : double
1631 2002513 : fujiwara_bound(GEN p)
1632 : {
1633 2002513 : pari_sp av = avma;
1634 2002513 : long i, n = degpol(p);
1635 : GEN cc;
1636 : double loglc, Lmax;
1637 :
1638 2002510 : if (n <= 0) pari_err_CONSTPOL("fujiwara_bound");
1639 2002510 : loglc = mydbllog2r( quicktofp(gel(p,n+2)) ); /* log_2 |lc(p)| */
1640 2002484 : cc = gel(p, 2);
1641 2002484 : if (gequal0(cc))
1642 612112 : Lmax = -pariINFINITY-1;
1643 : else
1644 1390383 : Lmax = (mydbllog2r(quicktofp(cc)) - loglc - 1) / n;
1645 6798347 : for (i = 1; i < n; i++)
1646 : {
1647 4795876 : GEN y = gel(p,i+2);
1648 : double L;
1649 4795876 : if (gequal0(y)) continue;
1650 3046310 : L = (mydbllog2r(quicktofp(y)) - loglc) / (n-i);
1651 3046303 : if (L > Lmax) Lmax = L;
1652 : }
1653 2002471 : return gc_double(av, Lmax+1);
1654 : }
1655 :
1656 : /* Fujiwara's bound, real roots. Based on the following remark: if
1657 : * p = x^n + sum a_i x^i and q = x^n + sum min(a_i,0)x^i
1658 : * then for all x >= 0, p(x) >= q(x). Thus any bound for the (positive) roots
1659 : * of q is a bound for the positive roots of p. */
1660 : double
1661 1147479 : fujiwara_bound_real(GEN p, long sign)
1662 : {
1663 1147479 : pari_sp av = avma;
1664 : GEN x;
1665 1147479 : long n = degpol(p), i, signodd, signeven;
1666 1147478 : if (n <= 0) pari_err_CONSTPOL("fujiwara_bound");
1667 1147478 : x = shallowcopy(p);
1668 1147483 : if (gsigne(gel(x, n+2)) > 0)
1669 1147461 : { signeven = 1; signodd = sign; }
1670 : else
1671 21 : { signeven = -1; signodd = -sign; }
1672 4938225 : for (i = 0; i < n; i++)
1673 : {
1674 3790741 : if ((n - i) % 2)
1675 2199814 : { if (gsigne(gel(x, i+2)) == signodd ) gel(x, i+2) = gen_0; }
1676 : else
1677 1590927 : { if (gsigne(gel(x, i+2)) == signeven) gel(x, i+2) = gen_0; }
1678 : }
1679 1147484 : return gc_double(av, fujiwara_bound(x));
1680 : }
1681 :
1682 : static GEN
1683 2152471 : mygprecrc_special(GEN x, long prec, long e)
1684 : {
1685 : GEN y;
1686 2152471 : switch(typ(x))
1687 : {
1688 34185 : case t_REAL:
1689 34185 : if (!signe(x)) return real_0_bit(minss(e, expo(x)));
1690 33205 : return (prec > realprec(x))? rtor(x, prec): x;
1691 12394 : case t_COMPLEX:
1692 12394 : y = cgetg(3,t_COMPLEX);
1693 12394 : gel(y,1) = mygprecrc_special(gel(x,1),prec,e);
1694 12394 : gel(y,2) = mygprecrc_special(gel(x,2),prec,e);
1695 12394 : return y;
1696 2105892 : default: return x;
1697 : }
1698 : }
1699 :
1700 : /* like mygprec but keep at least the same precision as before */
1701 : static GEN
1702 529229 : mygprec_special(GEN x, long bit)
1703 : {
1704 529229 : long e = gexpo(x) - bit, prec = nbits2prec(bit);
1705 529223 : switch(typ(x))
1706 : {
1707 2656906 : case t_POL: pari_APPLY_pol_normalized(mygprecrc_special(gel(x,i),prec,e));
1708 0 : default: return mygprecrc_special(x,prec,e);
1709 : }
1710 : }
1711 :
1712 : static GEN
1713 393219 : fix_roots1(GEN R)
1714 : {
1715 393219 : long i, l = lg(R);
1716 393219 : GEN v = cgetg(l, t_VEC);
1717 1747564 : for (i=1; i < l; i++) { GEN r = gel(R,i); gel(v,i) = gcopy(r); gunclone(r); }
1718 393224 : return v;
1719 : }
1720 : static GEN
1721 529223 : fix_roots(GEN R, long h, long bit)
1722 : {
1723 : long i, j, c, n, prec;
1724 : GEN v, Z, gh;
1725 :
1726 529223 : if (h == 1) return fix_roots1(R);
1727 136004 : prec = nbits2prec(bit); Z = grootsof1(h, prec); gh = utoipos(h);
1728 136003 : n = lg(R)-1; v = cgetg(h*n + 1, t_VEC);
1729 380115 : for (c = i = 1; i <= n; i++)
1730 : {
1731 244110 : GEN s, r = gel(R,i);
1732 244110 : s = (h == 2)? gsqrt(r, prec): gsqrtn(r, gh, NULL, prec);
1733 782639 : for (j = 1; j <= h; j++) gel(v, c++) = gmul(s, gel(Z,j));
1734 244090 : gunclone(r);
1735 : }
1736 136005 : return v;
1737 : }
1738 :
1739 : static GEN
1740 528443 : all_roots(GEN p, long bit)
1741 : {
1742 528443 : long bit2, i, e, h, n = degpol(p), elc = gexpo(leading_coeff(p));
1743 528443 : GEN q, R, m, pd = RgX_deflate_max(p, &h);
1744 528445 : double fb = fujiwara_bound(pd);
1745 : pari_sp av;
1746 :
1747 528442 : if (fb < 0) fb = 0;
1748 528442 : bit2 = bit + maxss(gexpo(p), 0) + (long)ceil(log2(n / h) + 2 * fb);
1749 529224 : for (av = avma, i = 1, e = 0;; i++, set_avma(av))
1750 : {
1751 529225 : R = vectrunc_init(n+1);
1752 529224 : bit2 += e + (n << i);
1753 529224 : q = RgX_gtofp_bit(mygprec(pd,bit2), bit2);
1754 529220 : q[1] = evalsigne(1)|evalvarn(0);
1755 529220 : m = split_complete(q, bit2, R);
1756 529223 : R = fix_roots(R, h, bit2);
1757 529229 : q = mygprec_special(pd,bit2);
1758 529225 : q[1] = evalsigne(1)|evalvarn(0);
1759 529225 : e = gexpo(RgX_sub(q, m)) - elc + (long)log2((double)n) + 1;
1760 529223 : if (e < 0)
1761 : {
1762 529223 : if (e < -2*bit2) e = -2*bit2; /* avoid e = -oo */
1763 529223 : e = bit + a_posteriori_errors(p, R, e);
1764 529225 : if (e < 0) return R;
1765 : }
1766 779 : if (DEBUGLEVEL)
1767 0 : err_printf("all_roots: restarting, i = %ld, e = %ld\n", i,e);
1768 : }
1769 : }
1770 :
1771 : INLINE int
1772 930942 : isexactscalar(GEN x) { long tx = typ(x); return is_rational_t(tx); }
1773 :
1774 : static int
1775 239296 : isexactpol(GEN p)
1776 : {
1777 239296 : long i,n = degpol(p);
1778 1161982 : for (i=0; i<=n; i++)
1779 930942 : if (!isexactscalar(gel(p,i+2))) return 0;
1780 231040 : return 1;
1781 : }
1782 :
1783 : /* p(0) != 0 [for efficiency] */
1784 : static GEN
1785 231040 : solve_exact_pol(GEN p, long bit)
1786 : {
1787 231040 : long i, j, k, m, n = degpol(p), iroots = 0;
1788 231040 : GEN ex, factors, v = zerovec(n);
1789 :
1790 231040 : factors = ZX_squff(Q_primpart(p), &ex);
1791 462080 : for (i=1; i<lg(factors); i++)
1792 : {
1793 231040 : GEN roots_fact = all_roots(gel(factors,i), bit);
1794 231040 : n = degpol(gel(factors,i));
1795 231040 : m = ex[i];
1796 922042 : for (j=1; j<=n; j++)
1797 1382004 : for (k=1; k<=m; k++) v[++iroots] = roots_fact[j];
1798 : }
1799 231040 : return v;
1800 : }
1801 :
1802 : /* return the roots of p with absolute error bit */
1803 : static GEN
1804 239296 : roots_com(GEN q, long bit)
1805 : {
1806 : GEN L, p;
1807 239296 : long v = RgX_valrem_inexact(q, &p);
1808 239296 : int ex = isexactpol(p);
1809 239296 : if (!ex) p = RgX_normalize1(p);
1810 239296 : if (lg(p) == 3)
1811 0 : L = cgetg(1,t_VEC); /* constant polynomial */
1812 : else
1813 239296 : L = ex? solve_exact_pol(p,bit): all_roots(p,bit);
1814 239296 : if (v)
1815 : {
1816 3935 : GEN M, z, t = gel(q,2);
1817 : long i, x, y, l, n;
1818 :
1819 3935 : if (isrationalzero(t)) x = -bit;
1820 : else
1821 : {
1822 7 : n = gexpo(t);
1823 7 : x = n / v; l = degpol(q);
1824 35 : for (i = v; i <= l; i++)
1825 : {
1826 28 : t = gel(q,i+2);
1827 28 : if (isrationalzero(t)) continue;
1828 28 : y = (n - gexpo(t)) / i;
1829 28 : if (y < x) x = y;
1830 : }
1831 : }
1832 3935 : z = real_0_bit(x); l = v + lg(L);
1833 3935 : M = cgetg(l, t_VEC); L -= v;
1834 7933 : for (i = 1; i <= v; i++) gel(M,i) = z;
1835 11826 : for ( ; i < l; i++) gel(M,i) = gel(L,i);
1836 3935 : L = M;
1837 : }
1838 239296 : return L;
1839 : }
1840 :
1841 : static GEN
1842 1196984 : tocomplex(GEN x, long l, long bit)
1843 : {
1844 : GEN y;
1845 1196984 : if (typ(x) == t_COMPLEX)
1846 : {
1847 1177577 : if (signe(gel(x,1))) return mygprecrc(x, l, -bit);
1848 137089 : x = gel(x,2);
1849 137089 : y = cgetg(3,t_COMPLEX);
1850 137091 : gel(y,1) = real_0_bit(-bit);
1851 137091 : gel(y,2) = mygprecrc(x, l, -bit);
1852 : }
1853 : else
1854 : {
1855 19407 : y = cgetg(3,t_COMPLEX);
1856 19407 : gel(y,1) = mygprecrc(x, l, -bit);
1857 19407 : gel(y,2) = real_0_bit(-bit);
1858 : }
1859 156498 : return y;
1860 : }
1861 :
1862 : /* x,y are t_COMPLEX of t_REALs or t_REAL, compare wrt |Im x| - |Im y|,
1863 : * then Re x - Re y, up to 2^-e absolute error */
1864 : static int
1865 2222579 : cmp_complex_appr(void *E, GEN x, GEN y)
1866 : {
1867 2222579 : long e = (long)E;
1868 : GEN z, xi, yi, xr, yr;
1869 : long sz, sxi, syi;
1870 2222579 : if (typ(x) == t_COMPLEX) { xr = gel(x,1); xi = gel(x,2); sxi = signe(xi); }
1871 835853 : else { xr = x; xi = NULL; sxi = 0; }
1872 2222579 : if (typ(y) == t_COMPLEX) { yr = gel(y,1); yi = gel(y,2); syi = signe(yi); }
1873 557835 : else { yr = y; yi = NULL; syi = 0; }
1874 : /* Compare absolute values of imaginary parts */
1875 2222579 : if (!sxi)
1876 : {
1877 855229 : if (syi && expo(yi) >= e) return -1;
1878 : /* |Im x| ~ |Im y| ~ 0 */
1879 : }
1880 1367350 : else if (!syi)
1881 : {
1882 49878 : if (sxi && expo(xi) >= e) return 1;
1883 : /* |Im x| ~ |Im y| ~ 0 */
1884 : }
1885 : else
1886 : {
1887 1317472 : z = addrr_sign(xi, 1, yi, -1); sz = signe(z);
1888 1317431 : if (sz && expo(z) >= e) return (int)sz;
1889 : }
1890 : /* |Im x| ~ |Im y|, sort according to real parts */
1891 1326414 : z = subrr(xr, yr); sz = signe(z);
1892 1326431 : if (sz && expo(z) >= e) return (int)sz;
1893 : /* Re x ~ Re y. Place negative imaginary part before positive */
1894 583872 : return (int) (sxi - syi);
1895 : }
1896 :
1897 : static GEN
1898 528488 : clean_roots(GEN L, long l, long bit, long clean)
1899 : {
1900 528488 : long i, n = lg(L), ex = 5 - bit;
1901 528488 : GEN res = cgetg(n,t_COL);
1902 2422165 : for (i=1; i<n; i++)
1903 : {
1904 1893677 : GEN c = gel(L,i);
1905 1893677 : if (clean && isrealappr(c,ex))
1906 : {
1907 696701 : if (typ(c) == t_COMPLEX) c = gel(c,1);
1908 696701 : c = mygprecrc(c, l, -bit);
1909 : }
1910 : else
1911 1196977 : c = tocomplex(c, l, bit);
1912 1893677 : gel(res,i) = c;
1913 : }
1914 528488 : gen_sort_inplace(res, (void*)ex, &cmp_complex_appr, NULL);
1915 528485 : return res;
1916 : }
1917 :
1918 : /* the vector of roots of p, with absolute error 2^(- prec2nbits(l)) */
1919 : static GEN
1920 239324 : roots_aux(GEN p, long l, long clean)
1921 : {
1922 239324 : pari_sp av = avma;
1923 : long bit;
1924 : GEN L;
1925 :
1926 239324 : if (typ(p) != t_POL)
1927 : {
1928 21 : if (gequal0(p)) pari_err_ROOTS0("roots");
1929 14 : if (!isvalidcoeff(p)) pari_err_TYPE("roots",p);
1930 7 : return cgetg(1,t_COL); /* constant polynomial */
1931 : }
1932 239303 : if (!signe(p)) pari_err_ROOTS0("roots");
1933 239303 : checkvalidpol(p,"roots");
1934 239296 : if (lg(p) == 3) return cgetg(1,t_COL); /* constant polynomial */
1935 239296 : if (l < LOWDEFAULTPREC) l = LOWDEFAULTPREC;
1936 239296 : bit = prec2nbits(l);
1937 239296 : L = roots_com(p, bit);
1938 239296 : return gerepilecopy(av, clean_roots(L, l, bit, clean));
1939 : }
1940 : GEN
1941 8018 : roots(GEN p, long l) { return roots_aux(p,l, 0); }
1942 : /* clean up roots. If root is real replace it by its real part */
1943 : GEN
1944 231306 : cleanroots(GEN p, long l) { return roots_aux(p,l, 1); }
1945 :
1946 : /* private variant of conjvec. Allow non rational coefficients, shallow
1947 : * function. */
1948 : GEN
1949 84 : polmod_to_embed(GEN x, long prec)
1950 : {
1951 84 : GEN v, T = gel(x,1), A = gel(x,2);
1952 : long i, l;
1953 84 : if (typ(A) != t_POL || varn(A) != varn(T))
1954 : {
1955 7 : checkvalidpol(T,"polmod_to_embed");
1956 7 : return const_col(degpol(T), A);
1957 : }
1958 77 : v = cleanroots(T,prec); l = lg(v);
1959 231 : for (i=1; i<l; i++) gel(v,i) = poleval(A,gel(v,i));
1960 77 : return v;
1961 : }
1962 :
1963 : GEN
1964 289189 : QX_complex_roots(GEN p, long l)
1965 : {
1966 289189 : pari_sp av = avma;
1967 : long bit, v;
1968 : GEN L;
1969 :
1970 289189 : if (!signe(p)) pari_err_ROOTS0("QX_complex_roots");
1971 289189 : if (lg(p) == 3) return cgetg(1,t_COL); /* constant polynomial */
1972 289189 : if (l < LOWDEFAULTPREC) l = LOWDEFAULTPREC;
1973 289189 : bit = prec2nbits(l);
1974 289189 : v = RgX_valrem(p, &p);
1975 289189 : L = lg(p) > 3? all_roots(Q_primpart(p), bit): cgetg(1,t_COL);
1976 289192 : if (v) L = shallowconcat(const_vec(v, real_0_bit(-bit)), L);
1977 289192 : return gerepilecopy(av, clean_roots(L, l, bit, 1));
1978 : }
1979 :
1980 : /********************************************************************/
1981 : /** **/
1982 : /** REAL ROOTS OF INTEGER POLYNOMIAL **/
1983 : /** **/
1984 : /********************************************************************/
1985 :
1986 : /* Count sign changes in the coefficients of (x+1)^deg(P)*P(1/(x+1)), P
1987 : * has no rational root. The inversion is implicit (we take coefficients
1988 : * backwards). */
1989 : static long
1990 5201911 : X2XP1(GEN P, GEN *Premapped)
1991 : {
1992 5201911 : const pari_sp av = avma;
1993 5201911 : GEN v = shallowcopy(P);
1994 5201994 : long i, j, nb, s, dP = degpol(P), vlim = dP+2;
1995 :
1996 31888816 : for (j = 2; j < vlim; j++) gel(v, j+1) = addii(gel(v, j), gel(v, j+1));
1997 5201575 : s = -signe(gel(v, vlim));
1998 5201575 : vlim--; nb = 0;
1999 15062579 : for (i = 1; i < dP; i++)
2000 : {
2001 13047108 : long s2 = -signe(gel(v, 2));
2002 13047108 : int flag = (s2 == s);
2003 87259930 : for (j = 2; j < vlim; j++)
2004 : {
2005 74212768 : gel(v, j+1) = addii(gel(v, j), gel(v, j+1));
2006 74212822 : if (flag) flag = (s2 != signe(gel(v, j+1)));
2007 : }
2008 13047162 : if (s == signe(gel(v, vlim)))
2009 : {
2010 4508948 : if (++nb >= 2) return gc_long(av,2);
2011 3318721 : s = -s;
2012 : }
2013 : /* if flag is set there will be no further sign changes */
2014 11856935 : if (flag && (!Premapped || !nb)) return gc_long(av, nb);
2015 9860463 : vlim--;
2016 9860463 : if (gc_needed(av, 3))
2017 : {
2018 0 : if (DEBUGMEM>1) pari_warn(warnmem, "X2XP1, i = %ld/%ld", i, dP-1);
2019 0 : if (!Premapped) setlg(v, vlim + 2);
2020 0 : v = gerepilecopy(av, v);
2021 : }
2022 : }
2023 2015471 : if (vlim >= 2 && s == signe(gel(v, vlim))) nb++;
2024 2015471 : if (Premapped && nb == 1) *Premapped = v; else set_avma(av);
2025 2015135 : return nb;
2026 : }
2027 :
2028 : static long
2029 0 : _intervalcmp(GEN x, GEN y)
2030 : {
2031 0 : if (typ(x) == t_VEC) x = gel(x, 1);
2032 0 : if (typ(y) == t_VEC) y = gel(y, 1);
2033 0 : return gcmp(x, y);
2034 : }
2035 :
2036 : static GEN
2037 11096055 : _gen_nored(void *E, GEN x) { (void)E; return x; }
2038 : static GEN
2039 24543673 : _mp_add(void *E, GEN x, GEN y) { (void)E; return mpadd(x, y); }
2040 : static GEN
2041 0 : _mp_sub(void *E, GEN x, GEN y) { (void)E; return mpsub(x, y); }
2042 : static GEN
2043 4353986 : _mp_mul(void *E, GEN x, GEN y) { (void)E; return mpmul(x, y); }
2044 : static GEN
2045 6264619 : _mp_sqr(void *E, GEN x) { (void)E; return mpsqr(x); }
2046 : static GEN
2047 14337928 : _gen_one(void *E) { (void)E; return gen_1; }
2048 : static GEN
2049 322087 : _gen_zero(void *E) { (void)E; return gen_0; }
2050 :
2051 : static struct bb_algebra mp_algebra = { _gen_nored, _mp_add, _mp_sub,
2052 : _mp_mul, _mp_sqr, _gen_one, _gen_zero };
2053 :
2054 : static GEN
2055 34543379 : _mp_cmul(void *E, GEN P, long a, GEN x) {(void)E; return mpmul(gel(P,a+2), x);}
2056 :
2057 : /* Split the polynom P in two parts, whose coeffs have constant sign:
2058 : * P(X) = X^D*Pp + Pm. Also compute the two parts of the derivative of P,
2059 : * Pprimem = Pm', Pprimep = X*Pp'+ D*Pp => P' = X^(D-1)*Pprimep + Pprimem;
2060 : * Pprimep[i] = (i+D) Pp[i]. Return D */
2061 : static long
2062 165286 : split_pols(GEN P, GEN *pPp, GEN *pPm, GEN *pPprimep, GEN *pPprimem)
2063 : {
2064 165286 : long i, D, dP = degpol(P), s0 = signe(gel(P,2));
2065 : GEN Pp, Pm, Pprimep, Pprimem;
2066 508962 : for(i=1; i <= dP; i++)
2067 508963 : if (signe(gel(P, i+2)) == -s0) break;
2068 165286 : D = i;
2069 165286 : Pm = cgetg(D + 2, t_POL);
2070 165293 : Pprimem = cgetg(D + 1, t_POL);
2071 165292 : Pp = cgetg(dP-D + 3, t_POL);
2072 165293 : Pprimep = cgetg(dP-D + 3, t_POL);
2073 165290 : Pm[1] = Pp[1] = Pprimem[1] = Pprimep[1] = P[1];
2074 674240 : for(i=0; i < D; i++)
2075 : {
2076 508953 : GEN c = gel(P, i+2);
2077 508953 : gel(Pm, i+2) = c;
2078 508953 : if (i) gel(Pprimem, i+1) = mului(i, c);
2079 : }
2080 688381 : for(; i <= dP; i++)
2081 : {
2082 523113 : GEN c = gel(P, i+2);
2083 523113 : gel(Pp, i+2-D) = c;
2084 523113 : gel(Pprimep, i+2-D) = mului(i, c);
2085 : }
2086 165268 : *pPm = normalizepol_lg(Pm, D+2);
2087 165286 : *pPprimem = normalizepol_lg(Pprimem, D+1);
2088 165288 : *pPp = normalizepol_lg(Pp, dP-D+3);
2089 165289 : *pPprimep = normalizepol_lg(Pprimep, dP-D+3);
2090 165291 : return dP - degpol(*pPp);
2091 : }
2092 :
2093 : static GEN
2094 5183319 : bkeval_single_power(long d, GEN V)
2095 : {
2096 5183319 : long mp = lg(V) - 2;
2097 5183319 : if (d > mp) return gmul(gpowgs(gel(V, mp+1), d/mp), gel(V, (d%mp)+1));
2098 5183319 : return gel(V, d+1);
2099 : }
2100 :
2101 : static GEN
2102 5183471 : splitpoleval(GEN Pp, GEN Pm, GEN pows, long D, long bitprec)
2103 : {
2104 5183471 : GEN vp = gen_bkeval_powers(Pp, degpol(Pp), pows, NULL, &mp_algebra, _mp_cmul);
2105 5183063 : GEN vm = gen_bkeval_powers(Pm, degpol(Pm), pows, NULL, &mp_algebra, _mp_cmul);
2106 5183207 : GEN xa = bkeval_single_power(D, pows);
2107 : GEN r;
2108 5183344 : if (!signe(vp)) return vm;
2109 5183344 : vp = gmul(vp, xa);
2110 5182155 : r = gadd(vp, vm);
2111 5178916 : if (gexpo(vp) - (signe(r)? gexpo(r): 0) > prec2nbits(realprec(vp)) - bitprec)
2112 339054 : return NULL;
2113 4840904 : return r;
2114 : }
2115 :
2116 : /* optimized Cauchy bound for P = X^D*Pp + Pm, D > deg(Pm) */
2117 : static GEN
2118 165292 : splitcauchy(GEN Pp, GEN Pm, long prec)
2119 : {
2120 165292 : GEN S = gel(Pp,2), A = gel(Pm,2);
2121 165292 : long i, lPm = lg(Pm), lPp = lg(Pp);
2122 505927 : for (i=3; i < lPm; i++) { GEN c = gel(Pm,i); if (abscmpii(A, c) < 0) A = c; }
2123 523120 : for (i=3; i < lPp; i++) S = addii(S, gel(Pp, i));
2124 165289 : return subsr(1, rdivii(A, S, prec)); /* 1 + |Pm|_oo / |Pp|_1 */
2125 : }
2126 :
2127 : static GEN
2128 15138 : ZX_deg1root(GEN P, long prec)
2129 : {
2130 15138 : GEN a = gel(P,3), b = gel(P,2);
2131 15138 : if (is_pm1(a))
2132 : {
2133 15138 : b = itor(b, prec); if (signe(a) > 0) togglesign(b);
2134 15138 : return b;
2135 : }
2136 0 : return rdivii(negi(b), a, prec);
2137 : }
2138 :
2139 : /* Newton for polynom P, P(0)!=0, with unique sign change => one root in ]0,oo[
2140 : * P' has also at most one zero there */
2141 : static GEN
2142 165285 : polsolve(GEN P, long bitprec)
2143 : {
2144 : pari_sp av;
2145 : GEN Pp, Pm, Pprimep, Pprimem, Pprime, Pprime2, ra, rb, rc, Pc;
2146 165285 : long dP = degpol(P), prec = nbits2prec(bitprec);
2147 : long expoold, iter, D, rt, s0, bitaddprec, cprec, PREC;
2148 :
2149 165289 : if (dP == 1) return ZX_deg1root(P, prec);
2150 165289 : Pprime = ZX_deriv(P);
2151 165285 : Pprime2 = ZX_deriv(Pprime);
2152 165287 : bitaddprec = 1 + 2*expu(dP); PREC = prec + nbits2prec(bitaddprec);
2153 165286 : D = split_pols(P, &Pp, &Pm, &Pprimep, &Pprimem); /* P = X^D*Pp + Pm */
2154 165291 : s0 = signe(gel(P, 2));
2155 165291 : rt = maxss(D, brent_kung_optpow(maxss(degpol(Pp), degpol(Pm)), 2, 1));
2156 165292 : rb = splitcauchy(Pp, Pm, DEFAULTPREC);
2157 165276 : for (cprec = DEFAULTPREC, expoold = LONG_MAX;;)
2158 0 : {
2159 165276 : GEN pows = gen_powers(rb, rt, 1, NULL, _mp_sqr, _mp_mul, _gen_one);
2160 165290 : Pc = splitpoleval(Pp, Pm, pows, D, bitaddprec);
2161 165282 : if (!Pc) { cprec += EXTRAPREC64; rb = rtor(rb, cprec); continue; }
2162 165282 : if (signe(Pc) != s0) break;
2163 0 : shiftr_inplace(rb,1);
2164 : }
2165 165282 : for (iter = 0, ra = NULL;;)
2166 1807813 : {
2167 : GEN wdth;
2168 1973095 : iter++;
2169 1973095 : if (ra)
2170 899124 : rc = shiftr(addrr(ra, rb), -1);
2171 : else
2172 1073971 : rc = shiftr(rb, -1);
2173 : for(;;)
2174 0 : {
2175 1973272 : GEN pows = gen_powers(rc, rt, 1, NULL, _mp_sqr, _mp_mul, _gen_one);
2176 1973086 : Pc = splitpoleval(Pp, Pm, pows, D, bitaddprec+2);
2177 1972857 : if (Pc) break;
2178 0 : cprec += EXTRAPREC64;
2179 0 : rc = rtor(rc, cprec);
2180 : }
2181 1972857 : if (signe(Pc) == s0)
2182 591321 : ra = rc;
2183 : else
2184 1381536 : rb = rc;
2185 1972857 : if (!ra) continue;
2186 1064144 : wdth = subrr(rb, ra);
2187 1064291 : if (!(iter % 8))
2188 : {
2189 166463 : GEN m1 = poleval(Pprime, ra), M2;
2190 166467 : if (signe(m1) == s0) continue;
2191 165312 : M2 = poleval(Pprime2, rb);
2192 165312 : if (abscmprr(gmul(M2, wdth), shiftr(m1, 1)) > 0) continue;
2193 162153 : break;
2194 : }
2195 897828 : else if (gexpo(wdth) <= -bitprec)
2196 3163 : break;
2197 : }
2198 165316 : rc = rb; av = avma;
2199 1358735 : for(;; rc = gerepileuptoleaf(av, rc))
2200 1358922 : {
2201 : long exponew;
2202 1524238 : GEN Ppc, dist, rcold = rc;
2203 1524238 : GEN pows = gen_powers(rc, rt, 1, NULL, _mp_sqr, _mp_mul, _gen_one);
2204 1523964 : Ppc = splitpoleval(Pprimep, Pprimem, pows, D-1, bitaddprec+4);
2205 1523662 : if (Ppc) Pc = splitpoleval(Pp, Pm, pows, D, bitaddprec+4);
2206 1523873 : if (!Ppc || !Pc)
2207 : {
2208 339065 : if (cprec >= PREC)
2209 44088 : cprec += EXTRAPREC64;
2210 : else
2211 294977 : cprec = minss(2*cprec, PREC);
2212 339070 : rc = rtor(rc, cprec); continue; /* backtrack one step */
2213 : }
2214 1184808 : dist = typ(Ppc) == t_REAL? divrr(Pc, Ppc): divri(Pc, Ppc);
2215 1185056 : rc = subrr(rc, dist);
2216 1184463 : if (cmprr(ra, rc) > 0 || cmprr(rb, rc) < 0)
2217 : {
2218 0 : if (cprec >= PREC) break;
2219 0 : cprec = minss(2*cprec, PREC);
2220 0 : rc = rtor(rcold, cprec); continue; /* backtrack one step */
2221 : }
2222 1184937 : if (expoold == LONG_MAX) { expoold = expo(dist); continue; }
2223 966638 : exponew = expo(dist);
2224 966638 : if (exponew < -bitprec - 1)
2225 : {
2226 230769 : if (cprec >= PREC) break;
2227 65483 : cprec = minss(2*cprec, PREC);
2228 65486 : rc = rtor(rc, cprec); continue;
2229 : }
2230 735869 : if (exponew > expoold - 2)
2231 : {
2232 53035 : if (cprec >= PREC) break;
2233 53035 : expoold = LONG_MAX;
2234 53035 : cprec = minss(2*cprec, PREC);
2235 53037 : rc = rtor(rc, cprec); continue;
2236 : }
2237 682834 : expoold = exponew;
2238 : }
2239 165286 : return rtor(rc, prec);
2240 : }
2241 :
2242 : /* Return primpart(P(x / 2)) */
2243 : static GEN
2244 1934983 : ZX_rescale2prim(GEN P)
2245 : {
2246 1934983 : long i, l = lg(P), v, n;
2247 : GEN Q;
2248 1934983 : if (l==2) return pol_0(varn(P));
2249 1934983 : Q = cgetg(l,t_POL); v = vali(gel(P,l-1));
2250 9700337 : for (i = l-2, n = 1; v > n && i >= 2; i--, n++)
2251 7765305 : v = minss(v, vali(gel(P,i)) + n);
2252 1935032 : gel(Q,l-1) = v? shifti(gel(P,l-1), -v): gel(P,l-1);
2253 11608811 : for (i = l-2, n = 1-v; i >= 2; i--, n++)
2254 9673845 : gel(Q,i) = shifti(gel(P,i), n);
2255 1934966 : Q[1] = P[1]; return Q;
2256 : }
2257 :
2258 : /* assume Q0 has no rational root */
2259 : static GEN
2260 939225 : usp(GEN Q0, long flag, long bitprec)
2261 : {
2262 939225 : const pari_sp av = avma;
2263 : GEN Qremapped, Q, c, Lc, Lk, sol;
2264 939225 : GEN *pQremapped = flag == 1? &Qremapped: NULL;
2265 939225 : const long prec = nbits2prec(bitprec), deg = degpol(Q0);
2266 939223 : long listsize = 64, nbr = 0, nb_todo, ind, indf, i, k, nb;
2267 :
2268 939223 : sol = zerocol(deg);
2269 939259 : Lc = zerovec(listsize);
2270 939278 : Lk = cgetg(listsize+1, t_VECSMALL);
2271 939278 : k = Lk[1] = 0;
2272 939278 : ind = 1; indf = 2;
2273 939278 : Q = Q0;
2274 939278 : c = gen_0;
2275 939278 : nb_todo = 1;
2276 6140960 : while (nb_todo)
2277 : {
2278 5201722 : GEN nc = gel(Lc, ind);
2279 : pari_sp av2;
2280 5201722 : if (Lk[ind] == k + 1)
2281 : {
2282 1934985 : Q = Q0 = ZX_rescale2prim(Q0);
2283 1934983 : c = gen_0;
2284 : }
2285 5201720 : if (!equalii(nc, c)) Q = ZX_translate(Q, subii(nc, c));
2286 5201796 : av2 = avma;
2287 5201796 : k = Lk[ind];
2288 5201796 : ind++;
2289 5201796 : c = nc;
2290 5201796 : nb_todo--;
2291 5201796 : nb = X2XP1(Q, pQremapped);
2292 :
2293 5201528 : if (nb == 1)
2294 : { /* exactly one root */
2295 1601994 : GEN s = gen_0;
2296 1601994 : if (flag == 0)
2297 : {
2298 0 : s = mkvec2(gmul2n(c,-k), gmul2n(addiu(c,1),-k));
2299 0 : s = gerepilecopy(av2, s);
2300 : }
2301 1601994 : else if (flag == 1) /* Caveat: Qremapped is the reciprocal polynomial */
2302 : {
2303 165286 : s = polsolve(*pQremapped, bitprec+1);
2304 165292 : s = addir(c, divrr(s, addsr(1, s)));
2305 165283 : shiftr_inplace(s, -k);
2306 165283 : if (realprec(s) != prec) s = rtor(s, prec);
2307 165292 : s = gerepileupto(av2, s);
2308 : }
2309 1436708 : else set_avma(av2);
2310 1602017 : gel(sol, ++nbr) = s;
2311 : }
2312 3599534 : else if (nb)
2313 : { /* unknown, add two nodes to refine */
2314 2131350 : if (indf + 2 > listsize)
2315 : {
2316 1628 : if (ind>1)
2317 : {
2318 4979 : for (i = ind; i < indf; i++)
2319 : {
2320 3351 : gel(Lc, i-ind+1) = gel(Lc, i);
2321 3351 : Lk[i-ind+1] = Lk[i];
2322 : }
2323 1628 : indf -= ind-1;
2324 1628 : ind = 1;
2325 : }
2326 1628 : if (indf + 2 > listsize)
2327 : {
2328 0 : listsize *= 2;
2329 0 : Lc = vec_lengthen(Lc, listsize);
2330 0 : Lk = vecsmall_lengthen(Lk, listsize);
2331 : }
2332 102469 : for (i = indf; i <= listsize; i++) gel(Lc, i) = gen_0;
2333 : }
2334 2131350 : gel(Lc, indf) = nc = shifti(c, 1);
2335 2131346 : gel(Lc, indf + 1) = addiu(nc, 1);
2336 2131350 : Lk[indf] = Lk[indf + 1] = k + 1;
2337 2131350 : indf += 2;
2338 2131350 : nb_todo += 2;
2339 : }
2340 5201551 : if (gc_needed(av, 2))
2341 : {
2342 0 : gerepileall(av, 6, &Q0, &Q, &c, &Lc, &Lk, &sol);
2343 0 : if (DEBUGMEM > 1) pari_warn(warnmem, "ZX_Uspensky", avma);
2344 : }
2345 : }
2346 939238 : setlg(sol, nbr+1);
2347 939240 : return gerepilecopy(av, sol);
2348 : }
2349 :
2350 : static GEN
2351 14 : ZX_Uspensky_equal_yes(GEN a, long flag, long bit)
2352 : {
2353 14 : if (flag == 2) return gen_1;
2354 7 : if (flag == 1 && typ(a) != t_REAL)
2355 : {
2356 7 : if (typ(a) == t_INT && !signe(a))
2357 0 : a = real_0_bit(bit);
2358 : else
2359 7 : a = gtofp(a, nbits2prec(bit));
2360 : }
2361 7 : return mkcol(a);
2362 : }
2363 : static GEN
2364 21 : ZX_Uspensky_no(long flag)
2365 21 : { return flag <= 1 ? cgetg(1, t_COL) : gen_0; }
2366 : /* ZX_Uspensky(P, [a,a], flag) */
2367 : static GEN
2368 28 : ZX_Uspensky_equal(GEN P, GEN a, long flag, long bit)
2369 : {
2370 28 : if (typ(a) != t_INFINITY && gequal0(poleval(P, a)))
2371 14 : return ZX_Uspensky_equal_yes(a, flag, bit);
2372 : else
2373 14 : return ZX_Uspensky_no(flag);
2374 : }
2375 : static int
2376 3350 : sol_ok(GEN r, GEN a, GEN b) { return gcmp(a, r) <= 0 && gcmp(r, b) <= 0; }
2377 :
2378 : /* P a ZX without real double roots; better if primitive and squarefree but
2379 : * caller should ensure that. If flag & 4 assume that P has no rational root
2380 : * (modest speedup) */
2381 : GEN
2382 1057566 : ZX_Uspensky(GEN P, GEN ab, long flag, long bitprec)
2383 : {
2384 1057566 : pari_sp av = avma;
2385 : GEN a, b, res, sol;
2386 : double fb;
2387 : long l, nbz, deg;
2388 :
2389 1057566 : if (ab)
2390 : {
2391 953069 : if (typ(ab) == t_VEC)
2392 : {
2393 925694 : if (lg(ab) != 3) pari_err_DIM("ZX_Uspensky");
2394 925694 : a = gel(ab, 1);
2395 925694 : b = gel(ab, 2);
2396 : }
2397 : else
2398 : {
2399 27375 : a = ab;
2400 27375 : b = mkoo();
2401 : }
2402 : }
2403 : else
2404 : {
2405 104497 : a = mkmoo();
2406 104497 : b = mkoo();
2407 : }
2408 1057566 : if (flag & 4)
2409 : {
2410 128325 : if (gcmp(a, b) >= 0) { set_avma(av); return ZX_Uspensky_no(flag); }
2411 128325 : flag &= ~4;
2412 128325 : sol = cgetg(1, t_COL);
2413 : }
2414 : else
2415 : {
2416 929241 : switch (gcmp(a, b))
2417 : {
2418 7 : case 1: set_avma(av); return ZX_Uspensky_no(flag);
2419 28 : case 0: return gerepilecopy(av, ZX_Uspensky_equal(P, a, flag, bitprec));
2420 : }
2421 929207 : sol = nfrootsQ(P);
2422 : }
2423 1057534 : nbz = 0; l = lg(sol);
2424 1057534 : if (l > 1)
2425 : {
2426 : long i, j;
2427 2699 : P = RgX_div(P, roots_to_pol(sol, varn(P)));
2428 2699 : if (!RgV_is_ZV(sol)) P = Q_primpart(P);
2429 6049 : for (i = j = 1; i < l; i++)
2430 3350 : if (sol_ok(gel(sol,i), a, b)) gel(sol,j++) = gel(sol,i);
2431 2699 : setlg(sol, j);
2432 2699 : if (flag == 2) { nbz = j-1; sol = utoi(nbz); }
2433 2552 : else if (flag == 1) sol = RgC_gtofp(sol, nbits2prec(bitprec));
2434 : }
2435 1054835 : else if (flag == 2) sol = gen_0;
2436 1057534 : deg = degpol(P);
2437 1057534 : if (deg == 0) return gerepilecopy(av, sol);
2438 1055594 : if (typ(a) == t_INFINITY && typ(b) != t_INFINITY && gsigne(b))
2439 : {
2440 28 : fb = fujiwara_bound_real(P, -1);
2441 28 : if (fb <= -pariINFINITY) a = gen_0;
2442 21 : else if (fb < 0) a = gen_m1;
2443 21 : else a = negi(int2n((long)ceil(fb)));
2444 : }
2445 1055594 : if (typ(b) == t_INFINITY && typ(a) != t_INFINITY && gsigne(a))
2446 : {
2447 21 : fb = fujiwara_bound_real(P, 1);
2448 21 : if (fb <= -pariINFINITY) b = gen_0;
2449 21 : else if (fb < 0) b = gen_1;
2450 7 : else b = int2n((long)ceil(fb));
2451 : }
2452 1055594 : if (typ(a) != t_INFINITY && typ(b) != t_INFINITY)
2453 : {
2454 : GEN d, ad, bd, diff;
2455 : long i;
2456 : /* can occur if one of a,b was initially a t_INFINITY */
2457 12582 : if (gequal(a,b)) return gerepilecopy(av, sol);
2458 12575 : d = lcmii(Q_denom(a), Q_denom(b));
2459 12575 : if (is_pm1(d)) { d = NULL; ad = a; bd = b; }
2460 : else
2461 14 : { P = ZX_rescale(P, d); ad = gmul(a, d); bd = gmul(b, d); }
2462 12575 : diff = subii(bd, ad);
2463 12575 : P = ZX_affine(P, diff, ad);
2464 12575 : res = usp(P, flag, bitprec);
2465 12575 : if (flag <= 1)
2466 : {
2467 34176 : for (i = 1; i < lg(res); i++)
2468 : {
2469 21916 : GEN z = gmul(diff, gel(res, i));
2470 21916 : if (typ(z) == t_VEC)
2471 : {
2472 0 : gel(z, 1) = gadd(ad, gel(z, 1));
2473 0 : gel(z, 2) = gadd(ad, gel(z, 2));
2474 : }
2475 : else
2476 21916 : z = gadd(ad, z);
2477 21916 : if (d) z = gdiv(z, d);
2478 21916 : gel(res, i) = z;
2479 : }
2480 12260 : sol = shallowconcat(sol, res);
2481 : }
2482 : else
2483 315 : nbz += lg(res) - 1;
2484 : }
2485 1055587 : if (typ(b) == t_INFINITY && (fb=fujiwara_bound_real(P, 1)) > -pariINFINITY)
2486 : {
2487 838643 : long bp = maxss((long)ceil(fb), 0);
2488 838642 : res = usp(ZX_unscale2n(P, bp), flag, bitprec);
2489 838643 : if (flag <= 1)
2490 70871 : sol = shallowconcat(sol, gmul2n(res, bp));
2491 : else
2492 767772 : nbz += lg(res)-1;
2493 : }
2494 1055591 : if (typ(a) == t_INFINITY && (fb=fujiwara_bound_real(P,-1)) > -pariINFINITY)
2495 : {
2496 88056 : long i, bm = maxss((long)ceil(fb), 0);
2497 88056 : res = usp(ZX_unscale2n(ZX_z_unscale(P, -1), bm), flag, bitprec);
2498 88056 : if (flag <= 1)
2499 : {
2500 74633 : for (i = 1; i < lg(res); i++)
2501 : {
2502 46748 : GEN z = gneg(gmul2n(gel(res, i), bm));
2503 46748 : if (typ(z) == t_VEC) swap(gel(z, 1), gel(z, 2));
2504 46748 : gel(res, i) = z;
2505 : }
2506 27885 : sol = shallowconcat(res, sol);
2507 : }
2508 : else
2509 60171 : nbz += lg(res)-1;
2510 : }
2511 1055589 : if (flag >= 2) return utoi(nbz);
2512 83215 : if (flag)
2513 83215 : sol = sort(sol);
2514 : else
2515 0 : sol = gen_sort(sol, (void *)_intervalcmp, cmp_nodata);
2516 83214 : return gerepileupto(av, sol);
2517 : }
2518 :
2519 : /* x a scalar */
2520 : static GEN
2521 42 : rootsdeg0(GEN x)
2522 : {
2523 42 : if (!is_real_t(typ(x))) pari_err_TYPE("realroots",x);
2524 35 : if (gequal0(x)) pari_err_ROOTS0("realroots");
2525 14 : return cgetg(1,t_COL); /* constant polynomial */
2526 : }
2527 : static void
2528 1849102 : checkbound(GEN a)
2529 : {
2530 1849102 : switch(typ(a))
2531 : {
2532 1849102 : case t_INT: case t_FRAC: case t_INFINITY: break;
2533 0 : default: pari_err_TYPE("polrealroots", a);
2534 : }
2535 1849102 : }
2536 : static GEN
2537 925986 : check_ab(GEN ab)
2538 : {
2539 : GEN a, b;
2540 925986 : if (!ab) return NULL;
2541 924551 : if (typ(ab) != t_VEC || lg(ab) != 3) pari_err_TYPE("polrootsreal",ab);
2542 924551 : a = gel(ab,1); checkbound(a);
2543 924551 : b = gel(ab,2); checkbound(b);
2544 924551 : if (typ(a) == t_INFINITY && inf_get_sign(a) < 0 &&
2545 448 : typ(b) == t_INFINITY && inf_get_sign(b) > 0) ab = NULL;
2546 924551 : return ab;
2547 : }
2548 : /* e^(1/h) assuming the h-th root is real, beware that sqrtnr assumes e >= 0 */
2549 : static GEN
2550 22596 : _sqrtnr(GEN e, long h)
2551 : {
2552 : long s;
2553 : GEN r;
2554 22596 : if (h == 2) return sqrtr(e);
2555 14 : s = signe(e); setsigne(e, 1); /* e < 0 is possible, implies h is odd */
2556 14 : r = sqrtnr(e, h); if (s < 0) setsigne(r, -1);
2557 14 : return r;
2558 : }
2559 : GEN
2560 50613 : realroots(GEN P, GEN ab, long prec)
2561 : {
2562 50613 : pari_sp av = avma;
2563 50613 : GEN sol = NULL, fa, ex;
2564 : long i, j, v, l;
2565 :
2566 50613 : ab = check_ab(ab);
2567 50613 : if (typ(P) != t_POL) return rootsdeg0(P);
2568 50592 : if (!RgX_is_ZX(P)) P = RgX_rescale_to_int(P);
2569 50592 : switch(degpol(P))
2570 : {
2571 14 : case -1: return rootsdeg0(gen_0);
2572 7 : case 0: return rootsdeg0(gel(P,2));
2573 : }
2574 50571 : v = ZX_valrem(Q_primpart(P), &P);
2575 50571 : fa = ZX_squff(P, &ex); l = lg(fa); sol = cgetg(l + 1, t_VEC);
2576 102683 : for (i = 1; i < l; i++)
2577 : {
2578 52112 : GEN Pi = gel(fa, i), soli, soli2;
2579 : long n, h;
2580 52112 : if (ab) h = 1; else Pi = ZX_deflate_max(Pi, &h);
2581 52112 : soli = ZX_Uspensky(Pi, odd(h)? ab: gen_0, 1, prec2nbits(prec));
2582 52112 : n = lg(soli); soli2 = odd(h)? NULL: cgetg(n, t_COL);
2583 119073 : for (j = 1; j < n; j++)
2584 : {
2585 66961 : GEN r = gel(soli, j); /* != 0 */
2586 66961 : if (typ(r) != t_REAL) gel(soli, j) = r = gtofp(r, prec);
2587 66961 : if (h > 1)
2588 : {
2589 77 : gel(soli, j) = r = _sqrtnr(r, h);
2590 77 : if (soli2) gel(soli2, j) = negr(r);
2591 : }
2592 : }
2593 52112 : if (soli2) soli = shallowconcat(soli, soli2);
2594 52112 : if (ex[i] > 1) soli = shallowconcat1( const_vec(ex[i], soli) );
2595 52112 : gel(sol, i) = soli;
2596 : }
2597 50571 : if (v && (!ab || (gsigne(gel(ab,1)) <= 0 && gsigne(gel(ab,2)) >= 0)))
2598 84 : gel(sol, i++) = const_col(v, real_0(prec));
2599 50571 : setlg(sol, i); if (i == 1) { set_avma(av); return cgetg(1,t_COL); }
2600 50557 : return gerepileupto(av, sort(shallowconcat1(sol)));
2601 : }
2602 : GEN
2603 48125 : ZX_realroots_irred(GEN P, long prec)
2604 : {
2605 48125 : long dP = degpol(P), j, n, h;
2606 : GEN sol, sol2;
2607 : pari_sp av;
2608 48125 : if (dP == 1) retmkvec(ZX_deg1root(P, prec));
2609 44894 : av = avma; P = ZX_deflate_max(P, &h);
2610 44894 : if (h == dP)
2611 : {
2612 11907 : GEN r = _sqrtnr(ZX_deg1root(P, prec), h);
2613 11906 : return gerepilecopy(av, odd(h)? mkvec(r): mkvec2(negr(r), r));
2614 : }
2615 32987 : sol = ZX_Uspensky(P, odd(h)? NULL: gen_0, 1 | 4, prec2nbits(prec));
2616 32987 : n = lg(sol); sol2 = odd(h)? NULL: cgetg(n, t_COL);
2617 132050 : for (j = 1; j < n; j++)
2618 : {
2619 99063 : GEN r = gel(sol, j);
2620 99063 : if (typ(r) != t_REAL) gel(sol, j) = r = gtofp(r, prec);
2621 99063 : if (h > 1)
2622 : {
2623 10612 : gel(sol, j) = r = _sqrtnr(r, h);
2624 10612 : if (sol2) gel(sol2, j) = negr(r);
2625 : }
2626 : }
2627 32987 : if (sol2) sol = shallowconcat(sol, sol2);
2628 32987 : return gerepileupto(av, sort(sol));
2629 : }
2630 :
2631 : static long
2632 119686 : ZX_sturm_i(GEN P, long flag)
2633 : {
2634 : pari_sp av;
2635 119686 : long h, r, dP = degpol(P);
2636 119686 : if (dP == 1) return 1;
2637 116421 : av = avma; P = ZX_deflate_max(P, &h);
2638 116421 : if (h == dP)
2639 : { /* now deg P = 1 */
2640 17921 : if (odd(h))
2641 658 : r = 1;
2642 : else
2643 17263 : r = (signe(gel(P,2)) != signe(gel(P,3)))? 2: 0;
2644 17921 : return gc_long(av, r);
2645 : }
2646 98500 : if (odd(h))
2647 76018 : r = itou(ZX_Uspensky(P, NULL, flag, 0));
2648 : else
2649 22482 : r = 2*itou(ZX_Uspensky(P, gen_0, flag, 0));
2650 98498 : return gc_long(av,r);
2651 : }
2652 : /* P nonconstant, squarefree ZX */
2653 : long
2654 875373 : ZX_sturmpart(GEN P, GEN ab)
2655 : {
2656 875373 : pari_sp av = avma;
2657 875373 : if (!check_ab(ab)) return ZX_sturm(P);
2658 873967 : return gc_long(av, itou(ZX_Uspensky(P, ab, 2, 0)));
2659 : }
2660 : /* P nonconstant, squarefree ZX */
2661 : long
2662 3736 : ZX_sturm(GEN P) { return ZX_sturm_i(P, 2); }
2663 : /* P irreducible ZX */
2664 : long
2665 115950 : ZX_sturm_irred(GEN P) { return ZX_sturm_i(P, 2 + 4); }
|