Code coverage tests

This page documents the degree to which the PARI/GP source code is tested by our public test suite, distributed with the source distribution in directory src/test/. This is measured by the gcov utility; we then process gcov output using the lcov frond-end.

We test a few variants depending on Configure flags on the pari.math.u-bordeaux.fr machine (x86_64 architecture), and agregate them in the final report:

The target is 90% coverage for all mathematical modules (given that branches depending on DEBUGLEVEL or DEBUGMEM are not covered). This script is run to produce the results below.

LCOV - code coverage report
Current view: top level - basemath - rootpol.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.10.0 lcov report (development 21188-db834f2) Lines: 1514 1658 91.3 %
Date: 2017-10-20 06:23:07 Functions: 109 116 94.0 %
Legend: Lines: hit not hit

          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. It is distributed in the hope that it will be useful, but WITHOUT
       8             : ANY WARRANTY WHATSOEVER.
       9             : 
      10             : Check the License for details. You should have received a copy of it, along
      11             : with the package; see the file 'COPYING'. If not, write to the Free Software
      12             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      13             : 
      14             : /*******************************************************************/
      15             : /*                                                                 */
      16             : /*                ROOTS OF COMPLEX POLYNOMIALS                     */
      17             : /*  (original code contributed by Xavier Gourdon, INRIA RR 1852)   */
      18             : /*                                                                 */
      19             : /*******************************************************************/
      20             : #include "pari.h"
      21             : #include "paripriv.h"
      22             : 
      23             : static const double pariINFINITY = 1./0.;
      24             : 
      25             : static long
      26       69527 : isvalidcoeff(GEN x)
      27             : {
      28       69527 :   switch (typ(x))
      29             :   {
      30       60301 :     case t_INT: case t_REAL: case t_FRAC: return 1;
      31        9212 :     case t_COMPLEX: return isvalidcoeff(gel(x,1)) && isvalidcoeff(gel(x,2));
      32             :   }
      33          14 :   return 0;
      34             : }
      35             : 
      36             : static void
      37        9564 : checkvalidpol(GEN p, const char *f)
      38             : {
      39        9564 :   long i,n = lg(p);
      40       60646 :   for (i=2; i<n; i++)
      41       51089 :     if (!isvalidcoeff(gel(p,i))) pari_err_TYPE(f, gel(p,i));
      42        9557 : }
      43             : 
      44             : /********************************************************************/
      45             : /**                                                                **/
      46             : /**                   FAST ARITHMETIC over Z[i]                    **/
      47             : /**                                                                **/
      48             : /********************************************************************/
      49             : static THREAD long KARASQUARE_LIMIT, COOKSQUARE_LIMIT;
      50             : 
      51             : /* fast sum of x,y: t_INT or t_COMPLEX(t_INT) */
      52             : static GEN
      53    13152449 : addCC(GEN x, GEN y)
      54             : {
      55             :   GEN z;
      56             : 
      57    13152449 :   if (typ(x) == t_INT)
      58             :   {
      59    11268823 :     if (typ(y) == t_INT) return addii(x,y);
      60             :     /* ty == t_COMPLEX */
      61        2044 :     z = cgetg(3,t_COMPLEX);
      62        2044 :     gel(z,1) = addii(x, gel(y,1));
      63        2044 :     gel(z,2) = icopy(gel(y,2)); return z;
      64             :   }
      65             :   /* tx == t_COMPLEX */
      66     1883626 :   z = cgetg(3,t_COMPLEX);
      67     1883626 :   if (typ(y) == t_INT)
      68             :   {
      69       22770 :     gel(z,1) = addii(gel(x,1),y);
      70       22770 :     gel(z,2) = icopy(gel(x,2)); return z;
      71             :   }
      72             :   /* ty == t_COMPLEX */
      73     1860856 :   gel(z,1) = addii(gel(x,1),gel(y,1));
      74     1860856 :   gel(z,2) = addii(gel(x,2),gel(y,2)); return z;
      75             : }
      76             : /* fast product of x,y: t_INT or t_COMPLEX(t_INT) */
      77             : static GEN
      78    25228910 : mulCC(GEN x, GEN y)
      79             : {
      80             :   GEN z;
      81             : 
      82    25228910 :   if (typ(x) == t_INT)
      83             :   {
      84    20755647 :     if (typ(y) == t_INT) return mulii(x,y);
      85             :     /* ty == t_COMPLEX */
      86       13035 :     z = cgetg(3,t_COMPLEX);
      87       13035 :     gel(z,1) = mulii(x, gel(y,1));
      88       13035 :     gel(z,2) = mulii(x, gel(y,2)); return z;
      89             :   }
      90             :   /* tx == t_COMPLEX */
      91     4473263 :   z = cgetg(3,t_COMPLEX);
      92     4473263 :   if (typ(y) == t_INT)
      93             :   {
      94     1265153 :     gel(z,1) = mulii(gel(x,1),y);
      95     1265153 :     gel(z,2) = mulii(gel(x,2),y); return z;
      96             :   }
      97             :   /* ty == t_COMPLEX */
      98             :   {
      99     3208110 :     pari_sp av = avma, tetpil;
     100             :     GEN p1, p2;
     101             : 
     102     3208110 :     p1 = mulii(gel(x,1),gel(y,1));
     103     3208110 :     p2 = mulii(gel(x,2),gel(y,2));
     104     6416220 :     y = mulii(addii(gel(x,1),gel(x,2)),
     105     6416220 :               addii(gel(y,1),gel(y,2)));
     106     3208110 :     x = addii(p1,p2); tetpil = avma;
     107     3208110 :     gel(z,1) = subii(p1,p2);
     108     3208110 :     gel(z,2) = subii(y,x); gerepilecoeffssp(av,tetpil,z+1,2);
     109     3208110 :     return z;
     110             :   }
     111             : }
     112             : /* fast squaring x: t_INT or t_COMPLEX(t_INT) */
     113             : static GEN
     114    21666222 : sqrCC(GEN x)
     115             : {
     116             :   GEN z;
     117             : 
     118    21666222 :   if (typ(x) == t_INT) return sqri(x);
     119             :   /* tx == t_COMPLEX */
     120     4354637 :   z = cgetg(3,t_COMPLEX);
     121             :   {
     122     4354637 :     pari_sp av = avma, tetpil;
     123             :     GEN y, p1, p2;
     124             : 
     125     4354637 :     p1 = sqri(gel(x,1));
     126     4354637 :     p2 = sqri(gel(x,2));
     127     4354637 :     y = sqri(addii(gel(x,1),gel(x,2)));
     128     4354637 :     x = addii(p1,p2); tetpil = avma;
     129     4354637 :     gel(z,1) = subii(p1,p2);
     130     4354637 :     gel(z,2) = subii(y,x); gerepilecoeffssp(av,tetpil,z+1,2);
     131     4354637 :     return z;
     132             :   }
     133             : }
     134             : 
     135             : static void
     136     3716384 : set_karasquare_limit(long bit)
     137             : {
     138     3716384 :   if (bit<600)       { KARASQUARE_LIMIT=8; COOKSQUARE_LIMIT=400; }
     139        1176 :   else if (bit<2000) { KARASQUARE_LIMIT=4; COOKSQUARE_LIMIT=200; }
     140           0 :   else if (bit<3000) { KARASQUARE_LIMIT=4; COOKSQUARE_LIMIT=125; }
     141           0 :   else if (bit<5000) { KARASQUARE_LIMIT=2; COOKSQUARE_LIMIT= 75; }
     142           0 :   else               { KARASQUARE_LIMIT=1; COOKSQUARE_LIMIT= 50; }
     143     3716384 : }
     144             : 
     145             : /* assume lP > 0, lP = lgpol(P) */
     146             : static GEN
     147     7713699 : CX_square_spec(GEN P, long lP)
     148             : {
     149             :   GEN s, t;
     150     7713699 :   long i, j, l, nn, n = lP - 1;
     151             :   pari_sp av;
     152             : 
     153     7713699 :   nn = n<<1; s = cgetg(nn+3,t_POL); s[1] = evalsigne(1)|evalvarn(0);
     154     7713699 :   gel(s,2) = sqrCC(gel(P,0)); /* i = 0 */
     155    19790160 :   for (i=1; i<=n; i++)
     156             :   {
     157    12076461 :     av = avma; l = (i+1)>>1;
     158    12076461 :     t = mulCC(gel(P,0), gel(P,i)); /* j = 0 */
     159    12076461 :     for (j=1; j<l; j++) t = addCC(t, mulCC(gel(P,j), gel(P,i-j)));
     160    12076461 :     t = gmul2n(t,1);
     161    12076461 :     if ((i & 1) == 0) t = addCC(t, sqrCC(gel(P,i>>1)));
     162    12076461 :     gel(s,i+2) = gerepileupto(av, t);
     163             :   }
     164     7713699 :   gel(s,nn+2) = sqrCC(gel(P,n)); /* i = nn */
     165    13952523 :   for (   ; i<nn; i++)
     166             :   {
     167     6238824 :     av = avma; l = (i+1)>>1;
     168     6238824 :     t = mulCC(gel(P,i-n),gel(P,n)); /* j = i-n */
     169     6238824 :     for (j=i-n+1; j<l; j++) t = addCC(t, mulCC(gel(P,j),gel(P,i-j)));
     170     6238824 :     t = gmul2n(t,1);
     171     6238824 :     if ((i & 1) == 0) t = addCC(t, sqrCC(gel(P,i>>1)));
     172     6238824 :     gel(s,i+2) = gerepileupto(av, t);
     173             :   }
     174     7713699 :   return normalizepol_lg(s, nn+3);
     175             : }
     176             : /* nx = lgpol(x) */
     177             : static GEN
     178           0 : RgX_s_mulspec(GEN x, long nx, long s)
     179             : {
     180             :   GEN z, t;
     181             :   long i;
     182           0 :   if (!s || !nx) return pol_0(0);
     183           0 :   z = cgetg(nx+2, t_POL); z[1] = evalsigne(1)|evalvarn(0); t = z + 2;
     184           0 :   for (i=0; i < nx; i++) gel(t,i) = gmulgs(gel(x,i), s);
     185           0 :   return z;
     186             : }
     187             : /* nx = lgpol(x), return x << s. Inefficient if s = 0... */
     188             : static GEN
     189           0 : RgX_shiftspec(GEN x, long nx, long s)
     190             : {
     191             :   GEN z, t;
     192             :   long i;
     193           0 :   if (!nx) return pol_0(0);
     194           0 :   z = cgetg(nx+2, t_POL); z[1] = evalsigne(1)|evalvarn(0); t = z + 2;
     195           0 :   for (i=0; i < nx; i++) gel(t,i) = gmul2n(gel(x,i), s);
     196           0 :   return z;
     197             : }
     198             : 
     199             : /* spec function. nP = lgpol(P) */
     200             : static GEN
     201     7855141 : karasquare(GEN P, long nP)
     202             : {
     203             :   GEN Q, s0, s1, s2, a, t;
     204     7855141 :   long n0, n1, i, l, N, N0, N1, n = nP - 1; /* degree(P) */
     205             :   pari_sp av;
     206             : 
     207     7855141 :   if (n <= KARASQUARE_LIMIT) return nP? CX_square_spec(P, nP): pol_0(0);
     208      140791 :   av = avma;
     209      140791 :   n0 = (n>>1) + 1; n1 = nP - n0;
     210      140791 :   s0 = karasquare(P, n0); Q = P + n0;
     211      140791 :   s2 = karasquare(Q, n1);
     212      140791 :   s1 = RgX_addspec_shallow(P, Q, n0, n1);
     213      140791 :   s1 = RgX_sub(karasquare(s1+2, lgpol(s1)), RgX_add(s0,s2));
     214      140791 :   N = (n<<1) + 1;
     215      140791 :   a = cgetg(N + 2, t_POL); a[1] = evalsigne(1)|evalvarn(0);
     216      140791 :   t = a+2; l = lgpol(s0); s0 += 2; N0 = n0<<1;
     217      140791 :   for (i=0; i < l;  i++) gel(t,i) = gel(s0,i);
     218      140791 :   for (   ; i < N0; i++) gel(t,i) = gen_0;
     219      140791 :   t = a+2 + N0; l = lgpol(s2); s2 += 2; N1 = N - N0;
     220      140791 :   for (i=0; i < l;  i++) gel(t,i) = gel(s2,i);
     221      140791 :   for (   ; i < N1; i++) gel(t,i) = gen_0;
     222      140791 :   t = a+2 + n0; l = lgpol(s1); s1 += 2;
     223      140791 :   for (i=0; i < l; i++)  gel(t,i) = gadd(gel(t,i), gel(s1,i));
     224      140791 :   return gerepilecopy(av, normalizepol_lg(a, N+2));
     225             : }
     226             : /* spec function. nP = lgpol(P) */
     227             : static GEN
     228     7432768 : cook_square(GEN P, long nP)
     229             : {
     230             :   GEN Q, p0, p1, p2, p3, q, r, t, vp, vm;
     231     7432768 :   long n0, n3, i, j, n = nP - 1;
     232             :   pari_sp av;
     233             : 
     234     7432768 :   if (n <= COOKSQUARE_LIMIT) return  nP? karasquare(P, nP): pol_0(0);
     235           0 :   av = avma;
     236             : 
     237           0 :   n0 = (n+1) >> 2; n3 = n+1 - 3*n0;
     238           0 :   p0 = P;
     239           0 :   p1 = p0+n0;
     240           0 :   p2 = p1+n0;
     241           0 :   p3 = p2+n0; /* lgpol(p0,p1,p2) = n0, lgpol(p3) = n3 */
     242             : 
     243           0 :   q = cgetg(8,t_VEC) + 4;
     244           0 :   Q = cook_square(p0, n0);
     245           0 :   r = RgX_addspec_shallow(p0,p2, n0,n0);
     246           0 :   t = RgX_addspec_shallow(p1,p3, n0,n3);
     247           0 :   gel(q,-1) = RgX_sub(r,t);
     248           0 :   gel(q,1)  = RgX_add(r,t);
     249           0 :   r = RgX_addspec_shallow(p0,RgX_shiftspec(p2,n0, 2)+2, n0,n0);
     250           0 :   t = gmul2n(RgX_addspec_shallow(p1,RgX_shiftspec(p3,n3, 2)+2, n0,n3), 1);
     251           0 :   gel(q,-2) = RgX_sub(r,t);
     252           0 :   gel(q,2)  = RgX_add(r,t);
     253           0 :   r = RgX_addspec_shallow(p0,RgX_s_mulspec(p2,n0, 9)+2, n0,n0);
     254           0 :   t = gmulsg(3, RgX_addspec_shallow(p1,RgX_s_mulspec(p3,n3, 9)+2, n0,n3));
     255           0 :   gel(q,-3) = RgX_sub(r,t);
     256           0 :   gel(q,3)  = RgX_add(r,t);
     257             : 
     258           0 :   r = new_chunk(7);
     259           0 :   vp = cgetg(4,t_VEC);
     260           0 :   vm = cgetg(4,t_VEC);
     261           0 :   for (i=1; i<=3; i++)
     262             :   {
     263           0 :     GEN a = gel(q,i), b = gel(q,-i);
     264           0 :     a = cook_square(a+2, lgpol(a));
     265           0 :     b = cook_square(b+2, lgpol(b));
     266           0 :     gel(vp,i) = RgX_add(b, a);
     267           0 :     gel(vm,i) = RgX_sub(b, a);
     268             :   }
     269           0 :   gel(r,0) = Q;
     270           0 :   gel(r,1) = gdivgs(gsub(gsub(gmulgs(gel(vm,2),9),gel(vm,3)),
     271           0 :                      gmulgs(gel(vm,1),45)),
     272             :                 60);
     273           0 :   gel(r,2) = gdivgs(gadd(gadd(gmulgs(gel(vp,1),270),gmulgs(Q,-490)),
     274           0 :                      gadd(gmulgs(gel(vp,2),-27),gmulgs(gel(vp,3),2))),
     275             :                 360);
     276           0 :   gel(r,3) = gdivgs(gadd(gadd(gmulgs(gel(vm,1),13),gmulgs(gel(vm,2),-8)),
     277           0 :                     gel(vm,3)),
     278             :                 48);
     279           0 :   gel(r,4) = gdivgs(gadd(gadd(gmulgs(Q,56),gmulgs(gel(vp,1),-39)),
     280           0 :                      gsub(gmulgs(gel(vp,2),12),gel(vp,3))),
     281             :                 144);
     282           0 :   gel(r,5) = gdivgs(gsub(gadd(gmulgs(gel(vm,1),-5),gmulgs(gel(vm,2),4)),
     283           0 :                      gel(vm,3)),
     284             :                 240);
     285           0 :   gel(r,6) = gdivgs(gadd(gadd(gmulgs(Q,-20),gmulgs(gel(vp,1),15)),
     286           0 :                      gadd(gmulgs(gel(vp,2),-6),gel(vp,3))),
     287             :                 720);
     288           0 :   q = cgetg(2*n+3,t_POL); q[1] = evalsigne(1)|evalvarn(0);
     289           0 :   t = q+2;
     290           0 :   for (i=0; i<=2*n; i++) gel(t,i) = gen_0;
     291           0 :   for (i=0; i<=6; i++,t += n0)
     292             :   {
     293           0 :     GEN h = gel(r,i);
     294           0 :     long d = lgpol(h);
     295           0 :     h += 2;
     296           0 :     for (j=0; j<d; j++) gel(t,j) = gadd(gel(t,j), gel(h,j));
     297             :   }
     298           0 :   return gerepilecopy(av, normalizepol_lg(q, 2*n+3));
     299             : }
     300             : 
     301             : static GEN
     302     3716384 : graeffe(GEN p)
     303             : {
     304             :   GEN p0, p1, s0, s1;
     305     3716384 :   long n = degpol(p), n0, n1, i;
     306             : 
     307     3716384 :   if (!n) return gcopy(p);
     308     3716384 :   n0 = (n>>1)+1; n1 = n+1 - n0; /* n1 <= n0 <= n1+1 */
     309     3716384 :   p0 = new_chunk(n0);
     310     3716384 :   p1 = new_chunk(n1);
     311    12255812 :   for (i=0; i<n1; i++)
     312             :   {
     313     8539428 :     p0[i] = p[2+(i<<1)];
     314     8539428 :     p1[i] = p[3+(i<<1)];
     315             :   }
     316     3716384 :   if (n1 != n0)
     317     1785676 :     p0[i] = p[2+(i<<1)];
     318     3716384 :   s0 = cook_square(p0, n0);
     319     3716384 :   s1 = cook_square(p1, n1);
     320     3716384 :   return RgX_sub(s0, RgX_shift_shallow(s1,1));
     321             : }
     322             : GEN
     323        6566 : ZX_graeffe(GEN p)
     324             : {
     325        6566 :   pari_sp av = avma;
     326             :   GEN p0, p1, s0, s1;
     327        6566 :   long n = degpol(p);
     328             : 
     329        6566 :   if (!n) return ZX_copy(p);
     330        6566 :   RgX_even_odd(p, &p0, &p1);
     331             :   /* p = p0(x^2) + x p1(x^2) */
     332        6566 :   s0 = ZX_sqr(p0);
     333        6566 :   s1 = ZX_sqr(p1);
     334        6566 :   return gerepileupto(av, ZX_sub(s0, RgX_shift_shallow(s1,1)));
     335             : }
     336             : GEN
     337          14 : polgraeffe(GEN p)
     338             : {
     339          14 :   pari_sp av = avma;
     340             :   GEN p0, p1, s0, s1;
     341          14 :   long n = degpol(p);
     342             : 
     343          14 :   if (typ(p) != t_POL) pari_err_TYPE("polgraeffe",p);
     344          14 :   n = degpol(p);
     345          14 :   if (!n) return gcopy(p);
     346          14 :   RgX_even_odd(p, &p0, &p1);
     347             :   /* p = p0(x^2) + x p1(x^2) */
     348          14 :   s0 = RgX_sqr(p0);
     349          14 :   s1 = RgX_sqr(p1);
     350          14 :   return gerepileupto(av, RgX_sub(s0, RgX_shift_shallow(s1,1)));
     351             : }
     352             : 
     353             : /********************************************************************/
     354             : /**                                                                **/
     355             : /**                       MODULUS OF ROOTS                         **/
     356             : /**                                                                **/
     357             : /********************************************************************/
     358             : 
     359             : /* Quick approximation to log2(|x|); first define y s.t. |y-x| < 2^-32 then
     360             :  * return y rounded to 2 ulp. In particular, if result < 2^21, absolute error
     361             :  * is bounded by 2^-31. If result > 2^21, it is correct to 2 ulp */
     362             : static double
     363    19847796 : mydbllog2i(GEN x)
     364             : {
     365             : #ifdef LONG_IS_64BIT
     366    17075640 :   const double W = 1/(4294967296. * 4294967296.); /* 2^-64 */
     367             : #else
     368     2772156 :   const double W = 1/4294967296.; /*2^-32*/
     369             : #endif
     370             :   GEN m;
     371    19847796 :   long lx = lgefint(x);
     372             :   double l;
     373    19847796 :   if (lx == 2) return -pariINFINITY;
     374    19711169 :   m = int_MSW(x);
     375    19711169 :   l = (double)(ulong)*m;
     376    19711169 :   if (lx == 3) return log2(l);
     377     9029125 :   l += ((double)(ulong)*int_precW(m)) * W;
     378             :   /* at least m = min(53,BIL) bits are correct in the mantissa, thus log2
     379             :    * is correct with error < log(1 + 2^-m) ~ 2^-m. Adding the correct
     380             :    * exponent BIL(lx-3) causes 1ulp further round-off error */
     381     9029125 :   return log2(l) + (double)(BITS_IN_LONG*(lx-3));
     382             : }
     383             : 
     384             : /* return log(|x|) or -pariINFINITY */
     385             : static double
     386     1115001 : mydbllogr(GEN x) {
     387     1115001 :   if (!signe(x)) return -pariINFINITY;
     388     1115001 :   return LOG2*dbllog2r(x);
     389             : }
     390             : 
     391             : /* return log2(|x|) or -pariINFINITY */
     392             : static double
     393    10227495 : mydbllog2r(GEN x) {
     394    10227495 :   if (!signe(x)) return -pariINFINITY;
     395    10174985 :   return dbllog2r(x);
     396             : }
     397             : static double
     398     6806962 : dbllog2mp(GEN x) { return typ(x) == t_INT? mydbllog2i(x): mydbllog2r(x); }
     399             : double
     400    26135066 : dbllog2(GEN z)
     401             : {
     402             :   double x, y;
     403    26135066 :   switch(typ(z))
     404             :   {
     405    14097197 :     case t_INT: return mydbllog2i(z);
     406         791 :     case t_FRAC: return mydbllog2i(gel(z,1))-mydbllog2i(gel(z,2));
     407     8633597 :     case t_REAL: return mydbllog2r(z);
     408             :     default: /*t_COMPLEX*/
     409     3403481 :       x = dbllog2mp(gel(z,1));
     410     3403481 :       y = dbllog2mp(gel(z,2));
     411     3403481 :       if (x == -pariINFINITY) return y;
     412     3375401 :       if (y == -pariINFINITY) return x;
     413     3291561 :       if (fabs(x-y) > 10) return maxdd(x,y);
     414     3212019 :       return x + 0.5*log2(1 + exp2(2*(y-x)));
     415             :   }
     416             : }
     417             : static GEN /* beware overflow */
     418      554146 : dblexp(double x) { return fabs(x) < 100.? dbltor(exp(x)): mpexp(dbltor(x)); }
     419             : 
     420             : /* find s such that  A_h <= 2^s <= 2 A_i  for one h and all i < n = deg(p),
     421             :  * with  A_i := (binom(n,i) lc(p) / p_i) ^ 1/(n-i), and  p = sum p_i X^i */
     422             : static long
     423     2892357 : findpower(GEN p)
     424             : {
     425     2892357 :   double x, L, mins = pariINFINITY;
     426     2892357 :   long n = degpol(p),i;
     427             : 
     428     2892357 :   L = dbllog2(gel(p,n+2)); /* log2(lc * binom(n,i)) */
     429    14263878 :   for (i=n-1; i>=0; i--)
     430             :   {
     431    11371521 :     L += log2((double)(i+1) / (double)(n-i));
     432    11371521 :     x = dbllog2(gel(p,i+2));
     433    11371521 :     if (x != -pariINFINITY)
     434             :     {
     435    11297373 :       double s = (L - x) / (double)(n-i);
     436    11297373 :       if (s < mins) mins = s;
     437             :     }
     438             :   }
     439     2892357 :   i = (long)ceil(mins);
     440     2892357 :   if (i - mins > 1 - 1e-12) i--;
     441     2892357 :   return i;
     442             : }
     443             : 
     444             : /* returns the exponent for logmodulus(), from the Newton diagram */
     445             : static long
     446      482641 : newton_polygon(GEN p, long k)
     447             : {
     448      482641 :   pari_sp av = avma;
     449             :   double *logcoef, slope;
     450      482641 :   long n = degpol(p), i, j, h, l, *vertex;
     451             : 
     452      482641 :   logcoef = (double*)stack_malloc_align((n+1)*sizeof(double), sizeof(double));
     453      482641 :   vertex = (long*)new_chunk(n+1);
     454             : 
     455             :   /* vertex[i] = 1 if i a vertex of convex hull, 0 otherwise */
     456      482641 :   for (i=0; i<=n; i++) { logcoef[i] = dbllog2(gel(p,2+i)); vertex[i] = 0; }
     457      482641 :   vertex[0] = 1; /* sentinel */
     458     2013814 :   for (i=0; i < n; i=h)
     459             :   {
     460     1531173 :     slope = logcoef[i+1]-logcoef[i];
     461     7219746 :     for (j = h = i+1; j<=n; j++)
     462             :     {
     463     5688573 :       double pij = (logcoef[j]-logcoef[i])/(double)(j-i);
     464     5688573 :       if (slope < pij) { slope = pij; h = j; }
     465             :     }
     466     1531173 :     vertex[h] = 1;
     467             :   }
     468      482641 :   h = k;   while (!vertex[h]) h++;
     469      482641 :   l = k-1; while (!vertex[l]) l--;
     470      482641 :   avma = av;
     471      482641 :   return (long)floor((logcoef[h]-logcoef[l])/(double)(h-l) + 0.5);
     472             : }
     473             : 
     474             : /* change z into z*2^e, where z is real or complex of real */
     475             : static void
     476     3261332 : myshiftrc(GEN z, long e)
     477             : {
     478     3261332 :   if (typ(z)==t_COMPLEX)
     479             :   {
     480      775158 :     if (signe(gel(z,1))) shiftr_inplace(gel(z,1), e);
     481      775158 :     if (signe(gel(z,2))) shiftr_inplace(gel(z,2), e);
     482             :   }
     483             :   else
     484     2486174 :     if (signe(z)) shiftr_inplace(z, e);
     485     3261332 : }
     486             : 
     487             : /* return z*2^e, where z is integer or complex of integer (destroy z) */
     488             : static GEN
     489    12221147 : myshiftic(GEN z, long e)
     490             : {
     491    12221147 :   if (typ(z)==t_COMPLEX)
     492             :   {
     493     2260412 :     gel(z,1) = signe(gel(z,1))? mpshift(gel(z,1),e): gen_0;
     494     2260412 :     gel(z,2) = mpshift(gel(z,2),e);
     495     2260412 :     return z;
     496             :   }
     497     9960735 :   return signe(z)? mpshift(z,e): gen_0;
     498             : }
     499             : 
     500             : static GEN
     501      622411 : RgX_gtofp_bit(GEN q, long bit)
     502             : {
     503      622411 :   if (bit < 0) bit = 0;
     504      622411 :   return RgX_gtofp(q, nbits2prec(bit));
     505             : }
     506             : 
     507             : static GEN
     508    23606562 : mygprecrc(GEN x, long prec, long e)
     509             : {
     510             :   GEN y;
     511    23606562 :   switch(typ(x))
     512             :   {
     513    17738259 :     case t_REAL: return signe(x)? rtor(x, prec): real_0_bit(e);
     514             :     case t_COMPLEX:
     515     4538534 :       y = cgetg(3,t_COMPLEX);
     516     4538534 :       gel(y,1) = mygprecrc(gel(x,1),prec,e);
     517     4538534 :       gel(y,2) = mygprecrc(gel(x,2),prec,e);
     518     4538534 :       return y;
     519     1329769 :     default: return gcopy(x);
     520             :   }
     521             : }
     522             : 
     523             : /* gprec behaves badly with the zero for polynomials.
     524             : The second parameter in mygprec is the precision in base 2 */
     525             : static GEN
     526     5885137 : mygprec(GEN x, long bit)
     527             : {
     528             :   long lx, i, e, prec;
     529             :   GEN y;
     530             : 
     531     5885137 :   if (bit < 0) bit = 0; /* should rarely happen */
     532     5885137 :   e = gexpo(x) - bit;
     533     5885137 :   prec = nbits2prec(bit);
     534     5885137 :   switch(typ(x))
     535             :   {
     536             :     case t_POL:
     537     3577944 :       y = cgetg_copy(x, &lx); y[1] = x[1];
     538     3577944 :       for (i=2; i<lx; i++) gel(y,i) = mygprecrc(gel(x,i),prec,e);
     539     3577944 :       break;
     540             : 
     541     2307193 :     default: y = mygprecrc(x,prec,e);
     542             :   }
     543     5885137 :   return y;
     544             : }
     545             : 
     546             : /* normalize a polynomial p, that is change it with coefficients in Z[i],
     547             : after making product by 2^shift */
     548             : static GEN
     549     1456724 : pol_to_gaussint(GEN p, long shift)
     550             : {
     551     1456724 :   long i, l = lg(p);
     552     1456724 :   GEN q = cgetg(l, t_POL); q[1] = p[1];
     553     1456724 :   for (i=2; i<l; i++) gel(q,i) = gtrunc2n(gel(p,i), shift);
     554     1456724 :   return q;
     555             : }
     556             : 
     557             : /* returns a polynomial q in Z[i][x] keeping bit bits of p */
     558             : static GEN
     559     1147009 : eval_rel_pol(GEN p, long bit)
     560             : {
     561             :   long i;
     562     8000871 :   for (i = 2; i < lg(p); i++)
     563     6853862 :     if (gequal0(gel(p,i))) gel(p,i) = gen_0; /* bad behavior of gexpo */
     564     1147009 :   return pol_to_gaussint(p, bit-gexpo(p)+1);
     565             : }
     566             : 
     567             : /* returns p(R*x)/R^n (in R or R[i]), R = exp(lrho), bit bits of precision */
     568             : static GEN
     569      129170 : homothetie(GEN p, double lrho, long bit)
     570             : {
     571             :   GEN q, r, t, iR;
     572      129170 :   long n = degpol(p), i;
     573             : 
     574      129170 :   iR = mygprec(dblexp(-lrho),bit);
     575      129170 :   q = mygprec(p, bit);
     576      129170 :   r = cgetg(n+3,t_POL); r[1] = p[1];
     577      129170 :   t = iR; r[n+2] = q[n+2];
     578      788293 :   for (i=n-1; i>0; i--)
     579             :   {
     580      659123 :     gel(r,i+2) = gmul(t, gel(q,i+2));
     581      659123 :     t = mulrr(t, iR);
     582             :   }
     583      129170 :   gel(r,2) = gmul(t, gel(q,2)); return r;
     584             : }
     585             : 
     586             : /* change q in 2^(n*e) p(x*2^(-e)), n=deg(q)  [ ~as above with R = 2^-e ]*/
     587             : static void
     588      792356 : homothetie2n(GEN p, long e)
     589             : {
     590      792356 :   if (e)
     591             :   {
     592      599398 :     long i,n = lg(p)-1;
     593      599398 :     for (i=2; i<=n; i++) myshiftrc(gel(p,i), (n-i)*e);
     594             :   }
     595      792356 : }
     596             : 
     597             : /* return 2^f * 2^(n*e) p(x*2^(-e)), n=deg(q) */
     598             : static void
     599     2582642 : homothetie_gauss(GEN p, long e, long f)
     600             : {
     601     2582642 :   if (e || f)
     602             :   {
     603     2353308 :     long i, n = lg(p)-1;
     604     2353308 :     for (i=2; i<=n; i++) gel(p,i) = myshiftic(gel(p,i), f+(n-i)*e);
     605             :   }
     606     2582642 : }
     607             : 
     608             : /* Lower bound on the modulus of the largest root z_0
     609             :  * k is set to an upper bound for #{z roots, |z-z_0| < eps} */
     610             : static double
     611     2892357 : lower_bound(GEN p, long *k, double eps)
     612             : {
     613     2892357 :   long n = degpol(p), i, j;
     614     2892357 :   pari_sp ltop = avma;
     615             :   GEN a, s, S, ilc;
     616             :   double r, R, rho;
     617             : 
     618     2892357 :   if (n < 4) { *k = n; return 0.; }
     619     1185032 :   S = cgetg(5,t_VEC);
     620     1185032 :   a = cgetg(5,t_VEC); ilc = gdiv(real_1(DEFAULTPREC), gel(p,n+2));
     621     1185032 :   for (i=1; i<=4; i++) gel(a,i) = gmul(ilc,gel(p,n+2-i));
     622             :   /* i = 1 split out from next loop for efficiency and initialization */
     623     1185032 :   s = gel(a,1);
     624     1185032 :   gel(S,1) = gneg(s); /* Newton sum S_i */
     625     1185032 :   rho = r = gtodouble(gabs(s,3));
     626     1185032 :   R = r / n;
     627     4740128 :   for (i=2; i<=4; i++)
     628             :   {
     629     3555096 :     s = gmulsg(i,gel(a,i));
     630     3555096 :     for (j=1; j<i; j++) s = gadd(s, gmul(gel(S,j),gel(a,i-j)));
     631     3555096 :     gel(S,i) = gneg(s); /* Newton sum S_i */
     632     3555096 :     r = gtodouble(gabs(s,3));
     633     3555096 :     if (r > 0.)
     634             :     {
     635     3547106 :       r = exp(log(r/n) / (double)i);
     636     3547106 :       if (r > R) R = r;
     637             :     }
     638             :   }
     639     1185032 :   if (R > 0. && eps < 1.2)
     640     1183799 :     *k = (long)floor((rho/R + n) / (1 + exp(-eps)*cos(eps)));
     641             :   else
     642        1233 :     *k = n;
     643     1185032 :   avma = ltop; return R;
     644             : }
     645             : 
     646             : /* log of modulus of the largest root of p with relative error tau. Assume
     647             :  * P(0) != 0 and P non constant */
     648             : static double
     649      309715 : logmax_modulus(GEN p, double tau)
     650             : {
     651             :   GEN r, q, aux, gunr;
     652      309715 :   pari_sp av, ltop = avma;
     653      309715 :   long i,k,n=degpol(p),nn,bit,M,e;
     654      309715 :   double rho,eps, tau2 = (tau > 3.0)? 0.5: tau/6.;
     655             : 
     656      309715 :   r = cgeti(BIGDEFAULTPREC);
     657      309715 :   av = avma;
     658             : 
     659      309715 :   eps = - 1/log(1.5*tau2); /* > 0 */
     660      309715 :   bit = (long) ((double) n*log2(1./tau2)+3*log2((double) n))+1;
     661      309715 :   gunr = real_1_bit(bit+2*n);
     662      309715 :   aux = gdiv(gunr, gel(p,2+n));
     663      309715 :   q = RgX_Rg_mul(p, aux); gel(q,2+n) = gunr;
     664      309715 :   e = findpower(q);
     665      309715 :   homothetie2n(q,e);
     666      309715 :   affsi(e, r);
     667      309715 :   q = pol_to_gaussint(q, bit);
     668      309715 :   M = (long) (log2( log(4.*n) / (2*tau2) )) + 2;
     669      309715 :   nn = n;
     670      309715 :   for (i=0,e=0;;)
     671             :   { /* nn = deg(q) */
     672     2892357 :     rho = lower_bound(q, &k, eps);
     673     2892357 :     if (rho > exp2(-(double)e)) e = (long)-floor(log2(rho));
     674     2892357 :     affii(shifti(addis(r,e), 1), r);
     675     2892357 :     if (++i == M) break;
     676             : 
     677     7747926 :     bit = (long) ((double)k * log2(1./tau2) +
     678     5165284 :                      (double)(nn-k)*log2(1./eps) + 3*log2((double)nn)) + 1;
     679     2582642 :     homothetie_gauss(q, e, bit-(long)floor(dbllog2(gel(q,2+nn))+0.5));
     680     2582642 :     nn -= RgX_valrem(q, &q);
     681     2582642 :     set_karasquare_limit(gexpo(q));
     682     2582642 :     q = gerepileupto(av, graeffe(q));
     683     2582642 :     tau2 *= 1.5; if (tau2 > 0.9) tau2 = 0.5;
     684     2582642 :     eps = -1/log(tau2); /* > 0 */
     685     2582642 :     e = findpower(q);
     686     2582642 :   }
     687      309715 :   if (!signe(r)) { avma = ltop; return 0.; }
     688      290417 :   r = itor(r, DEFAULTPREC); shiftr_inplace(r, -M);
     689      290417 :   avma = ltop; return -rtodbl(r) * LOG2; /* -log(2) sum e_i 2^-i */
     690             : }
     691             : 
     692             : static GEN
     693        8667 : RgX_normalize1(GEN x)
     694             : {
     695        8667 :   long i, n = lg(x)-1;
     696             :   GEN y;
     697        8681 :   for (i = n; i > 1; i--)
     698        8674 :     if (!gequal0( gel(x,i) )) break;
     699        8667 :   if (i == n) return x;
     700          14 :   pari_warn(warner,"normalizing a polynomial with 0 leading term");
     701          14 :   if (i == 1) pari_err_ROOTS0("roots");
     702          14 :   y = cgetg(i+1, t_POL); y[1] = x[1];
     703          14 :   for (; i > 1; i--) gel(y,i) = gel(x,i);
     704          14 :   return y;
     705             : }
     706             : 
     707             : static GEN
     708        4936 : polrootsbound_i(GEN P)
     709             : {
     710        4936 :   pari_sp av = avma;
     711             :   double d;
     712        4936 :   (void)RgX_valrem_inexact(P,&P);
     713        4936 :   P = RgX_normalize1(P);
     714        4936 :   switch(degpol(P))
     715             :   {
     716           7 :     case -1: pari_err_ROOTS0("roots");
     717          63 :     case 0:  avma = av; return gen_0;
     718             :   }
     719        4866 :   d = logmax_modulus(P, 0.01) + 0.01;
     720        4866 :   avma = av; return dblexp(d);
     721             : }
     722             : GEN
     723        4943 : polrootsbound(GEN P)
     724             : {
     725        4943 :   if (typ(P) != t_POL) pari_err_TYPE("polrootsbound",P);
     726        4936 :   checkvalidpol(P, "polrootsbound");
     727        4936 :   return polrootsbound_i(P);
     728             : }
     729             : 
     730             : /* log of modulus of the smallest root of p, with relative error tau */
     731             : static double
     732      108823 : logmin_modulus(GEN p, double tau)
     733             : {
     734      108823 :   pari_sp av = avma;
     735             :   double r;
     736             : 
     737      108823 :   if (gequal0(gel(p,2))) return -pariINFINITY;
     738      108823 :   r = - logmax_modulus(RgX_recip_shallow(p),tau);
     739      108823 :   avma = av; return r;
     740             : }
     741             : 
     742             : /* return the log of the k-th modulus (ascending order) of p, rel. error tau*/
     743             : static double
     744       52843 : logmodulus(GEN p, long k, double tau)
     745             : {
     746             :   GEN q;
     747       52843 :   long i, kk = k, imax, n = degpol(p), nn, bit, e;
     748       52843 :   pari_sp av, ltop=avma;
     749       52843 :   double r, tau2 = tau/6;
     750             : 
     751       52843 :   bit = (long)(n * (2. + log2(3.*n/tau2)));
     752       52843 :   av = avma;
     753       52843 :   q = gprec_w(p, nbits2prec(bit));
     754       52843 :   q = RgX_gtofp_bit(q, bit);
     755       52843 :   e = newton_polygon(q,k);
     756       52843 :   r = (double)e;
     757       52843 :   homothetie2n(q,e);
     758       52843 :   imax = (long)(log2(3./tau) + log2(log(4.*n)))+1;
     759      482641 :   for (i=1; i<imax; i++)
     760             :   {
     761      429798 :     q = eval_rel_pol(q,bit);
     762      429798 :     kk -= RgX_valrem(q, &q);
     763      429798 :     nn = degpol(q);
     764             : 
     765      429798 :     set_karasquare_limit(bit);
     766      429798 :     q = gerepileupto(av, graeffe(q));
     767      429798 :     e = newton_polygon(q,kk);
     768      429798 :     r += e / exp2((double)i);
     769      429798 :     q = RgX_gtofp_bit(q, bit);
     770      429798 :     homothetie2n(q,e);
     771             : 
     772      429798 :     tau2 *= 1.5; if (tau2 > 1.) tau2 = 1.;
     773      429798 :     bit = 1 + (long)(nn*(2. + log2(3.*nn/tau2)));
     774             :   }
     775       52843 :   avma = ltop; return -r * LOG2;
     776             : }
     777             : 
     778             : /* return the log of the k-th modulus r_k of p, rel. error tau, knowing that
     779             :  * rmin < r_k < rmax. This information helps because we may reduce precision
     780             :  * quicker */
     781             : static double
     782       52843 : logpre_modulus(GEN p, long k, double tau, double lrmin, double lrmax)
     783             : {
     784             :   GEN q;
     785       52843 :   long n = degpol(p), i, imax, imax2, bit;
     786       52843 :   pari_sp ltop = avma, av;
     787       52843 :   double lrho, aux, tau2 = tau/6.;
     788             : 
     789       52843 :   aux = (lrmax - lrmin) / 2. + 4*tau2;
     790       52843 :   imax = (long) log2(log((double)n)/ aux);
     791       52843 :   if (imax <= 0) return logmodulus(p,k,tau);
     792             : 
     793       51786 :   lrho  = (lrmin + lrmax) / 2;
     794       51786 :   av = avma;
     795       51786 :   bit = (long)(n*(2. + aux / LOG2 - log2(tau2)));
     796       51786 :   q = homothetie(p, lrho, bit);
     797       51786 :   imax2 = (long)(log2(3./tau * log(4.*n))) + 1;
     798       51786 :   if (imax > imax2) imax = imax2;
     799             : 
     800      159660 :   for (i=0; i<imax; i++)
     801             :   {
     802      107874 :     q = eval_rel_pol(q,bit);
     803      107874 :     set_karasquare_limit(bit);
     804      107874 :     q = gerepileupto(av, graeffe(q));
     805      107874 :     aux = 2*aux + 2*tau2;
     806      107874 :     tau2 *= 1.5;
     807      107874 :     bit = (long)(n*(2. + aux / LOG2 - log2(1-exp(-tau2))));
     808      107874 :     q = RgX_gtofp_bit(q, bit);
     809             :   }
     810       51786 :   aux = exp2((double)imax);
     811       51786 :   aux = logmodulus(q,k, aux*tau/3.) / aux;
     812       51786 :   avma = ltop; return lrho + aux;
     813             : }
     814             : 
     815             : static double
     816       64117 : ind_maxlog2(GEN q)
     817             : {
     818       64117 :   long i, k = -1;
     819       64117 :   double L = - pariINFINITY;
     820      175000 :   for (i=0; i<=degpol(q); i++)
     821             :   {
     822      110883 :     double d = dbllog2(gel(q,2+i));
     823      110883 :     if (d > L) { L = d; k = i; }
     824             :   }
     825       64117 :   return k;
     826             : }
     827             : 
     828             : /* Returns k such that r_k e^(-tau) < R < r_{k+1} e^tau.
     829             :  * Assume that l <= k <= n-l */
     830             : static long
     831       77384 : dual_modulus(GEN p, double lrho, double tau, long l)
     832             : {
     833       77384 :   long i, imax, delta_k = 0, n = degpol(p), nn, v2, v, bit, ll = l;
     834       77384 :   double tau2 = tau * 7./8.;
     835       77384 :   pari_sp av = avma;
     836             :   GEN q;
     837             : 
     838       77384 :   bit = 6*n - 5*l + (long)(n*(-log2(tau2) + tau2 * 8./7.));
     839       77384 :   q = homothetie(p, lrho, bit);
     840       77384 :   imax = (long)(log(log(2.*n)/tau2)/log(7./4.)+1);
     841             : 
     842      673454 :   for (i=0; i<imax; i++)
     843             :   {
     844      609337 :     q = eval_rel_pol(q,bit); v2 = n - degpol(q);
     845      609337 :     v = RgX_valrem(q, &q);
     846      609337 :     ll -= maxss(v, v2); if (ll < 0) ll = 0;
     847             : 
     848      609337 :     nn = degpol(q); delta_k += v;
     849      609337 :     if (!nn) return delta_k;
     850             : 
     851      596070 :     set_karasquare_limit(bit);
     852      596070 :     q = gerepileupto(av, graeffe(q));
     853      596070 :     tau2 *= 7./4.;
     854      596070 :     bit = 6*nn - 5*ll + (long)(nn*(-log2(tau2) + tau2 * 8./7.));
     855             :   }
     856       64117 :   avma = av; return delta_k + (long)ind_maxlog2(q);
     857             : }
     858             : 
     859             : /********************************************************************/
     860             : /**                                                                **/
     861             : /**              FACTORS THROUGH CIRCLE INTEGRATION                **/
     862             : /**                                                                **/
     863             : /********************************************************************/
     864             : /* l power of 2 */
     865             : static void
     866     2158173 : fft(GEN Omega, GEN p, GEN f, long step, long l)
     867             : {
     868             :   pari_sp ltop;
     869             :   long i, l1, l2, l3, rapi, step4;
     870             :   GEN f1, f2, f3, f02, f13, g02, g13, ff;
     871             : 
     872     2158173 :   if (l == 2)
     873             :   {
     874     1198492 :     gel(f,0) = gadd(gel(p,0),gel(p,step));
     875     1198492 :     gel(f,1) = gsub(gel(p,0),gel(p,step)); return;
     876             :   }
     877      959681 :   if (l == 4)
     878             :   {
     879      528671 :     f1 = gadd(gel(p,0),   gel(p,step<<1));
     880      528671 :     f2 = gsub(gel(p,0),   gel(p,step<<1));
     881      528671 :     f3 = gadd(gel(p,step),gel(p,3*step));
     882      528671 :     f02= gsub(gel(p,step),gel(p,3*step));
     883      528671 :     f02 = mulcxI(f02);
     884      528671 :     gel(f,0) = gadd(f1, f3);
     885      528671 :     gel(f,1) = gadd(f2, f02);
     886      528671 :     gel(f,2) = gsub(f1, f3);
     887      528671 :     gel(f,3) = gsub(f2, f02); return;
     888             :   }
     889             : 
     890      431010 :   ltop = avma;
     891      431010 :   l1 = l>>2; l2 = 2*l1; l3 = l1+l2; step4 = step<<2;
     892      431010 :   fft(Omega,p,          f,   step4,l1);
     893      431010 :   fft(Omega,p+step,     f+l1,step4,l1);
     894      431010 :   fft(Omega,p+(step<<1),f+l2,step4,l1);
     895      431010 :   fft(Omega,p+3*step,   f+l3,step4,l1);
     896             : 
     897      431010 :   ff = cgetg(l+1,t_VEC);
     898     1750532 :   for (i=0; i<l1; i++)
     899             :   {
     900     1319522 :     rapi = step*i;
     901     1319522 :     f1 = gmul(gel(Omega,rapi),    gel(f,i+l1));
     902     1319522 :     f2 = gmul(gel(Omega,rapi<<1), gel(f,i+l2));
     903     1319522 :     f3 = gmul(gel(Omega,3*rapi),  gel(f,i+l3));
     904             : 
     905     1319522 :     f02 = gadd(gel(f,i),f2);
     906     1319522 :     g02 = gsub(gel(f,i),f2);
     907     1319522 :     f13 = gadd(f1,f3);
     908     1319522 :     g13 = mulcxI(gsub(f1,f3));
     909             : 
     910     1319522 :     gel(ff,i+1)    = gadd(f02, f13);
     911     1319522 :     gel(ff,i+l1+1) = gadd(g02, g13);
     912     1319522 :     gel(ff,i+l2+1) = gsub(f02, f13);
     913     1319522 :     gel(ff,i+l3+1) = gsub(g02, g13);
     914             :   }
     915      431010 :   ff = gerepilecopy(ltop,ff);
     916      431010 :   for (i=0; i<l; i++) f[i] = ff[i+1];
     917             : }
     918             : 
     919             : GEN
     920           0 : FFTinit(long k, long prec)
     921             : {
     922           0 :   if (k <= 0) pari_err_DOMAIN("FFTinit", "k", "<=", gen_0, stoi(k));
     923           0 :   return grootsof1(1L << k, prec);
     924             : }
     925             : 
     926             : GEN
     927           0 : FFT(GEN x, GEN Omega)
     928             : {
     929           0 :   long i, l = lg(Omega), n = lg(x);
     930             :   GEN y, z;
     931           0 :   if (!is_vec_t(typ(x))) pari_err_TYPE("FFT",x);
     932           0 :   if (typ(Omega) != t_VEC) pari_err_TYPE("FFT",Omega);
     933           0 :   if (n > l) pari_err_DIM("FFT");
     934             : 
     935           0 :   if (n < l) {
     936           0 :     z = cgetg(l, t_VECSMALL); /* cf stackdummy */
     937           0 :     for (i = 1; i < n; i++) z[i] = x[i];
     938           0 :     for (     ; i < l; i++) gel(z,i) = gen_0;
     939             :   }
     940           0 :   else z = x;
     941           0 :   y = cgetg(l, t_VEC);
     942           0 :   fft(Omega+1, z+1, y+1, 1, l-1);
     943           0 :   return y;
     944             : }
     945             : 
     946             : /* returns 1 if p has only real coefficients, 0 else */
     947             : static int
     948       67979 : isreal(GEN p)
     949             : {
     950             :   long i;
     951      400666 :   for (i = lg(p)-1; i > 1; i--)
     952      353194 :     if (typ(gel(p,i)) == t_COMPLEX) return 0;
     953       47472 :   return 1;
     954             : }
     955             : 
     956             : /* x non complex */
     957             : static GEN
     958       47814 : abs_update_r(GEN x, double *mu) {
     959       47814 :   GEN y = gtofp(x, DEFAULTPREC);
     960       47814 :   double ly = mydbllogr(y); if (ly < *mu) *mu = ly;
     961       47814 :   setabssign(y); return y;
     962             : }
     963             : /* return |x|, low accuracy. Set *mu = min(log(y), *mu) */
     964             : static GEN
     965      995652 : abs_update(GEN x, double *mu) {
     966             :   GEN y, xr, yr;
     967             :   double ly;
     968      995652 :   if (typ(x) != t_COMPLEX) return abs_update_r(x, mu);
     969      950022 :   xr = gel(x,1);
     970      950022 :   yr = gel(x,2);
     971      950022 :   if (gequal0(xr)) return abs_update_r(yr,mu);
     972      949896 :   if (gequal0(yr)) return abs_update_r(xr,mu);
     973             :   /* have to treat 0 specially: 0E-10 + 1e-20 = 0E-10 */
     974      947838 :   xr = gtofp(xr, DEFAULTPREC);
     975      947838 :   yr = gtofp(yr, DEFAULTPREC);
     976      947838 :   y = sqrtr(addrr(sqrr(xr), sqrr(yr)));
     977      947838 :   ly = mydbllogr(y); if (ly < *mu) *mu = ly;
     978      947838 :   return y;
     979             : }
     980             : 
     981             : static void
     982       73396 : initdft(GEN *Omega, GEN *prim, long N, long Lmax, long bit)
     983             : {
     984       73396 :   long prec = nbits2prec(bit);
     985       73396 :   *Omega = grootsof1(Lmax, prec) + 1;
     986       73396 :   *prim = rootsof1u_cx(N, prec);
     987       73396 : }
     988             : 
     989             : static void
     990       35521 : parameters(GEN p, long *LMAX, double *mu, double *gamma,
     991             :            int polreal, double param, double param2)
     992             : {
     993             :   GEN q, pc, Omega, A, RU, prim, g, TWO;
     994       35521 :   long n = degpol(p), bit, NN, K, i, j, Lmax;
     995       35521 :   pari_sp av2, av = avma;
     996             : 
     997       35521 :   bit = gexpo(p) + (long)param2+8;
     998       35521 :   Lmax = 4; while (Lmax <= n) Lmax <<= 1;
     999       35521 :   NN = (long)(param*3.14)+1; if (NN < Lmax) NN = Lmax;
    1000       35521 :   K = NN/Lmax; if (K & 1) K++;
    1001       35521 :   NN = Lmax*K;
    1002       35521 :   if (polreal) K = K/2+1;
    1003             : 
    1004       35521 :   initdft(&Omega, &prim, NN, Lmax, bit);
    1005       35521 :   q = mygprec(p,bit) + 2;
    1006       35521 :   A = cgetg(Lmax+1,t_VEC); A++;
    1007       35521 :   pc= cgetg(Lmax+1,t_VEC); pc++;
    1008       35521 :   for (i=0; i <= n; i++) gel(pc,i)= gel(q,i);
    1009       35521 :   for (   ; i<Lmax; i++) gel(pc,i) = gen_0;
    1010             : 
    1011       35521 :   *mu = pariINFINITY;
    1012       35521 :   g = real_0_bit(-bit);
    1013       35521 :   TWO = real2n(1, DEFAULTPREC);
    1014       35521 :   av2 = avma;
    1015       35521 :   RU = gen_1;
    1016      138586 :   for (i=0; i<K; i++)
    1017             :   {
    1018      103065 :     if (i) {
    1019       67544 :       GEN z = RU;
    1020      455636 :       for (j=1; j<n; j++)
    1021             :       {
    1022      388092 :         gel(pc,j) = gmul(gel(q,j),z);
    1023      388092 :         z = gmul(z,RU); /* RU = prim^i, z=prim^(ij) */
    1024             :       }
    1025       67544 :       gel(pc,n) = gmul(gel(q,n),z);
    1026             :     }
    1027             : 
    1028      103065 :     fft(Omega,pc,A,1,Lmax);
    1029      117460 :     if (polreal && i>0 && i<K-1)
    1030       14395 :       for (j=0; j<Lmax; j++) g = addrr(g, divrr(TWO, abs_update(gel(A,j),mu)));
    1031             :     else
    1032       88670 :       for (j=0; j<Lmax; j++) g = addrr(g, invr(abs_update(gel(A,j),mu)));
    1033      103065 :     RU = gmul(RU, prim);
    1034      103065 :     if (gc_needed(av,1))
    1035             :     {
    1036           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"parameters");
    1037           0 :       gerepileall(av2,2, &g,&RU);
    1038             :     }
    1039             :   }
    1040       35521 :   *gamma = mydbllog2r(divru(g,NN));
    1041       35521 :   *LMAX = Lmax; avma = av;
    1042       35521 : }
    1043             : 
    1044             : /* NN is a multiple of Lmax */
    1045             : static void
    1046       37875 : dft(GEN p, long k, long NN, long Lmax, long bit, GEN F, GEN H, long polreal)
    1047             : {
    1048             :   GEN Omega, q, qd, pc, pd, A, B, C, RU, aux, U, W, prim, prim2;
    1049       37875 :   long n = degpol(p), i, j, K;
    1050             :   pari_sp ltop;
    1051             : 
    1052       37875 :   initdft(&Omega, &prim, NN, Lmax, bit);
    1053       37875 :   RU = cgetg(n+2,t_VEC) + 1;
    1054             : 
    1055       37875 :   K = NN/Lmax; if (polreal) K = K/2+1;
    1056       37875 :   q = mygprec(p,bit);
    1057       37875 :   qd = RgX_deriv(q);
    1058             : 
    1059       37875 :   A = cgetg(Lmax+1,t_VEC); A++;
    1060       37875 :   B = cgetg(Lmax+1,t_VEC); B++;
    1061       37875 :   C = cgetg(Lmax+1,t_VEC); C++;
    1062       37875 :   pc = cgetg(Lmax+1,t_VEC); pc++;
    1063       37875 :   pd = cgetg(Lmax+1,t_VEC); pd++;
    1064       37875 :   pc[0] = q[2];  for (i=n+1; i<Lmax; i++) gel(pc,i) = gen_0;
    1065       37875 :   pd[0] = qd[2]; for (i=n;   i<Lmax; i++) gel(pd,i) = gen_0;
    1066             : 
    1067       37875 :   ltop = avma;
    1068       37875 :   W = cgetg(k+1,t_VEC);
    1069       37875 :   U = cgetg(k+1,t_VEC);
    1070       37875 :   for (i=1; i<=k; i++) gel(W,i) = gel(U,i) = gen_0;
    1071             : 
    1072       37875 :   gel(RU,0) = gen_1;
    1073       37875 :   prim2 = gen_1;
    1074      120642 :   for (i=0; i<K; i++)
    1075             :   {
    1076       82767 :     gel(RU,1) = prim2;
    1077       82767 :     for (j=1; j<n; j++) gel(RU,j+1) = gmul(gel(RU,j),prim2);
    1078             :     /* RU[j] = prim^(ij)= prim2^j */
    1079             : 
    1080       82767 :     for (j=1; j<n; j++) gel(pd,j) = gmul(gel(qd,j+2),gel(RU,j));
    1081       82767 :     fft(Omega,pd,A,1,Lmax);
    1082       82767 :     for (j=1; j<=n; j++) gel(pc,j) = gmul(gel(q,j+2),gel(RU,j));
    1083       82767 :     fft(Omega,pc,B,1,Lmax);
    1084       82767 :     for (j=0; j<Lmax; j++) gel(C,j) = ginv(gel(B,j));
    1085       82767 :     for (j=0; j<Lmax; j++) gel(B,j) = gmul(gel(A,j),gel(C,j));
    1086       82767 :     fft(Omega,B,A,1,Lmax);
    1087       82767 :     fft(Omega,C,B,1,Lmax);
    1088             : 
    1089       82767 :     if (polreal) /* p has real coefficients */
    1090             :     {
    1091       61734 :       if (i>0 && i<K-1)
    1092             :       {
    1093       34524 :         for (j=1; j<=k; j++)
    1094             :         {
    1095       28659 :           gel(W,j) = gadd(gel(W,j), gshift(mulreal(gel(A,j+1),gel(RU,j+1)),1));
    1096       28659 :           gel(U,j) = gadd(gel(U,j), gshift(mulreal(gel(B,j),gel(RU,j)),1));
    1097             :         }
    1098             :       }
    1099             :       else
    1100             :       {
    1101      157076 :         for (j=1; j<=k; j++)
    1102             :         {
    1103      107072 :           gel(W,j) = gadd(gel(W,j), mulreal(gel(A,j+1),gel(RU,j+1)));
    1104      107072 :           gel(U,j) = gadd(gel(U,j), mulreal(gel(B,j),gel(RU,j)));
    1105             :         }
    1106             :       }
    1107             :     }
    1108             :     else
    1109             :     {
    1110       75565 :       for (j=1; j<=k; j++)
    1111             :       {
    1112       48667 :         gel(W,j) = gadd(gel(W,j), gmul(gel(A,j+1),gel(RU,j+1)));
    1113       48667 :         gel(U,j) = gadd(gel(U,j), gmul(gel(B,j),gel(RU,j)));
    1114             :       }
    1115             :     }
    1116       82767 :     prim2 = gmul(prim2,prim);
    1117       82767 :     gerepileall(ltop,3, &W,&U,&prim2);
    1118             :   }
    1119             : 
    1120      112753 :   for (i=1; i<=k; i++)
    1121             :   {
    1122       74878 :     aux=gel(W,i);
    1123       74878 :     for (j=1; j<i; j++) aux = gadd(aux, gmul(gel(W,i-j),gel(F,k+2-j)));
    1124       74878 :     gel(F,k+2-i) = gdivgs(aux,-i*NN);
    1125             :   }
    1126      112753 :   for (i=0; i<k; i++)
    1127             :   {
    1128       74878 :     aux=gel(U,k-i);
    1129       74878 :     for (j=1+i; j<k; j++) aux = gadd(aux,gmul(gel(F,2+j),gel(U,j-i)));
    1130       74878 :     gel(H,i+2) = gdivgs(aux,NN);
    1131             :   }
    1132       37875 : }
    1133             : 
    1134             : #define NEWTON_MAX 10
    1135             : static GEN
    1136      191303 : refine_H(GEN F, GEN G, GEN HH, long bit, long Sbit)
    1137             : {
    1138      191303 :   GEN H = HH, D, aux;
    1139      191303 :   pari_sp ltop = avma;
    1140             :   long error, i, bit1, bit2;
    1141             : 
    1142      191303 :   D = Rg_RgX_sub(gen_1, RgX_rem(RgX_mul(H,G),F)); error = gexpo(D);
    1143      191303 :   bit2 = bit + Sbit;
    1144      357625 :   for (i=0; error>-bit && i<NEWTON_MAX && error<=0; i++)
    1145             :   {
    1146      166322 :     if (gc_needed(ltop,1))
    1147             :     {
    1148           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"refine_H");
    1149           0 :       gerepileall(ltop,2, &D,&H);
    1150             :     }
    1151      166322 :     bit1 = -error + Sbit;
    1152      166322 :     aux = RgX_mul(mygprec(H,bit1), mygprec(D,bit1));
    1153      166322 :     aux = RgX_rem(mygprec(aux,bit1), mygprec(F,bit1));
    1154             : 
    1155      166322 :     bit1 = -error*2 + Sbit; if (bit1 > bit2) bit1 = bit2;
    1156      166322 :     H = RgX_add(mygprec(H,bit1), aux);
    1157      166322 :     D = Rg_RgX_sub(gen_1, RgX_rem(RgX_mul(H,G),F));
    1158      166322 :     error = gexpo(D); if (error < -bit1) error = -bit1;
    1159             :   }
    1160      191303 :   if (error > -bit/2) return NULL; /* FAIL */
    1161      191162 :   return gerepilecopy(ltop,H);
    1162             : }
    1163             : 
    1164             : /* return 0 if fails, 1 else */
    1165             : static long
    1166       37875 : refine_F(GEN p, GEN *F, GEN *G, GEN H, long bit, double gamma)
    1167             : {
    1168       37875 :   GEN f0, FF, GG, r, HH = H;
    1169       37875 :   long error, i, bit1 = 0, bit2, Sbit, Sbit2,  enh, normF, normG, n = degpol(p);
    1170       37875 :   pari_sp av = avma;
    1171             : 
    1172       37875 :   FF = *F; GG = RgX_divrem(p, FF, &r);
    1173       37875 :   error = gexpo(r); if (error <= -bit) error = 1-bit;
    1174       37875 :   normF = gexpo(FF);
    1175       37875 :   normG = gexpo(GG);
    1176       37875 :   enh = gexpo(H); if (enh < 0) enh = 0;
    1177       37875 :   Sbit = normF + 2*normG + enh + (long)(4.*log2((double)n)+gamma) + 1;
    1178       37875 :   Sbit2 = enh + 2*(normF+normG) + (long)(2.*gamma+5.*log2((double)n)) + 1;
    1179       37875 :   bit2 = bit + Sbit;
    1180      229037 :   for (i=0; error>-bit && i<NEWTON_MAX && error<=0; i++)
    1181             :   {
    1182      191303 :     if (bit1 == bit2 && i >= 2) { Sbit += n; Sbit2 += n; bit2 += n; }
    1183      191303 :     if (gc_needed(av,1))
    1184             :     {
    1185           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"refine_F");
    1186           0 :       gerepileall(av,4, &FF,&GG,&r,&HH);
    1187             :     }
    1188             : 
    1189      191303 :     bit1 = -error + Sbit2;
    1190      191303 :     HH = refine_H(mygprec(FF,bit1), mygprec(GG,bit1), mygprec(HH,bit1),
    1191             :                   1-error, Sbit2);
    1192      191303 :     if (!HH) return 0; /* FAIL */
    1193             : 
    1194      191162 :     bit1 = -error + Sbit;
    1195      191162 :     r = RgX_mul(mygprec(HH,bit1), mygprec(r,bit1));
    1196      191162 :     f0 = RgX_rem(mygprec(r,bit1), mygprec(FF,bit1));
    1197             : 
    1198      191162 :     bit1 = -2*error + Sbit; if (bit1 > bit2) bit1 = bit2;
    1199      191162 :     FF = gadd(mygprec(FF,bit1),f0);
    1200             : 
    1201      191162 :     bit1 = -3*error + Sbit; if (bit1 > bit2) bit1 = bit2;
    1202      191162 :     GG = RgX_divrem(mygprec(p,bit1), mygprec(FF,bit1), &r);
    1203      191162 :     error = gexpo(r); if (error < -bit1) error = -bit1;
    1204             :   }
    1205       37734 :   if (error>-bit) return 0; /* FAIL */
    1206       35521 :   *F = FF; *G = GG; return 1;
    1207             : }
    1208             : 
    1209             : /* returns F and G from the unit circle U such that |p-FG|<2^(-bit) |cd|,
    1210             : where cd is the leading coefficient of p */
    1211             : static void
    1212       35521 : split_fromU(GEN p, long k, double delta, long bit,
    1213             :             GEN *F, GEN *G, double param, double param2)
    1214             : {
    1215             :   GEN pp, FF, GG, H;
    1216       35521 :   long n = degpol(p), NN, bit2, Lmax;
    1217       35521 :   int polreal = isreal(p);
    1218             :   pari_sp ltop;
    1219             :   double mu, gamma;
    1220             : 
    1221       35521 :   pp = gdiv(p, gel(p,2+n));
    1222       35521 :   parameters(pp, &Lmax,&mu,&gamma, polreal,param,param2);
    1223             : 
    1224       35521 :   H  = cgetg(k+2,t_POL); H[1] = p[1];
    1225       35521 :   FF = cgetg(k+3,t_POL); FF[1]= p[1];
    1226       35521 :   gel(FF,k+2) = gen_1;
    1227             : 
    1228       35521 :   NN = (long)(0.5/delta); NN |= 1; if (NN < 2) NN = 2;
    1229       35521 :   NN *= Lmax; ltop = avma;
    1230             :   for(;;)
    1231             :   {
    1232       37875 :     bit2 = (long)(((double)NN*delta-mu)/LOG2) + gexpo(pp) + 8;
    1233       37875 :     dft(pp, k, NN, Lmax, bit2, FF, H, polreal);
    1234       37875 :     if (refine_F(pp,&FF,&GG,H,bit,gamma)) break;
    1235        2354 :     NN <<= 1; avma = ltop;
    1236        2354 :   }
    1237       35521 :   *G = gmul(GG,gel(p,2+n)); *F = FF;
    1238       35521 : }
    1239             : 
    1240             : static void
    1241       35521 : optimize_split(GEN p, long k, double delta, long bit,
    1242             :             GEN *F, GEN *G, double param, double param2)
    1243             : {
    1244       35521 :   long n = degpol(p);
    1245             :   GEN FF, GG;
    1246             : 
    1247       35521 :   if (k <= n/2)
    1248       26733 :     split_fromU(p,k,delta,bit,F,G,param,param2);
    1249             :   else
    1250             :   {
    1251        8788 :     split_fromU(RgX_recip_shallow(p),n-k,delta,bit,&FF,&GG,param,param2);
    1252        8788 :     *F = RgX_recip_shallow(GG);
    1253        8788 :     *G = RgX_recip_shallow(FF);
    1254             :   }
    1255       35521 : }
    1256             : 
    1257             : /********************************************************************/
    1258             : /**                                                                **/
    1259             : /**               SEARCH FOR SEPARATING CIRCLE                     **/
    1260             : /**                                                                **/
    1261             : /********************************************************************/
    1262             : 
    1263             : /* return p(2^e*x) *2^(-n*e) */
    1264             : static void
    1265           0 : scalepol2n(GEN p, long e)
    1266             : {
    1267           0 :   long i,n=lg(p)-1;
    1268           0 :   for (i=2; i<=n; i++) gel(p,i) = gmul2n(gel(p,i),(i-n)*e);
    1269           0 : }
    1270             : 
    1271             : /* returns p(x/R)*R^n */
    1272             : static GEN
    1273      301311 : scalepol(GEN p, GEN R, long bit)
    1274             : {
    1275             :   GEN q,aux,gR;
    1276             :   long i;
    1277             : 
    1278      301311 :   aux = gR = mygprec(R,bit); q = mygprec(p,bit);
    1279     1385147 :   for (i=lg(p)-2; i>=2; i--)
    1280             :   {
    1281     1083836 :     gel(q,i) = gmul(aux,gel(q,i));
    1282     1083836 :     aux = gmul(aux,gR);
    1283             :   }
    1284      301311 :   return q;
    1285             : }
    1286             : 
    1287             : /* return (conj(a)X-1)^n * p[ (X-a) / (conj(a)X-1) ] */
    1288             : static GEN
    1289       97374 : conformal_pol(GEN p, GEN a)
    1290             : {
    1291       97374 :   GEN z, r, ma = gneg(a), ca = gconj(a);
    1292       97374 :   long n = degpol(p), i;
    1293       97374 :   pari_sp av = avma;
    1294             : 
    1295       97374 :   z = mkpoln(2, ca, gen_m1);
    1296       97374 :   r = scalarpol(gel(p,2+n), 0);
    1297      341594 :   for (i=n-1; ; i--)
    1298             :   {
    1299      341594 :     r = RgX_addmulXn_shallow(r, gmul(ma,r), 1); /* r *= (X - a) */
    1300      341594 :     r = gadd(r, gmul(z, gel(p,2+i)));
    1301      438968 :     if (i == 0) return gerepileupto(av, r);
    1302      244220 :     z = RgX_addmulXn_shallow(gmul(z,ca), gneg(z), 1); /* z *= conj(a)X - 1 */
    1303      244220 :     if (gc_needed(av,2))
    1304             :     {
    1305           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"conformal_pol");
    1306           0 :       gerepileall(av,2, &r,&z);
    1307             :     }
    1308      244220 :   }
    1309             : }
    1310             : 
    1311             : static const double UNDEF = -100000.;
    1312             : 
    1313             : static double
    1314       35521 : logradius(double *radii, GEN p, long k, double aux, double *delta)
    1315             : {
    1316       35521 :   long i, n = degpol(p);
    1317             :   double lrho, lrmin, lrmax;
    1318       35521 :   if (k > 1)
    1319             :   {
    1320       23252 :     i = k-1; while (i>0 && radii[i] == UNDEF) i--;
    1321       23252 :     lrmin = logpre_modulus(p,k,aux, radii[i], radii[k]);
    1322             :   }
    1323             :   else /* k=1 */
    1324       12269 :     lrmin = logmin_modulus(p,aux);
    1325       35521 :   radii[k] = lrmin;
    1326             : 
    1327       35521 :   if (k+1<n)
    1328             :   {
    1329       29591 :     i = k+2; while (i<=n && radii[i] == UNDEF) i++;
    1330       29591 :     lrmax = logpre_modulus(p,k+1,aux, radii[k+1], radii[i]);
    1331             :   }
    1332             :   else /* k+1=n */
    1333        5930 :     lrmax = logmax_modulus(p,aux);
    1334       35521 :   radii[k+1] = lrmax;
    1335             : 
    1336       35521 :   lrho = radii[k];
    1337       82591 :   for (i=k-1; i>=1; i--)
    1338             :   {
    1339       47070 :     if (radii[i] == UNDEF || radii[i] > lrho)
    1340       32794 :       radii[i] = lrho;
    1341             :     else
    1342       14276 :       lrho = radii[i];
    1343             :   }
    1344       35521 :   lrho = radii[k+1];
    1345      153254 :   for (i=k+1; i<=n; i++)
    1346             :   {
    1347      117733 :     if (radii[i] == UNDEF || radii[i] < lrho)
    1348       69877 :       radii[i] = lrho;
    1349             :     else
    1350       47856 :       lrho = radii[i];
    1351             :   }
    1352       35521 :   *delta = (lrmax - lrmin) / 2;
    1353       35521 :   if (*delta > 1.) *delta = 1.;
    1354       35521 :   return (lrmin + lrmax) / 2;
    1355             : }
    1356             : 
    1357             : static void
    1358       35521 : update_radius(long n, double *radii, double lrho, double *par, double *par2)
    1359             : {
    1360       35521 :   double t, param = 0., param2 = 0.;
    1361             :   long i;
    1362      235845 :   for (i=1; i<=n; i++)
    1363             :   {
    1364      200324 :     radii[i] -= lrho;
    1365      200324 :     t = fabs(rtodbl( invr(subsr(1, dblexp(radii[i]))) ));
    1366      200324 :     param += t; if (t > 1.) param2 += log2(t);
    1367             :   }
    1368       35521 :   *par = param; *par2 = param2;
    1369       35521 : }
    1370             : 
    1371             : /* apply the conformal mapping then split from U */
    1372             : static void
    1373       32458 : conformal_mapping(double *radii, GEN ctr, GEN p, long k, long bit,
    1374             :                   double aux, GEN *F,GEN *G)
    1375             : {
    1376       32458 :   long bit2, n = degpol(p), i;
    1377       32458 :   pari_sp ltop = avma, av;
    1378             :   GEN q, FF, GG, a, R;
    1379             :   double lrho, delta, param, param2;
    1380             :   /* n * (2.*log2(2.732)+log2(1.5)) + 1 */
    1381       32458 :   bit2 = bit + (long)(n*3.4848775) + 1;
    1382       32458 :   a = sqrtr_abs( stor(3, 2*MEDDEFAULTPREC - 2) );
    1383       32458 :   a = divrs(a, -6);
    1384       32458 :   a = gmul(mygprec(a,bit2), mygprec(ctr,bit2)); /* a = -ctr/2sqrt(3) */
    1385             : 
    1386       32458 :   av = avma;
    1387       32458 :   q = conformal_pol(mygprec(p,bit2), a);
    1388      203255 :   for (i=1; i<=n; i++)
    1389      170797 :     if (radii[i] != UNDEF) /* update array radii */
    1390             :     {
    1391      119349 :       pari_sp av2 = avma;
    1392      119349 :       GEN t, r = dblexp(radii[i]), r2 = sqrr(r);
    1393             :       /* 2(r^2 - 1) / (r^2 - 3(r-1)) */
    1394      119349 :       t = divrr(shiftr((subrs(r2,1)),1), subrr(r2, mulur(3,subrs(r,1))));
    1395      119349 :       radii[i] = mydbllogr(addsr(1,t)) / 2;
    1396      119349 :       avma = av2;
    1397             :     }
    1398       32458 :   lrho = logradius(radii, q,k,aux/10., &delta);
    1399       32458 :   update_radius(n, radii, lrho, &param, &param2);
    1400             : 
    1401       32458 :   bit2 += (long)(n * fabs(lrho)/LOG2 + 1.);
    1402       32458 :   R = mygprec(dblexp(-lrho), bit2);
    1403       32458 :   q = scalepol(q,R,bit2);
    1404       32458 :   gerepileall(av,2, &q,&R);
    1405             : 
    1406       32458 :   optimize_split(q,k,delta,bit2,&FF,&GG,param,param2);
    1407       32458 :   bit2 += n; R = invr(R);
    1408       32458 :   FF = scalepol(FF,R,bit2);
    1409       32458 :   GG = scalepol(GG,R,bit2);
    1410             : 
    1411       32458 :   a = mygprec(a,bit2);
    1412       32458 :   FF = conformal_pol(FF,a);
    1413       32458 :   GG = conformal_pol(GG,a);
    1414             : 
    1415       32458 :   a = invr(subsr(1, gnorm(a)));
    1416       32458 :   FF = RgX_Rg_mul(FF, powru(a,k));
    1417       32458 :   GG = RgX_Rg_mul(GG, powru(a,n-k));
    1418             : 
    1419       32458 :   *F = mygprec(FF,bit+n);
    1420       32458 :   *G = mygprec(GG,bit+n); gerepileall(ltop,2, F,G);
    1421       32458 : }
    1422             : 
    1423             : /* split p, this time without scaling. returns in F and G two polynomials
    1424             :  * such that |p-FG|< 2^(-bit)|p| */
    1425             : static void
    1426       35521 : split_2(GEN p, long bit, GEN ctr, double thickness, GEN *F, GEN *G)
    1427             : {
    1428             :   GEN q, FF, GG, R;
    1429             :   double aux, delta, param, param2;
    1430       35521 :   long n = degpol(p), i, j, k, bit2;
    1431             :   double lrmin, lrmax, lrho, *radii;
    1432             : 
    1433       35521 :   radii = (double*) stack_malloc_align((n+1) * sizeof(double), sizeof(double));
    1434             : 
    1435       35521 :   for (i=2; i<n; i++) radii[i] = UNDEF;
    1436       35521 :   aux = thickness/(double)(4 * n);
    1437       35521 :   lrmin = logmin_modulus(p, aux);
    1438       35521 :   lrmax = logmax_modulus(p, aux);
    1439       35521 :   radii[1] = lrmin;
    1440       35521 :   radii[n] = lrmax;
    1441       35521 :   i = 1; j = n;
    1442       35521 :   lrho = (lrmin + lrmax) / 2;
    1443       35521 :   k = dual_modulus(p, lrho, aux, 1);
    1444       35521 :   if (5*k < n || (n < 2*k && 5*k < 4*n))
    1445        7122 :     { lrmax = lrho; j=k+1; radii[j] = lrho; }
    1446             :   else
    1447       28399 :     { lrmin = lrho; i=k;   radii[i] = lrho; }
    1448      112905 :   while (j > i+1)
    1449             :   {
    1450       41863 :     if (i+j == n+1)
    1451       17377 :       lrho = (lrmin + lrmax) / 2;
    1452             :     else
    1453             :     {
    1454       24486 :       double kappa = 2. - log(1. + minss(i,n-j)) / log(1. + minss(j,n-i));
    1455       24486 :       if (i+j < n+1) lrho = lrmax * kappa + lrmin;
    1456       18062 :       else           lrho = lrmin * kappa + lrmax;
    1457       24486 :       lrho /= 1+kappa;
    1458             :     }
    1459       41863 :     aux = (lrmax - lrmin) / (4*(j-i));
    1460       41863 :     k = dual_modulus(p, lrho, aux, minss(i,n+1-j));
    1461       41863 :     if (k-i < j-k-1 || (k-i == j-k-1 && 2*k > n))
    1462       28808 :       { lrmax = lrho; j=k+1; radii[j] = lrho - aux; }
    1463             :     else
    1464       13055 :       { lrmin = lrho; i=k;   radii[i] = lrho + aux; }
    1465             :   }
    1466       35521 :   aux = lrmax - lrmin;
    1467             : 
    1468       35521 :   if (ctr)
    1469             :   {
    1470       32458 :     lrho = (lrmax + lrmin) / 2;
    1471      203255 :     for (i=1; i<=n; i++)
    1472      170797 :       if (radii[i] != UNDEF) radii[i] -= lrho;
    1473             : 
    1474       32458 :     bit2 = bit + (long)(n * fabs(lrho)/LOG2 + 1.);
    1475       32458 :     R = mygprec(dblexp(-lrho), bit2);
    1476       32458 :     q = scalepol(p,R,bit2);
    1477       32458 :     conformal_mapping(radii, ctr, q, k, bit2, aux, &FF, &GG);
    1478             :   }
    1479             :   else
    1480             :   {
    1481        3063 :     lrho = logradius(radii, p, k, aux/10., &delta);
    1482        3063 :     update_radius(n, radii, lrho, &param, &param2);
    1483             : 
    1484        3063 :     bit2 = bit + (long)(n * fabs(lrho)/LOG2 + 1.);
    1485        3063 :     R = mygprec(dblexp(-lrho), bit2);
    1486        3063 :     q = scalepol(p,R,bit2);
    1487        3063 :     optimize_split(q, k, delta, bit2, &FF, &GG, param, param2);
    1488             :   }
    1489       35521 :   bit  += n;
    1490       35521 :   bit2 += n; R = invr(mygprec(R,bit2));
    1491       35521 :   *F = mygprec(scalepol(FF,R,bit2), bit);
    1492       35521 :   *G = mygprec(scalepol(GG,R,bit2), bit);
    1493       35521 : }
    1494             : 
    1495             : /* procedure corresponding to steps 5,6,.. page 44 in RR n. 1852 */
    1496             : /* put in F and G two polynomial such that |p-FG|<2^(-bit)|p|
    1497             :  * where the maximum modulus of the roots of p is <=1.
    1498             :  * Assume sum of roots is 0. */
    1499             : static void
    1500       32458 : split_1(GEN p, long bit, GEN *F, GEN *G)
    1501             : {
    1502       32458 :   long i, imax, n = degpol(p), polreal = isreal(p), ep = gexpo(p), bit2 = bit+n;
    1503             :   GEN ctr, q, qq, FF, GG, v, gr, r, newq;
    1504             :   double lrmin, lrmax, lthick;
    1505       32458 :   const double LOG3 = 1.098613;
    1506             : 
    1507       32458 :   lrmax = logmax_modulus(p, 0.01);
    1508       32458 :   gr = mygprec(dblexp(-lrmax), bit2);
    1509       32458 :   q = scalepol(p,gr,bit2);
    1510             : 
    1511       32458 :   bit2 = bit + gexpo(q) - ep + (long)((double)n*2.*log2(3.)+1);
    1512       32458 :   v = cgetg(5,t_VEC);
    1513       32458 :   gel(v,1) = gen_2;
    1514       32458 :   gel(v,2) = gen_m2;
    1515       32458 :   gel(v,3) = mkcomplex(gen_0, gel(v,1));
    1516       32458 :   gel(v,4) = mkcomplex(gen_0, gel(v,2));
    1517       32458 :   q = mygprec(q,bit2); lthick = 0;
    1518       32458 :   newq = ctr = NULL; /* -Wall */
    1519       32458 :   imax = polreal? 3: 4;
    1520       61935 :   for (i=1; i<=imax; i++)
    1521             :   {
    1522       61033 :     qq = RgX_translate(q, gel(v,i));
    1523       61033 :     lrmin = logmin_modulus(qq,0.05);
    1524       61033 :     if (LOG3 > lrmin + lthick)
    1525             :     {
    1526       59447 :       double lquo = logmax_modulus(qq,0.05) - lrmin;
    1527       59447 :       if (lquo > lthick) { lthick = lquo; newq = qq; ctr = gel(v,i); }
    1528             :     }
    1529       61033 :     if (lthick > LOG2) break;
    1530       34057 :     if (polreal && i==2 && lthick > LOG3 - LOG2) break;
    1531             :   }
    1532       32458 :   bit2 = bit + gexpo(newq) - ep + (long)(n*LOG3/LOG2 + 1);
    1533       32458 :   split_2(newq, bit2, ctr, lthick, &FF, &GG);
    1534       32458 :   r = gneg(mygprec(ctr,bit2));
    1535       32458 :   FF = RgX_translate(FF,r);
    1536       32458 :   GG = RgX_translate(GG,r);
    1537             : 
    1538       32458 :   gr = invr(gr); bit2 = bit - ep + gexpo(FF)+gexpo(GG);
    1539       32458 :   *F = scalepol(FF,gr,bit2);
    1540       32458 :   *G = scalepol(GG,gr,bit2);
    1541       32458 : }
    1542             : 
    1543             : /* put in F and G two polynomials such that |P-FG|<2^(-bit)|P|,
    1544             : where the maximum modulus of the roots of p is < 0.5 */
    1545             : static int
    1546       32593 : split_0_2(GEN p, long bit, GEN *F, GEN *G)
    1547             : {
    1548             :   GEN q, b;
    1549       32593 :   long n = degpol(p), k, bit2, eq;
    1550       32593 :   double aux0 = dbllog2(gel(p,n+2)); /* != -oo */
    1551       32593 :   double aux1 = dbllog2(gel(p,n+1)), aux;
    1552             : 
    1553       32593 :   if (aux1 == -pariINFINITY) /* p1 = 0 */
    1554        1013 :     aux = 0;
    1555             :   else
    1556             :   {
    1557       31580 :     aux = aux1 - aux0; /* log2(p1/p0) */
    1558             :     /* beware double overflow */
    1559       31580 :     if (aux >= 0 && (aux > 1e4 || exp2(aux) > 2.5*n)) return 0;
    1560       31580 :     aux = (aux < -300)? 0.: n*log2(1 + exp2(aux)/(double)n);
    1561             :   }
    1562       32593 :   bit2 = bit+1 + (long)(log2((double)n) + aux);
    1563       32593 :   q = mygprec(p,bit2);
    1564       32593 :   if (aux1 == -pariINFINITY) b = NULL;
    1565             :   else
    1566             :   {
    1567       31580 :     b = gdivgs(gdiv(gel(q,n+1),gel(q,n+2)),-n);
    1568       31580 :     q = RgX_translate(q,b);
    1569             :   }
    1570       32593 :   gel(q,n+1) = gen_0; eq = gexpo(q);
    1571       32593 :   k = 0;
    1572       65551 :   while (k <= n/2 && (- gexpo(gel(q,k+2)) > bit2 + 2*(n-k) + eq
    1573       32836 :                       || gequal0(gel(q,k+2)))) k++;
    1574       32593 :   if (k > 0)
    1575             :   {
    1576         135 :     if (k > n/2) k = n/2;
    1577         135 :     bit2 += k<<1;
    1578         135 :     *F = pol_xn(k, 0);
    1579         135 :     *G = RgX_shift_shallow(q, -k);
    1580             :   }
    1581             :   else
    1582             :   {
    1583       32458 :     split_1(q,bit2,F,G);
    1584       32458 :     bit2 = bit + gexpo(*F) + gexpo(*G) - gexpo(p) + (long)aux+1;
    1585       32458 :     *F = mygprec(*F,bit2);
    1586             :   }
    1587       32593 :   *G = mygprec(*G,bit2);
    1588       32593 :   if (b)
    1589             :   {
    1590       31580 :     GEN mb = mygprec(gneg(b), bit2);
    1591       31580 :     *F = RgX_translate(*F, mb);
    1592       31580 :     *G = RgX_translate(*G, mb);
    1593             :   }
    1594       32593 :   return 1;
    1595             : }
    1596             : 
    1597             : /* put in F and G two polynomials such that |P-FG|<2^(-bit)|P|.
    1598             :  * Assume max_modulus(p) < 2 */
    1599             : static void
    1600       32593 : split_0_1(GEN p, long bit, GEN *F, GEN *G)
    1601             : {
    1602             :   GEN FF, GG;
    1603             :   long n, bit2, normp;
    1604             : 
    1605       65186 :   if  (split_0_2(p,bit,F,G)) return;
    1606             : 
    1607           0 :   normp = gexpo(p);
    1608           0 :   scalepol2n(p,2); /* p := 4^(-n) p(4*x) */
    1609           0 :   n = degpol(p); bit2 = bit + 2*n + gexpo(p) - normp;
    1610           0 :   split_1(mygprec(p,bit2), bit2,&FF,&GG);
    1611           0 :   scalepol2n(FF,-2);
    1612           0 :   scalepol2n(GG,-2); bit2 = bit + gexpo(FF) + gexpo(GG) - normp;
    1613           0 :   *F = mygprec(FF,bit2);
    1614           0 :   *G = mygprec(GG,bit2);
    1615             : }
    1616             : 
    1617             : /* put in F and G two polynomials such that |P-FG|<2^(-bit)|P| */
    1618             : static void
    1619       35656 : split_0(GEN p, long bit, GEN *F, GEN *G)
    1620             : {
    1621       35656 :   const double LOG1_9 = 0.6418539;
    1622       35656 :   long n = degpol(p), k = 0;
    1623             :   GEN q;
    1624             : 
    1625       35656 :   while (gexpo(gel(p,k+2)) < -bit && k <= n/2) k++;
    1626       35656 :   if (k > 0)
    1627             :   {
    1628           0 :     if (k > n/2) k = n/2;
    1629           0 :     *F = pol_xn(k, 0);
    1630           0 :     *G = RgX_shift_shallow(p, -k);
    1631             :   }
    1632             :   else
    1633             :   {
    1634       35656 :     double lr = logmax_modulus(p, 0.05);
    1635       35656 :     if (lr < LOG1_9) split_0_1(p, bit, F, G);
    1636             :     else
    1637             :     {
    1638       27014 :       q = RgX_recip_shallow(p);
    1639       27014 :       lr = logmax_modulus(q,0.05);
    1640       27014 :       if (lr < LOG1_9)
    1641             :       {
    1642       23951 :         split_0_1(q, bit, F, G);
    1643       23951 :         *F = RgX_recip_shallow(*F);
    1644       23951 :         *G = RgX_recip_shallow(*G);
    1645             :       }
    1646             :       else
    1647        3063 :         split_2(p,bit,NULL, 1.2837,F,G);
    1648             :     }
    1649             :   }
    1650       35656 : }
    1651             : 
    1652             : /********************************************************************/
    1653             : /**                                                                **/
    1654             : /**                ERROR ESTIMATE FOR THE ROOTS                    **/
    1655             : /**                                                                **/
    1656             : /********************************************************************/
    1657             : 
    1658             : static GEN
    1659      139128 : root_error(long n, long k, GEN roots_pol, long err, GEN shatzle)
    1660             : {
    1661      139128 :   GEN rho, d, eps, epsbis, eps2, aux, rap = NULL;
    1662             :   long i, j;
    1663             : 
    1664      139128 :   d = cgetg(n+1,t_VEC);
    1665     1635368 :   for (i=1; i<=n; i++)
    1666             :   {
    1667     1496240 :     if (i!=k)
    1668             :     {
    1669     1357112 :       aux = gsub(gel(roots_pol,i), gel(roots_pol,k));
    1670     1357112 :       gel(d,i) = gabs(mygprec(aux,31), DEFAULTPREC);
    1671             :     }
    1672             :   }
    1673      139128 :   rho = gabs(mygprec(gel(roots_pol,k),31), DEFAULTPREC);
    1674      139128 :   if (expo(rho) < 0) rho = real_1(DEFAULTPREC);
    1675      139128 :   eps = mulrr(rho, shatzle);
    1676      139128 :   aux = shiftr(powru(rho,n), err);
    1677             : 
    1678      434217 :   for (j=1; j<=2 || (j<=5 && cmprr(rap, dbltor(1.2)) > 0); j++)
    1679             :   {
    1680      295089 :     GEN prod = NULL; /* 1. */
    1681      295089 :     long m = n;
    1682      295089 :     epsbis = mulrr(eps, dbltor(1.25));
    1683     3962843 :     for (i=1; i<=n; i++)
    1684             :     {
    1685     3667754 :       if (i != k && cmprr(gel(d,i),epsbis) > 0)
    1686             :       {
    1687     3341831 :         GEN dif = subrr(gel(d,i),eps);
    1688     3341831 :         prod = prod? mulrr(prod, dif): dif;
    1689     3341831 :         m--;
    1690             :       }
    1691             :     }
    1692      295089 :     eps2 = prod? divrr(aux, prod): aux;
    1693      295089 :     if (m > 1) eps2 = sqrtnr(shiftr(eps2, 2*m-2), m);
    1694      295089 :     rap = divrr(eps,eps2); eps = eps2;
    1695             :   }
    1696      139128 :   return eps;
    1697             : }
    1698             : 
    1699             : /* round a complex or real number x to an absolute value of 2^(-bit) */
    1700             : static GEN
    1701      337272 : mygprec_absolute(GEN x, long bit)
    1702             : {
    1703             :   long e;
    1704             :   GEN y;
    1705             : 
    1706      337272 :   switch(typ(x))
    1707             :   {
    1708             :     case t_REAL:
    1709      221486 :       e = expo(x) + bit;
    1710      221486 :       return (e <= 0 || !signe(x))? real_0_bit(-bit): rtor(x, nbits2prec(e));
    1711             :     case t_COMPLEX:
    1712      100474 :       if (gexpo(gel(x,2)) < -bit) return mygprec_absolute(gel(x,1),bit);
    1713       97670 :       y = cgetg(3,t_COMPLEX);
    1714       97670 :       gel(y,1) = mygprec_absolute(gel(x,1),bit);
    1715       97670 :       gel(y,2) = mygprec_absolute(gel(x,2),bit);
    1716       97670 :       return y;
    1717       15312 :     default: return x;
    1718             :   }
    1719             : }
    1720             : 
    1721             : static long
    1722       31896 : a_posteriori_errors(GEN p, GEN roots_pol, long err)
    1723             : {
    1724       31896 :   long i, n = degpol(p), e_max = -(long)EXPOBITS;
    1725             :   GEN sigma, shatzle;
    1726             : 
    1727       31896 :   err += (long)log2((double)n) + 1;
    1728       31896 :   if (err > -2) return 0;
    1729       31896 :   sigma = real2n(-err, LOWDEFAULTPREC);
    1730             :   /*  2 / ((s - 1)^(1/n) - 1) */
    1731       31896 :   shatzle = divur(2, subrs(sqrtnr(subrs(sigma,1),n), 1));
    1732      171024 :   for (i=1; i<=n; i++)
    1733             :   {
    1734      139128 :     pari_sp av = avma;
    1735      139128 :     GEN x = root_error(n,i,roots_pol,err,shatzle);
    1736      139128 :     long e = gexpo(x);
    1737      139128 :     avma = av; if (e > e_max) e_max = e;
    1738      139128 :     gel(roots_pol,i) = mygprec_absolute(gel(roots_pol,i), -e);
    1739             :   }
    1740       31896 :   return e_max;
    1741             : }
    1742             : 
    1743             : /********************************************************************/
    1744             : /**                                                                **/
    1745             : /**                           MAIN                                 **/
    1746             : /**                                                                **/
    1747             : /********************************************************************/
    1748             : static GEN
    1749      109103 : append_clone(GEN r, GEN a) { a = gclone(a); vectrunc_append(r, a); return a; }
    1750             : 
    1751             : /* put roots in placeholder roots_pol so that |P - L_1...L_n| < 2^(-bit)|P|
    1752             :  * returns prod (x-roots_pol[i]) */
    1753             : static GEN
    1754      103208 : split_complete(GEN p, long bit, GEN roots_pol)
    1755             : {
    1756      103208 :   long n = degpol(p);
    1757             :   pari_sp ltop;
    1758             :   GEN p1, F, G, a, b, m1, m2;
    1759             : 
    1760      103208 :   if (n == 1)
    1761             :   {
    1762       26001 :     a = gneg_i(gdiv(gel(p,2), gel(p,3)));
    1763       26001 :     (void)append_clone(roots_pol,a); return p;
    1764             :   }
    1765       77207 :   ltop = avma;
    1766       77207 :   if (n == 2)
    1767             :   {
    1768       41551 :     F = gsub(gsqr(gel(p,3)), gmul2n(gmul(gel(p,2),gel(p,4)), 2));
    1769       41551 :     F = gsqrt(F, nbits2prec(bit));
    1770       41551 :     p1 = ginv(gmul2n(gel(p,4),1));
    1771       41551 :     a = gneg_i(gmul(gadd(F,gel(p,3)), p1));
    1772       41551 :     b =        gmul(gsub(F,gel(p,3)), p1);
    1773       41551 :     a = append_clone(roots_pol,a);
    1774       41551 :     b = append_clone(roots_pol,b); avma = ltop;
    1775       41551 :     a = mygprec(a, 3*bit);
    1776       41551 :     b = mygprec(b, 3*bit);
    1777       41551 :     return gmul(gel(p,4), mkpoln(3, gen_1, gneg(gadd(a,b)), gmul(a,b)));
    1778             :   }
    1779       35656 :   split_0(p,bit,&F,&G);
    1780       35656 :   m1 = split_complete(F,bit,roots_pol);
    1781       35656 :   m2 = split_complete(G,bit,roots_pol);
    1782       35656 :   return gerepileupto(ltop, gmul(m1,m2));
    1783             : }
    1784             : 
    1785             : static GEN
    1786      500432 : quicktofp(GEN x)
    1787             : {
    1788      500432 :   const long prec = DEFAULTPREC;
    1789      500432 :   switch(typ(x))
    1790             :   {
    1791      490856 :     case t_INT: return itor(x, prec);
    1792        4095 :     case t_REAL: return rtor(x, prec);
    1793           0 :     case t_FRAC: return fractor(x, prec);
    1794             :     case t_COMPLEX: {
    1795        5481 :       GEN a = gel(x,1), b = gel(x,2);
    1796             :       /* avoid problem with 0, e.g. x = 0 + I*1e-100. We don't want |x| = 0. */
    1797        5481 :       if (isintzero(a)) return cxcompotor(b, prec);
    1798        5481 :       if (isintzero(b)) return cxcompotor(a, prec);
    1799        5481 :       a = cxcompotor(a, prec);
    1800        5481 :       b = cxcompotor(b, prec); return sqrtr(addrr(sqrr(a), sqrr(b)));
    1801             :     }
    1802           0 :     default: pari_err_TYPE("quicktofp",x);
    1803             :       return NULL;/*LCOV_EXCL_LINE*/
    1804             :   }
    1805             : 
    1806             : }
    1807             : 
    1808             : /* bound log_2 |largest root of p| (Fujiwara's bound) */
    1809             : double
    1810      120171 : fujiwara_bound(GEN p)
    1811             : {
    1812      120171 :   pari_sp av = avma;
    1813      120171 :   long i, n = degpol(p);
    1814             :   GEN cc;
    1815             :   double loglc, Lmax;
    1816             : 
    1817      120171 :   if (n <= 0) pari_err_CONSTPOL("fujiwara_bound");
    1818      120171 :   loglc = mydbllog2r( quicktofp(gel(p,n+2)) ); /* log_2 |lc(p)| */
    1819      120171 :   cc = gel(p, 2);
    1820      120171 :   if (gequal0(cc))
    1821       16962 :     Lmax = -pariINFINITY-1;
    1822             :   else
    1823      103209 :     Lmax = (mydbllog2r(quicktofp(cc)) - loglc - 1) / n;
    1824      478586 :   for (i = 1; i < n; i++)
    1825             :   {
    1826      358415 :     GEN y = gel(p,i+2);
    1827             :     double L;
    1828      358415 :     if (gequal0(y)) continue;
    1829      277052 :     L = (mydbllog2r(quicktofp(y)) - loglc) / (n-i);
    1830      277052 :     if (L > Lmax) Lmax = L;
    1831             :   }
    1832      120171 :   avma = av; return Lmax + 1;
    1833             : }
    1834             : 
    1835             : /* Fujiwara's bound, real roots. Based on the following remark: if
    1836             :  *   p = x^n + sum a_i x^i and q = x^n + sum min(a_i,0)x^i
    1837             :  * then for all x >= 0, p(x) >= q(x). Thus any bound for the (positive) roots
    1838             :  * of q is a bound for the positive roots of p. */
    1839             : double
    1840       24114 : fujiwara_bound_real(GEN p, long sign)
    1841             : {
    1842       24114 :   pari_sp av = avma;
    1843             :   GEN x;
    1844       24114 :   long n = degpol(p), i, signodd, signeven;
    1845             :   double fb;
    1846       24114 :   if (n <= 0) pari_err_CONSTPOL("fujiwara_bound");
    1847       24114 :   x = shallowcopy(p);
    1848       24114 :   if (gsigne(gel(x, n+2)) > 0)
    1849             :   {
    1850       24114 :     signeven = 1;
    1851       24114 :     signodd = sign;
    1852             :   }
    1853             :   else
    1854             :   {
    1855           0 :     signeven = -1;
    1856           0 :     signodd = -sign;
    1857             :   }
    1858      119072 :   for (i = 0; i < n; i++)
    1859             :   {
    1860       94958 :     if ((n - i) % 2)
    1861             :     {
    1862       49823 :       if (gsigne(gel(x, i+2)) == signodd ) gel(x, i+2) = gen_0;
    1863             :     }
    1864             :     else
    1865             :     {
    1866       45135 :       if (gsigne(gel(x, i+2)) == signeven) gel(x, i+2) = gen_0;
    1867             :     }
    1868             :   }
    1869       24114 :   fb = fujiwara_bound(x);
    1870       24114 :   avma = av; return fb;
    1871             : }
    1872             : 
    1873             : static GEN
    1874      182602 : mygprecrc_special(GEN x, long prec, long e)
    1875             : {
    1876             :   GEN y;
    1877      182602 :   switch(typ(x))
    1878             :   {
    1879             :     case t_REAL:
    1880       16259 :       if (!signe(x)) return real_0_bit(minss(e, expo(x)));
    1881       15838 :       return (prec > realprec(x))? rtor(x, prec): x;
    1882             :     case t_COMPLEX:
    1883        5789 :       y = cgetg(3,t_COMPLEX);
    1884        5789 :       gel(y,1) = mygprecrc_special(gel(x,1),prec,e);
    1885        5789 :       gel(y,2) = mygprecrc_special(gel(x,2),prec,e);
    1886        5789 :       return y;
    1887      160554 :     default: return x;
    1888             :   }
    1889             : }
    1890             : 
    1891             : /* like mygprec but keep at least the same precision as before */
    1892             : static GEN
    1893       31896 : mygprec_special(GEN x, long bit)
    1894             : {
    1895             :   long lx, i, e, prec;
    1896             :   GEN y;
    1897             : 
    1898       31896 :   if (bit < 0) bit = 0; /* should not happen */
    1899       31896 :   e = gexpo(x) - bit;
    1900       31896 :   prec = nbits2prec(bit);
    1901       31896 :   switch(typ(x))
    1902             :   {
    1903             :     case t_POL:
    1904       31896 :       y = cgetg_copy(x, &lx); y[1] = x[1];
    1905       31896 :       for (i=2; i<lx; i++) gel(y,i) = mygprecrc_special(gel(x,i),prec,e);
    1906       31896 :       break;
    1907             : 
    1908           0 :     default: y = mygprecrc_special(x,prec,e);
    1909             :   }
    1910       31896 :   return y;
    1911             : }
    1912             : 
    1913             : static GEN
    1914       20655 : fix_roots1(GEN r)
    1915             : {
    1916       20655 :   long i, l = lg(r);
    1917       20655 :   GEN allr = cgetg(l, t_VEC);
    1918      107057 :   for (i=1; i<l; i++)
    1919             :   {
    1920       86402 :     GEN t = gel(r,i);
    1921       86402 :     gel(allr,i) = gcopy(t); gunclone(t);
    1922             :   }
    1923       20655 :   return allr;
    1924             : }
    1925             : static GEN
    1926       31896 : fix_roots(GEN r, GEN *m, long h, long bit)
    1927             : {
    1928             :   long i, j, k, l, prec;
    1929             :   GEN allr, ro1;
    1930       31896 :   if (h == 1) return fix_roots1(r);
    1931       11241 :   prec = nbits2prec(bit);
    1932       11241 :   ro1 = grootsof1(h, prec) + 1;
    1933       11241 :   l = lg(r)-1;
    1934       11241 :   allr = cgetg(h*l+1, t_VEC);
    1935       33942 :   for (k=1,i=1; i<=l; i++)
    1936             :   {
    1937       22701 :     GEN p2, p1 = gel(r,i);
    1938       22701 :     p2 = (h == 2)? gsqrt(p1, prec): gsqrtn(p1, utoipos(h), NULL, prec);
    1939       22701 :     for (j=0; j<h; j++) gel(allr,k++) = gmul(p2, gel(ro1,j));
    1940       22701 :     gunclone(p1);
    1941             :   }
    1942       11241 :   *m = roots_to_pol(allr, 0);
    1943       11241 :   return allr;
    1944             : }
    1945             : 
    1946             : static GEN
    1947       31542 : all_roots(GEN p, long bit)
    1948             : {
    1949             :   GEN lc, pd, q, roots_pol, m;
    1950       31542 :   long bit0,  bit2, i, e, h, n = degpol(p);
    1951             :   pari_sp av;
    1952             : 
    1953       31542 :   pd = RgX_deflate_max(p, &h); lc = leading_coeff(pd);
    1954       31542 :   e = (long)(2 * fujiwara_bound(pd)); if (e < 0) e = 0;
    1955       31542 :   bit0 = bit + gexpo(pd) - gexpo(lc) + (long)log2(n/h)+1+e;
    1956       31542 :   bit2 = bit0; e = 0;
    1957       31896 :   for (av=avma,i=1;; i++,avma=av)
    1958             :   {
    1959       31896 :     roots_pol = vectrunc_init(n+1);
    1960       31896 :     bit2 += e + (n << i);
    1961       31896 :     q = RgX_gtofp_bit(mygprec(pd,bit2), bit2);
    1962       31896 :     q[1] = evalsigne(1)|evalvarn(0);
    1963       31896 :     m = split_complete(q,bit2,roots_pol);
    1964       31896 :     roots_pol = fix_roots(roots_pol, &m, h, bit2);
    1965       31896 :     q = mygprec_special(p,bit2); lc = leading_coeff(q);
    1966       31896 :     q[1] = evalsigne(1)|evalvarn(0);
    1967       31896 :     if (h > 1) m = gmul(m,lc);
    1968             : 
    1969       31896 :     e = gexpo(gsub(q, m)) - gexpo(lc) + (long)log2((double)n) + 1;
    1970       31896 :     if (e < -2*bit2) e = -2*bit2; /* avoid e = -pariINFINITY */
    1971       31896 :     if (e < 0)
    1972             :     {
    1973       31896 :       e = bit + a_posteriori_errors(p,roots_pol,e);
    1974       63438 :       if (e < 0) return roots_pol;
    1975             :     }
    1976         354 :     if (DEBUGLEVEL > 7)
    1977           0 :       err_printf("all_roots: restarting, i = %ld, e = %ld\n", i,e);
    1978         354 :   }
    1979             : }
    1980             : 
    1981             : INLINE int
    1982        7585 : isexactscalar(GEN x) { long tx = typ(x); return is_rational_t(tx); }
    1983             : 
    1984             : static int
    1985        4621 : isexactpol(GEN p)
    1986             : {
    1987        4621 :   long i,n = degpol(p);
    1988        8475 :   for (i=0; i<=n; i++)
    1989        7585 :     if (!isexactscalar(gel(p,i+2))) return 0;
    1990         890 :   return 1;
    1991             : }
    1992             : 
    1993             : static GEN
    1994         890 : solve_exact_pol(GEN p, long bit)
    1995             : {
    1996         890 :   long i, j, k, m, n = degpol(p), iroots = 0;
    1997         890 :   GEN ex, factors, v = zerovec(n);
    1998             : 
    1999         890 :   factors = ZX_squff(Q_primpart(p), &ex);
    2000        1780 :   for (i=1; i<lg(factors); i++)
    2001             :   {
    2002         890 :     GEN roots_fact = all_roots(gel(factors,i), bit);
    2003         890 :     n = degpol(gel(factors,i));
    2004         890 :     m = ex[i];
    2005        3546 :     for (j=1; j<=n; j++)
    2006        2656 :       for (k=1; k<=m; k++) v[++iroots] = roots_fact[j];
    2007             :   }
    2008         890 :   return v;
    2009             : }
    2010             : 
    2011             : /* return the roots of p with absolute error bit */
    2012             : static GEN
    2013        4621 : roots_com(GEN q, long bit)
    2014             : {
    2015             :   GEN L, p;
    2016        4621 :   long v = RgX_valrem_inexact(q, &p);
    2017        4621 :   int ex = isexactpol(p);
    2018        4621 :   if (!ex) p = RgX_normalize1(p);
    2019        4621 :   if (lg(p) == 3)
    2020           7 :     L = cgetg(1,t_VEC); /* constant polynomial */
    2021             :   else
    2022        4614 :     L = ex? solve_exact_pol(p,bit): all_roots(p,bit);
    2023        4621 :   if (v)
    2024             :   {
    2025         112 :     GEN M, z, t = gel(q,2);
    2026             :     long i, x, y, l, n;
    2027             : 
    2028         112 :     if (isrationalzero(t)) x = -bit;
    2029             :     else
    2030             :     {
    2031          14 :       n = gexpo(t);
    2032          14 :       x = n / v; l = degpol(q);
    2033          56 :       for (i = v; i <= l; i++)
    2034             :       {
    2035          42 :         t  = gel(q,i+2);
    2036          42 :         if (isrationalzero(t)) continue;
    2037          42 :         y = (n - gexpo(t)) / i;
    2038          42 :         if (y < x) x = y;
    2039             :       }
    2040             :     }
    2041         112 :     z = real_0_bit(x); l = v + lg(L);
    2042         112 :     M = cgetg(l, t_VEC); L -= v;
    2043         112 :     for (i = 1; i <= v; i++) gel(M,i) = z;
    2044         112 :     for (     ; i <  l; i++) gel(M,i) = gel(L,i);
    2045         112 :     L = M;
    2046             :   }
    2047        4621 :   return L;
    2048             : }
    2049             : 
    2050             : static GEN
    2051      102019 : tocomplex(GEN x, long l, long bit)
    2052             : {
    2053             :   GEN y;
    2054      102019 :   if (typ(x) == t_COMPLEX)
    2055             :   {
    2056       96556 :     if (signe(gel(x,1))) return mygprecrc(x, l, -bit);
    2057       15152 :     x = gel(x,2);
    2058       15152 :     y = cgetg(3,t_COMPLEX);
    2059       15152 :     gel(y,1) = real_0_bit(-bit);
    2060       15152 :     gel(y,2) = mygprecrc(x, l, -bit);
    2061             :   }
    2062             :   else
    2063             :   {
    2064        5463 :     y = cgetg(3,t_COMPLEX);
    2065        5463 :     gel(y,1) = mygprecrc(x, l, -bit);
    2066        5463 :     gel(y,2) = real_0_bit(-bit);
    2067             :   }
    2068       20615 :   return y;
    2069             : }
    2070             : 
    2071             : /* x,y are t_COMPLEX of t_REALs or t_REAL, compare lexicographically,
    2072             :  * up to 2^-e absolute error */
    2073             : static int
    2074      226310 : cmp_complex_appr(void *E, GEN x, GEN y)
    2075             : {
    2076      226310 :   long e = (long)E;
    2077             :   GEN z, xi, yi, xr, yr;
    2078             :   long sxi, syi;
    2079      226310 :   if (typ(x) == t_COMPLEX) { xr = gel(x,1); xi = gel(x,2); sxi = signe(xi); }
    2080       79466 :   else { xr = x; xi = NULL; sxi = 0; }
    2081      226310 :   if (typ(y) == t_COMPLEX) { yr = gel(y,1); yi = gel(y,2); syi = signe(yi); }
    2082       72202 :   else { yr = y; yi = NULL; syi = 0; }
    2083             :   /* Compare absolute values of imaginary parts */
    2084      226310 :   if (!sxi)
    2085             :   {
    2086       85201 :     if (syi && expo(yi) >= e) return -1;
    2087             :     /* |Im x| ~ |Im y| ~ 0 */
    2088             :   }
    2089      141109 :   else if (!syi)
    2090             :   {
    2091        3587 :     if (sxi && expo(xi) >= e) return 1;
    2092             :     /* |Im x| ~ |Im y| ~ 0 */
    2093             :   }
    2094             :   else
    2095             :   {
    2096             :     long sz;
    2097      137522 :     z = addrr_sign(xi, 1, yi, -1);
    2098      137522 :     sz = signe(z);
    2099      137522 :     if (sz && expo(z) >= e) return (int)sz;
    2100             :   }
    2101             :   /* |Im x| ~ |Im y|, sort according to real parts */
    2102      133921 :   z = subrr(xr, yr);
    2103      133921 :   if (expo(z) >= e) return (int)signe(z);
    2104             :   /* Re x ~ Re y. Place negative absolute value before positive */
    2105       41714 :   return (int) (sxi - syi);
    2106             : }
    2107             : 
    2108             : static GEN
    2109       31549 : clean_roots(GEN L, long l, long bit, long clean)
    2110             : {
    2111       31549 :   long i, n = lg(L), ex = 5 - bit;
    2112       31549 :   GEN res = cgetg(n,t_COL);
    2113      168419 :   for (i=1; i<n; i++)
    2114             :   {
    2115      136870 :     GEN c = gel(L,i);
    2116      136870 :     if (clean && isrealappr(c,ex))
    2117             :     {
    2118       34851 :       if (typ(c) == t_COMPLEX) c = gel(c,1);
    2119       34851 :       c = mygprecrc(c, l, -bit);
    2120             :     }
    2121             :     else
    2122      102019 :       c = tocomplex(c, l, bit);
    2123      136870 :     gel(res,i) = c;
    2124             :   }
    2125       31549 :   gen_sort_inplace(res, (void*)ex, &cmp_complex_appr, NULL);
    2126       31549 :   return res;
    2127             : }
    2128             : 
    2129             : /* the vector of roots of p, with absolute error 2^(- prec2nbits(l)) */
    2130             : static GEN
    2131        4649 : roots_aux(GEN p, long l, long clean)
    2132             : {
    2133        4649 :   pari_sp av = avma;
    2134             :   long bit;
    2135             :   GEN L;
    2136             : 
    2137        4649 :   if (typ(p) != t_POL)
    2138             :   {
    2139          21 :     if (gequal0(p)) pari_err_ROOTS0("roots");
    2140          14 :     if (!isvalidcoeff(p)) pari_err_TYPE("roots",p);
    2141           7 :     return cgetg(1,t_COL); /* constant polynomial */
    2142             :   }
    2143        4628 :   if (!signe(p)) pari_err_ROOTS0("roots");
    2144        4628 :   checkvalidpol(p,"roots");
    2145        4621 :   if (lg(p) == 3) return cgetg(1,t_COL); /* constant polynomial */
    2146        4621 :   if (l < LOWDEFAULTPREC) l = LOWDEFAULTPREC;
    2147        4621 :   bit = prec2nbits(l);
    2148        4621 :   L = roots_com(p, bit);
    2149        4621 :   return gerepileupto(av, clean_roots(L, l, bit, clean));
    2150             : }
    2151             : GEN
    2152        4502 : roots(GEN p, long l) { return roots_aux(p,l, 0); }
    2153             : /* clean up roots. If root is real replace it by its real part */
    2154             : GEN
    2155         147 : cleanroots(GEN p, long l) { return roots_aux(p,l, 1); }
    2156             : 
    2157             : /* private variant of conjvec. Allow non rational coefficients, shallow
    2158             :  * function. */
    2159             : GEN
    2160          77 : polmod_to_embed(GEN x, long prec)
    2161             : {
    2162          77 :   GEN v, T = gel(x,1), A = gel(x,2);
    2163             :   long i, l;
    2164          77 :   if (typ(A) != t_POL || varn(A) != varn(T))
    2165             :   {
    2166           0 :     checkvalidpol(T,"polmod_to_embed");
    2167           0 :     return const_col(degpol(T), A);
    2168             :   }
    2169          77 :   v = cleanroots(T,prec); l = lg(v);
    2170          77 :   for (i=1; i<l; i++) gel(v,i) = poleval(A,gel(v,i));
    2171          77 :   return v;
    2172             : }
    2173             : 
    2174             : GEN
    2175       26928 : QX_complex_roots(GEN p, long l)
    2176             : {
    2177       26928 :   pari_sp av = avma;
    2178             :   long bit;
    2179             :   GEN L;
    2180             : 
    2181       26928 :   if (!signe(p)) pari_err_ROOTS0("QX_complex_roots");
    2182       26928 :   if (lg(p) == 3) return cgetg(1,t_COL); /* constant polynomial */
    2183       26928 :   if (l < LOWDEFAULTPREC) l = LOWDEFAULTPREC;
    2184       26928 :   bit = prec2nbits(l);
    2185       26928 :   L = all_roots(Q_primpart(p), bit);
    2186       26928 :   return gerepileupto(av, clean_roots(L, l, bit, 1));
    2187             : }
    2188             : 
    2189             : /********************************************************************/
    2190             : /**                                                                **/
    2191             : /**                REAL ROOTS OF INTEGER POLYNOMIAL                **/
    2192             : /**                                                                **/
    2193             : /********************************************************************/
    2194             : 
    2195             : /* Count sign changes in the coefficients of (x+1)^deg(P)*P(1/(x+1))
    2196             :  * The inversion is implicit (we take coefficients backwards). Roots of P
    2197             :  * at 0 and 1 (mapped to oo and 0) are ignored here and must be dealt with
    2198             :  * by the caller */
    2199             : static long
    2200      437258 : X2XP1(GEN P, long deg, int *root1, GEN *Premapped)
    2201             : {
    2202      437258 :   const pari_sp av = avma;
    2203      437258 :   GEN v = shallowcopy(P);
    2204             :   long i, j, vlim, nb, s;
    2205             : 
    2206      437258 :   for (i = 0, vlim = deg+2;;)
    2207             :   {
    2208      437328 :     for (j = 2; j < vlim; j++) gel(v, j+1) = addii(gel(v, j), gel(v, j+1));
    2209      437328 :     s = -signe(gel(v, vlim));
    2210      437328 :     vlim--; i++; if (s) break;
    2211          70 :   }
    2212      437258 :   if (vlim == deg+1) *root1 = 0;
    2213             :   else
    2214             :   {
    2215          70 :     *root1 = 1;
    2216          70 :     if (Premapped) setlg(v, vlim + 2);
    2217             :   }
    2218             : 
    2219      437258 :   nb = 0;
    2220     1633813 :   for (; i < deg; i++)
    2221             :   {
    2222     1566651 :     long s2 = -signe(gel(v, 2));
    2223     1566651 :     int flag = (s2 == s);
    2224    20558138 :     for (j = 2; j < vlim; j++)
    2225             :     {
    2226    18991487 :       gel(v, j+1) = addii(gel(v, j), gel(v, j+1));
    2227    18991487 :       if (flag) flag = (s2 != signe(gel(v, j+1)));
    2228             :     }
    2229     1566651 :     if (s == signe(gel(v, vlim)))
    2230             :     {
    2231      414693 :       if (++nb >= 2) { avma = av; return 2; }
    2232      240494 :       s = -s;
    2233             :     }
    2234             :     /* if flag is set there will be no further sign changes */
    2235     1392452 :     if (flag && (!Premapped || !nb)) goto END;
    2236     1196555 :     vlim--;
    2237     1196555 :     if (gc_needed(av, 3))
    2238             :     {
    2239           0 :       if (!Premapped) setlg(v, vlim + 2);
    2240           0 :       v = gerepileupto(av, v);
    2241           0 :       if (DEBUGMEM > 1) pari_warn(warnmem, "X2XP1");
    2242             :     }
    2243             :   }
    2244       67162 :   if (vlim >= 2 && s == signe(gel(v, vlim))) nb++;
    2245             : END:
    2246      263059 :   if (Premapped && nb == 1) *Premapped = v; else avma = av;
    2247      263059 :   return nb;
    2248             : }
    2249             : 
    2250             : static long
    2251           0 : _intervalcmp(GEN x, GEN y)
    2252             : {
    2253           0 :   if (typ(x) == t_VEC) x = gel(x, 1);
    2254           0 :   if (typ(y) == t_VEC) y = gel(y, 1);
    2255           0 :   return gcmp(x, y);
    2256             : }
    2257             : 
    2258             : static GEN
    2259     2332716 : _gen_nored(void *E, GEN x) { (void)E; return x; }
    2260             : static GEN
    2261    13284888 : _mp_add(void *E, GEN x, GEN y) { (void)E; return mpadd(x, y); }
    2262             : static GEN
    2263           0 : _mp_sub(void *E, GEN x, GEN y) { (void)E; return mpsub(x, y); }
    2264             : static GEN
    2265     2968144 : _mp_mul(void *E, GEN x, GEN y) { (void)E; return mpmul(x, y); }
    2266             : static GEN
    2267     2833704 : _mp_sqr(void *E, GEN x) { (void)E; return mpsqr(x); }
    2268             : static GEN
    2269     2666072 : _gen_one(void *E) { (void)E; return gen_1; }
    2270             : static GEN
    2271       28879 : _gen_zero(void *E) { (void)E; return gen_0; }
    2272             : 
    2273             : static struct bb_algebra mp_algebra = { _gen_nored, _mp_add, _mp_sub,
    2274             :                          _mp_mul, _mp_sqr, _gen_one, _gen_zero };
    2275             : 
    2276             : static GEN
    2277    14919361 : _mp_cmul(void *E, GEN P, long a, GEN x) {(void)E; return mpmul(gel(P,a+2), x);}
    2278             : 
    2279             : /* Split the polynom P in two parts, whose coeffs have constant sign:
    2280             :  * P(X) = X^D*Pp + Pm. Also compute the two parts of the derivative of P,
    2281             :  * Pprimem = Pm', Pprimep = X*Pp'+ D*Pp => P' = X^(D-1)*Pprimep + Pprimem;
    2282             :  * Pprimep[i] = (i+D) Pp[i]. Return D */
    2283             : static long
    2284       28538 : split_pols(GEN P, GEN *pPp, GEN *pPm, GEN *pPprimep, GEN *pPprimem)
    2285             : {
    2286       28538 :   long i, D, dP = degpol(P), s0 = signe(gel(P,2));
    2287             :   GEN Pp, Pm, Pprimep, Pprimem;
    2288      194840 :   for(i=1; i <= dP; i++)
    2289      194840 :     if (signe(gel(P, i+2)) == -s0) break;
    2290       28538 :   D = i;
    2291       28538 :   Pm = cgetg(D + 2, t_POL);
    2292       28538 :   Pprimem = cgetg(D + 1, t_POL);
    2293       28538 :   Pp = cgetg(dP-D + 3, t_POL);
    2294       28538 :   Pprimep = cgetg(dP-D + 3, t_POL);
    2295       28538 :   Pm[1] = Pp[1] = Pprimem[1] = Pprimep[1] = P[1];
    2296      223378 :   for(i=0; i < D; i++)
    2297             :   {
    2298      194840 :     GEN c = gel(P, i+2);
    2299      194840 :     gel(Pm, i+2) = c;
    2300      194840 :     if (i) gel(Pprimem, i+1) = mului(i, c);
    2301             :   }
    2302      235870 :   for(; i <= dP; i++)
    2303             :   {
    2304      207332 :     GEN c = gel(P, i+2);
    2305      207332 :     gel(Pp, i+2-D) = c;
    2306      207332 :     gel(Pprimep, i+2-D) = mului(i, c);
    2307             :   }
    2308       28538 :   *pPm = normalizepol_lg(Pm, D+2);
    2309       28538 :   *pPprimem = normalizepol_lg(Pprimem, D+1);
    2310       28538 :   *pPp = normalizepol_lg(Pp, dP-D+3);
    2311       28538 :   *pPprimep = normalizepol_lg(Pprimep, dP-D+3);
    2312       28538 :   return dP - degpol(*pPp);
    2313             : }
    2314             : 
    2315             : static GEN
    2316      831676 : bkeval_single_power(long d, GEN V)
    2317             : {
    2318      831676 :   long mp = lg(V) - 2;
    2319      831676 :   if (d > mp) return gmul(gpowgs(gel(V, mp+1), d/mp), gel(V, (d%mp)+1));
    2320      831676 :   return gel(V, d+1);
    2321             : }
    2322             : 
    2323             : static GEN
    2324      831676 : splitpoleval(GEN Pp, GEN Pm, GEN pows, long D, long bitprec)
    2325             : {
    2326      831676 :   GEN vp = gen_bkeval_powers(Pp, degpol(Pp), pows, NULL, &mp_algebra, _mp_cmul);
    2327      831676 :   GEN vm = gen_bkeval_powers(Pm, degpol(Pm), pows, NULL, &mp_algebra, _mp_cmul);
    2328      831676 :   GEN xa = bkeval_single_power(D, pows);
    2329             :   GEN r;
    2330      831676 :   if (!signe(vp)) return vm;
    2331      831676 :   vp = gmul(vp, xa);
    2332      831676 :   r = gadd(vp, vm);
    2333      831676 :   if (gexpo(vp) - (signe(r)? gexpo(r): 0) > prec2nbits(realprec(vp)) - bitprec)
    2334       42403 :     return NULL;
    2335      789273 :   return r;
    2336             : }
    2337             : 
    2338             : /* optimized Cauchy bound for P = X^D*Pp + Pm, D > deg(Pm) */
    2339             : static GEN
    2340       28538 : splitcauchy(GEN Pp, GEN Pm, long prec)
    2341             : {
    2342       28538 :   GEN S = gel(Pp,2), A = gel(Pm,2);
    2343       28538 :   long i, lPm = lg(Pm), lPp = lg(Pp);
    2344       28538 :   for (i=3; i < lPm; i++) { GEN c = gel(Pm,i); if (abscmpii(A, c) < 0) A = c; }
    2345       28538 :   for (i=3; i < lPp; i++) S = addii(S, gel(Pp, i));
    2346       28538 :   return subsr(1, rdivii(A, S, prec)); /* 1 + |Pm|_oo / |Pp|_1 */
    2347             : }
    2348             : 
    2349             : /* Newton for polynom P, P(0)!=0, with unique sign change => one root in ]0,oo[
    2350             :  * P' has also at most one zero there */
    2351             : static GEN
    2352       28538 : polsolve(GEN P, long bitprec)
    2353             : {
    2354       28538 :   pari_sp av = avma, av2;
    2355             :   GEN Pp, Pm, Pprimep, Pprimem, Pprime, Pprime2, ra, rb, rc, Pc;
    2356       28538 :   long deg = degpol(P);
    2357       28538 :   long expoold = LONG_MAX, cprec = DEFAULTPREC, prec = nbits2prec(bitprec);
    2358             :   long iter, D, rt, s0, bitaddprec, addprec;
    2359       28538 :   if (deg == 1)
    2360           0 :     return gerepileuptoleaf(av, rdivii(negi(gel(P,2)), gel(P,3), prec));
    2361       28538 :   Pprime = ZX_deriv(P);
    2362       28538 :   Pprime2 = ZX_deriv(Pprime);
    2363       28538 :   bitaddprec = 1 + 2*expu(deg); addprec = nbits2prec(bitaddprec);
    2364       28538 :   D = split_pols(P, &Pp, &Pm, &Pprimep, &Pprimem); /* P = X^D*Pp + Pm */
    2365       28538 :   s0 = signe(gel(P, 2));
    2366       28538 :   rt = maxss(D, brent_kung_optpow(maxss(degpol(Pp), degpol(Pm)), 2, 1));
    2367       28538 :   rb = splitcauchy(Pp, Pm, DEFAULTPREC);
    2368             :   for(;;)
    2369             :   {
    2370       28538 :     GEN pows = gen_powers(rb, rt, 1, NULL, _mp_sqr, _mp_mul, _gen_one);
    2371       28538 :     Pc = splitpoleval(Pp, Pm, pows, D, bitaddprec);
    2372       28538 :     if (!Pc) { cprec++; rb = rtor(rb, cprec); continue; }
    2373       28538 :     if (signe(Pc) != s0) break;
    2374           0 :     shiftr_inplace(rb,1);
    2375           0 :   }
    2376       28538 :   ra = NULL;
    2377       28538 :   iter = 0;
    2378             :   for(;;)
    2379             :   {
    2380             :     GEN wdth;
    2381      326592 :     iter++;
    2382      326592 :     if (ra)
    2383      178389 :       rc = shiftr(addrr(ra, rb), -1);
    2384             :     else
    2385      148203 :       rc = shiftr(rb, -1);
    2386             :     do
    2387             :     {
    2388      326592 :       GEN pows = gen_powers(rc, rt, 1, NULL, _mp_sqr, _mp_mul, _gen_one);
    2389      326592 :       Pc = splitpoleval(Pp, Pm, pows, D, bitaddprec+2);
    2390      326592 :       if (Pc) break;
    2391           0 :       cprec++; rc = rtor(rc, cprec);
    2392           0 :     } while (1);
    2393      326592 :     if (signe(Pc) == s0)
    2394      113709 :       ra = rc;
    2395             :     else
    2396      212883 :       rb = rc;
    2397      326592 :     if (!ra) continue;
    2398      206927 :     wdth = subrr(rb, ra);
    2399      206927 :     if (!(iter % 8))
    2400             :     {
    2401       31828 :       GEN m1 = poleval(Pprime, ra), M2;
    2402       31828 :       if (signe(m1) == s0) continue;
    2403       31002 :       M2 = poleval(Pprime2, rb);
    2404       31002 :       if (abscmprr(gmul(M2, wdth), shiftr(m1, 1)) > 0) continue;
    2405       28538 :       break;
    2406             :     }
    2407      175099 :     else if (gexpo(wdth) <= -bitprec)
    2408           0 :       break;
    2409      298054 :   }
    2410       28538 :   rc = rb;
    2411       28538 :   av2 = avma;
    2412      209735 :   for(;; rc = gerepileuptoleaf(av2, rc))
    2413             :   {
    2414             :     long exponew;
    2415      238273 :     GEN Ppc, dist, rcold = rc;
    2416      238273 :     GEN pows = gen_powers(rc, rt, 1, NULL, _mp_sqr, _mp_mul, _gen_one);
    2417      238273 :     Ppc = splitpoleval(Pprimep, Pprimem, pows, D-1, bitaddprec+4);
    2418      238273 :     if (Ppc)
    2419      238273 :       Pc = splitpoleval(Pp, Pm, pows, D, bitaddprec+4);
    2420      238273 :     if (!Ppc || !Pc)
    2421             :     {
    2422       42403 :       if (cprec >= prec+addprec)
    2423        3407 :         cprec += EXTRAPRECWORD;
    2424             :       else
    2425       38996 :         cprec = minss(2*cprec, prec+addprec);
    2426       42403 :       rc = rtor(rc, cprec); continue; /* backtrack one step */
    2427             :     }
    2428      195870 :     dist = typ(Ppc) == t_REAL? divrr(Pc, Ppc): divri(Pc, Ppc);
    2429      195870 :     rc = subrr(rc, dist);
    2430      195870 :     if (cmprr(ra, rc) > 0 || cmprr(rb, rc) < 0)
    2431             :     {
    2432           0 :       if (cprec >= prec+addprec) break;
    2433           0 :       cprec = minss(2*cprec, prec+addprec);
    2434           0 :       rc = rtor(rcold, cprec); continue; /* backtrack one step */
    2435             :     }
    2436      195870 :     if (expoold == LONG_MAX) { expoold = expo(dist); continue; }
    2437      164439 :     exponew = expo(dist);
    2438      164439 :     if (exponew < -bitprec - 1)
    2439             :     {
    2440       37876 :       if (cprec >= prec+addprec) break;
    2441        9338 :       cprec = minss(2*cprec, prec+addprec);
    2442        9338 :       rc = rtor(rc, cprec); continue;
    2443             :     }
    2444      126563 :     if (exponew > expoold - 2)
    2445             :     {
    2446        2893 :       if (cprec >= prec+addprec) break;
    2447        2893 :       expoold = LONG_MAX;
    2448        2893 :       cprec = minss(2*cprec, prec+addprec);
    2449        2893 :       rc = rtor(rc, cprec); continue;
    2450             :     }
    2451      123670 :     expoold = exponew;
    2452      209735 :   }
    2453       28538 :   return gerepileuptoleaf(av, rtor(rc, prec));
    2454             : }
    2455             : 
    2456             : static GEN
    2457       25490 : usp(GEN Q0, long deg, long flag, long bitprec)
    2458             : {
    2459       25490 :   const pari_sp av = avma;
    2460             :   GEN Q, sol, c, Lc, Lk;
    2461       25490 :   long listsize = 64, nbr = 0, nb_todo, ind, deg0, indf, i, k, nb;
    2462             : 
    2463             : 
    2464       25490 :   sol = const_col(deg, gen_0);
    2465       25490 :   deg0 = deg;
    2466       25490 :   Lc = const_vec(listsize, gen_0);
    2467       25490 :   Lk = cgetg(listsize+1, t_VECSMALL);
    2468       25490 :   c = gen_0;
    2469       25490 :   k = Lk[1] = 0;
    2470       25490 :   ind = 1; indf = 2;
    2471       25490 :   Q = leafcopy(Q0);
    2472             : 
    2473       25490 :   nb_todo = 1;
    2474      488238 :   while (nb_todo)
    2475             :   {
    2476      437258 :     GEN nc = gel(Lc, ind), Qremapped;
    2477             :     pari_sp av2;
    2478             :     int root1;
    2479      437258 :     if (Lk[ind] == k + 1)
    2480             :     {
    2481      136584 :       deg0 = deg;
    2482      136584 :       setlg(Q0, deg + 3);
    2483      136584 :       Q0 = ZX_rescale2n(Q0, 1);
    2484      136584 :       Q = Q_primpart(Q0);
    2485      136584 :       c = gen_0;
    2486             :     }
    2487             : 
    2488      437258 :     if (!equalii(nc, c)) Q = ZX_translate(Q, subii(nc, c));
    2489             : 
    2490      437258 :     k = Lk[ind];
    2491      437258 :     c = nc;
    2492      437258 :     ind++;
    2493      437258 :     nb_todo--;
    2494             : 
    2495      437258 :     if (equalii(gel(Q, 2), gen_0))
    2496             :     { /* Q(0) = 0 */
    2497         168 :       GEN s = gmul2n(c, -k);
    2498             :       long j;
    2499         210 :       for (j = 1; j <= nbr; j++)
    2500         126 :         if (gequal(gel(sol, j), s)) break;
    2501         168 :       if (j > nbr) gel(sol, ++nbr) = s;
    2502             : 
    2503         168 :       deg0--;
    2504         168 :       for (j = 2; j <= deg0 + 2; j++) gel(Q, j) = gel(Q, j+1);
    2505         168 :       setlg(Q, j);
    2506             :     }
    2507             : 
    2508      437258 :     av2 = avma;
    2509      437258 :     nb = X2XP1(Q, deg0, &root1, flag == 1 ? &Qremapped : NULL);
    2510             : 
    2511      437258 :     if      (nb == 0) /* no root in this open interval */;
    2512      247325 :     else if (nb == 1) /* exactly one root */
    2513             :     {
    2514       41441 :       GEN s = gen_0;
    2515       41441 :       if (flag == 0)
    2516           0 :         s = mkvec2(gmul2n(c,-k), gmul2n(addiu(c,1),-k));
    2517       41441 :       else if (flag == 1) /* Caveat: Qremapped is the reciprocal polynomial */
    2518             :       {
    2519       28538 :         s = polsolve(Qremapped, bitprec+1);
    2520       28538 :         s = divrr(s, addsr(1, s));
    2521       28538 :         s = gmul2n(addir(c, s), -k);
    2522       28538 :         s = rtor(s, nbits2prec(bitprec));
    2523             :       }
    2524       41441 :       gel(sol, ++nbr) = gerepileupto(av2, s);
    2525             :     }
    2526             :     else
    2527             :     { /* unknown, add two nodes to refine */
    2528      205884 :       if (indf + 2 > listsize)
    2529             :       {
    2530         294 :         if (ind>1)
    2531             :         {
    2532        2702 :           for (i = ind; i < indf; i++)
    2533             :           {
    2534        2408 :             gel(Lc, i-ind+1) = gel(Lc, i);
    2535        2408 :             Lk[i-ind+1] = Lk[i];
    2536             :           }
    2537         294 :           indf -= ind-1;
    2538         294 :           ind = 1;
    2539             :         }
    2540         294 :         if (indf + 2 > listsize)
    2541             :         {
    2542           0 :           listsize *= 2;
    2543           0 :           Lc = vec_lengthen(Lc, listsize);
    2544           0 :           Lk = vecsmall_lengthen(Lk, listsize);
    2545             :         }
    2546         294 :         for (i = indf; i <= listsize; i++) gel(Lc, i) = gen_0;
    2547             :       }
    2548      205884 :       nc = shifti(c, 1);
    2549      205884 :       gel(Lc, indf) = nc;
    2550      205884 :       gel(Lc, indf + 1) = addiu(nc, 1);
    2551      205884 :       Lk[indf] = Lk[indf + 1] = k + 1;
    2552      205884 :       indf += 2;
    2553      205884 :       nb_todo += 2;
    2554             :     }
    2555      437258 :     if (root1)
    2556             :     { /* Q(1) = 0 */
    2557          70 :       GEN s = gmul2n(addiu(c,1), -k);
    2558             :       long j;
    2559         112 :       for (j = 1; j <= nbr; j++)
    2560          49 :         if (gequal(gel(sol, j), s)) break;
    2561          70 :       if (j > nbr) gel(sol, ++nbr) = s;
    2562             :     }
    2563             : 
    2564      437258 :     if (gc_needed(av, 2))
    2565             :     {
    2566           0 :       gerepileall(av, 6, &Q0, &Q, &c, &Lc, &Lk, &sol);
    2567           0 :       if (DEBUGMEM > 1) pari_warn(warnmem, "ZX_Uspensky", avma);
    2568             :     }
    2569             :   }
    2570             : 
    2571       25490 :   setlg(sol, nbr+1);
    2572       25490 :   return gerepilecopy(av, sol);
    2573             : }
    2574             : 
    2575             : static GEN
    2576        4249 : ZX_Uspensky_cst_pol(long nbz, long flag, long bitprec)
    2577             : {
    2578        4249 :   switch(flag)
    2579             :   {
    2580           0 :     case 0:  return zerocol(nbz);
    2581           0 :     case 1:  retconst_col(nbz, real_0_bit(-bitprec));
    2582        4249 :     default: return utoi(nbz);
    2583             :   }
    2584             : }
    2585             : 
    2586             : /* ZX_Uspensky(P, [a,a], flag) */
    2587             : static GEN
    2588          28 : ZX_Uspensky_equal(GEN P, GEN a, long flag)
    2589             : {
    2590          28 :   if (typ(a) != t_INFINITY && gequal0(poleval(P, a)))
    2591          21 :     return flag <= 1 ? mkcol(a): gen_1;
    2592             :   else
    2593           7 :     return flag <= 1 ? cgetg(1, t_COL) : gen_0;
    2594             : }
    2595             : 
    2596             : GEN
    2597       29282 : ZX_Uspensky(GEN P, GEN ab, long flag, long bitprec)
    2598             : {
    2599       29282 :   pari_sp av = avma;
    2600       29282 :   GEN a, b, sol = NULL, Pcur;
    2601             :   double fb;
    2602             :   long nbz, deg;
    2603             : 
    2604       29282 :   deg = degpol(P);
    2605       29282 :   if (deg == 0) return flag <= 1 ? cgetg(1, t_COL) : gen_0;
    2606       29282 :   if (ab)
    2607             :   {
    2608       16731 :     if (typ(ab) == t_VEC)
    2609             :     {
    2610        8756 :       if (lg(ab) != 3) pari_err_DIM("ZX_Uspensky");
    2611        8756 :       a = gel(ab, 1);
    2612        8756 :       b = gel(ab, 2);
    2613             :     }
    2614             :     else
    2615             :     {
    2616        7975 :       a = ab;
    2617        7975 :       b = mkoo();
    2618             :     }
    2619             :   }
    2620             :   else
    2621             :   {
    2622       12551 :     a = mkmoo();
    2623       12551 :     b = mkoo();
    2624             :   }
    2625       29282 :   switch (gcmp(a, b))
    2626             :   {
    2627           7 :     case 1: avma = av; return flag <= 1 ? cgetg(1, t_COL) : gen_0;
    2628          21 :     case 0: return gerepilecopy(av, ZX_Uspensky_equal(P, a, flag));
    2629             :   }
    2630       29254 :   nbz = ZX_valrem(P, &Pcur);
    2631       29254 :   deg -= nbz;
    2632       29254 :   if (!nbz) Pcur = P;
    2633       29254 :   if (nbz && (gsigne(a) > 0 || gsigne(b) < 0)) nbz = 0;
    2634       29254 :   if (deg == 0) { avma = av; return ZX_Uspensky_cst_pol(nbz, flag, bitprec); }
    2635       27959 :   if (deg == 1)
    2636             :   {
    2637        7275 :     sol = gdiv(gneg(gel(Pcur, 2)), pollead(Pcur, -1));
    2638        7275 :     if (gcmp(a, sol) > 0 || gcmp(sol, b) > 0)
    2639             :     {
    2640        2954 :       avma = av;
    2641        2954 :       return ZX_Uspensky_cst_pol(nbz, flag, bitprec);
    2642             :     }
    2643        4321 :     if (flag >= 2) { avma = av; return utoi(nbz+1); }
    2644        1983 :     sol = gconcat(zerocol(nbz), mkcol(sol));
    2645        1983 :     if (flag == 1) sol = RgC_gtofp(sol, nbits2prec(bitprec));
    2646        1983 :     return gerepilecopy(av, sol);
    2647             :   }
    2648       20684 :   switch(flag)
    2649             :   {
    2650             :     case 0:
    2651           0 :       sol = zerocol(nbz);
    2652           0 :       break;
    2653             :     case 1:
    2654       10108 :       sol = const_col(nbz, real_0_bit(-bitprec));
    2655       10108 :       break;
    2656             :     /* case 2: nothing */
    2657             :   }
    2658             : 
    2659       20684 :   if (typ(a) == t_INFINITY && typ(b) != t_INFINITY && gsigne(b))
    2660             :   {
    2661          35 :     fb = fujiwara_bound_real(Pcur, -1);
    2662          35 :     if (fb > -pariINFINITY) a = negi(int2n((long)ceil(fb))); else a = gen_0;
    2663             :   }
    2664       20684 :   if (typ(b) == t_INFINITY && typ(a) != t_INFINITY && gsigne(a))
    2665             :   {
    2666          21 :     fb = fujiwara_bound_real(Pcur, 1);
    2667          21 :     if (fb > -pariINFINITY) b = int2n((long)ceil(fb)); else b = gen_0;
    2668             :   }
    2669             : 
    2670       20684 :   if (typ(a) != t_INFINITY && typ(b) != t_INFINITY)
    2671             :   {
    2672             :     GEN den, diff, unscaledres, co, Pdiv, ascaled;
    2673             :     pari_sp av1;
    2674             :     long i;
    2675        7287 :     if (gequal(a,b)) /* can occur if one of a,b was initially a t_INFINITY */
    2676           7 :       return gerepilecopy(av, ZX_Uspensky_equal(P, a, flag));
    2677        7280 :     den = lcmii(Q_denom(a), Q_denom(b));
    2678        7280 :     if (!is_pm1(den))
    2679             :     {
    2680          42 :       Pcur = ZX_rescale(Pcur, den);
    2681          42 :       ascaled = gmul(a, den);
    2682             :     }
    2683             :     else
    2684             :     {
    2685        7238 :       den = NULL;
    2686        7238 :       ascaled = a;
    2687             :     }
    2688        7280 :     diff = subii(den ? gmul(b,den) : b, ascaled);
    2689        7280 :     Pcur = ZX_unscale(ZX_translate(Pcur, ascaled), diff);
    2690        7280 :     av1 = avma;
    2691        7280 :     Pdiv = cgetg(deg+2, t_POL);
    2692        7280 :     Pdiv[1] = Pcur[1];
    2693        7280 :     co = gel(Pcur, deg+2);
    2694       56959 :     for (i = deg; --i >= 0; )
    2695             :     {
    2696       42399 :       gel(Pdiv, i+2) = co;
    2697       42399 :       co = addii(co, gel(Pcur, i+2));
    2698             :     }
    2699        7280 :     if (!signe(co))
    2700             :     {
    2701          77 :       Pcur = Pdiv;
    2702          77 :       deg--;
    2703          77 :       if (flag <= 1)
    2704          21 :         sol = gconcat(sol, b);
    2705             :       else
    2706          56 :         nbz++;
    2707             :     }
    2708             :     else
    2709        7203 :       avma = av1;
    2710        7280 :     unscaledres = usp(Pcur, deg, flag, bitprec);
    2711        7280 :     if (flag <= 1)
    2712             :     {
    2713       19719 :       for (i = 1; i < lg(unscaledres); i++)
    2714             :       {
    2715       12565 :         GEN z = gmul(diff, gel(unscaledres, i));
    2716       12565 :         if (typ(z) == t_VEC)
    2717             :         {
    2718           0 :           gel(z, 1) = gadd(ascaled, gel(z, 1));
    2719           0 :           gel(z, 2) = gadd(ascaled, gel(z, 2));
    2720             :         }
    2721             :         else
    2722       12565 :           z = gadd(ascaled, z);
    2723       12565 :         if (den) z = gdiv(z, den);
    2724       12565 :         gel(unscaledres, i) = z;
    2725             :       }
    2726        7154 :       sol = gconcat(sol, unscaledres);
    2727             :     }
    2728             :     else
    2729         126 :       nbz += lg(unscaledres) - 1;
    2730             :   }
    2731       20677 :   if (typ(b) == t_INFINITY && (fb=fujiwara_bound_real(Pcur, 1)) > -pariINFINITY)
    2732             :   {
    2733             :     GEN Pcurp, unscaledres;
    2734       11973 :     long bp = (long)ceil(fb);
    2735       11973 :     if (bp < 0) bp = 0;
    2736       11973 :     Pcurp = ZX_unscale2n(Pcur, bp);
    2737       11973 :     unscaledres = usp(Pcurp, deg, flag, bitprec);
    2738       11973 :     if (flag <= 1)
    2739        2947 :       sol = shallowconcat(sol, gmul2n(unscaledres, bp));
    2740             :     else
    2741        9026 :       nbz += lg(unscaledres)-1;
    2742             :   }
    2743       20677 :   if (typ(a) == t_INFINITY && (fb=fujiwara_bound_real(Pcur,-1)) > -pariINFINITY)
    2744             :   {
    2745             :     GEN Pcurm, unscaledres;
    2746        6237 :     long i, bm = (long)ceil(fb);
    2747        6237 :     if (bm < 0) bm = 0;
    2748        6237 :     Pcurm = ZX_unscale2n(ZX_z_unscale(Pcur, -1), bm);
    2749        6237 :     unscaledres = usp(Pcurm, deg, flag, bitprec);
    2750        6237 :     if (flag <= 1)
    2751             :     {
    2752        6594 :       for (i = 1; i < lg(unscaledres); i++)
    2753             :       {
    2754        4242 :         GEN z = gneg(gmul2n(gel(unscaledres, i), bm));
    2755        4242 :         if (typ(z) == t_VEC) swap(gel(z, 1), gel(z, 2));
    2756        4242 :         gel(unscaledres, i) = z;
    2757             :       }
    2758        2352 :       sol = shallowconcat(unscaledres, sol);
    2759             :     }
    2760             :     else
    2761        3885 :       nbz += lg(unscaledres)-1;
    2762             :   }
    2763       20677 :   if (flag >= 2) return utoi(nbz);
    2764       10108 :   if (flag)
    2765       10108 :     sol = sort(sol);
    2766             :   else
    2767           0 :     sol = gen_sort(sol, (void *)_intervalcmp, cmp_nodata);
    2768       10108 :   return gerepileupto(av, sol);
    2769             : }
    2770             : 
    2771             : /* x a scalar */
    2772             : static GEN
    2773          35 : rootsdeg0(GEN x)
    2774             : {
    2775          35 :   if (!is_real_t(typ(x))) pari_err_TYPE("realroots",x);
    2776          28 :   if (gequal0(x)) pari_err_ROOTS0("realroots");
    2777          14 :   return cgetg(1,t_COL); /* constant polynomial */
    2778             : }
    2779             : static void
    2780        4674 : checkbound(GEN a)
    2781             : {
    2782        4674 :   switch(typ(a))
    2783             :   {
    2784        4674 :     case t_INT: case t_FRAC: case t_INFINITY: break;
    2785           0 :     default: pari_err_TYPE("polrealroots", a);
    2786             :   }
    2787        4674 : }
    2788             : static GEN
    2789        8450 : check_ab(GEN ab)
    2790             : {
    2791             :   GEN a, b;
    2792        8450 :   if (!ab) return NULL;
    2793        2337 :   if (typ(ab) != t_VEC || lg(ab) != 3) pari_err_TYPE("polrootsreal",ab);
    2794        2337 :   a = gel(ab,1); checkbound(a);
    2795        2337 :   b = gel(ab,2); checkbound(b);
    2796        3051 :   if (typ(a) == t_INFINITY && inf_get_sign(a) < 0 &&
    2797        1393 :       typ(b) == t_INFINITY && inf_get_sign(b) > 0) ab = NULL;
    2798        2337 :   return ab;
    2799             : }
    2800             : GEN
    2801        6169 : realroots(GEN P, GEN ab, long prec)
    2802             : {
    2803        6169 :   pari_sp av = avma;
    2804        6169 :   long nrr = 0;
    2805        6169 :   GEN sol = NULL, fa, ex;
    2806             :   long i, j, v;
    2807             : 
    2808        6169 :   ab = check_ab(ab);
    2809        6169 :   if (typ(P) != t_POL) return rootsdeg0(P);
    2810        6148 :   switch(degpol(P))
    2811             :   {
    2812           7 :     case -1: return rootsdeg0(gen_0);
    2813           7 :     case 0: return rootsdeg0(gel(P,2));
    2814             :   }
    2815        6134 :   if (!RgX_is_ZX(P)) P = RgX_rescale_to_int(P);
    2816        6134 :   P = Q_primpart(P);
    2817        6134 :   v = ZX_valrem(P,&P);
    2818        6134 :   if (v && (!ab || (gsigne(gel(ab,1)) <= 0 && gsigne(gel(ab,2)) >= 0)))
    2819        1204 :     sol = const_col(v, real_0(prec));
    2820        6134 :   fa = ZX_squff(P, &ex);
    2821       11085 :   for (i = 1; i < lg(fa); i++)
    2822             :   {
    2823        4951 :     GEN Pi = gel(fa, i), soli, soli2 = NULL;
    2824        4951 :     long n, nrri = 0, h;
    2825        4951 :     if (ab)
    2826          49 :       h = 1;
    2827             :     else
    2828        4902 :       Pi = ZX_deflate_max(Pi, &h);
    2829        4951 :     soli = ZX_Uspensky(Pi, h%2 ? ab: gen_0, 1, prec2nbits(prec));
    2830        4951 :     n = lg(soli);
    2831        4951 :     if (!(h % 2)) soli2 = cgetg(n, t_COL);
    2832       22095 :     for (j = 1; j < n; j++)
    2833             :     {
    2834       17144 :       GEN elt = gel(soli, j); /* != 0 */
    2835       17144 :       if (typ(elt) != t_REAL)
    2836             :       {
    2837          84 :         nrri++; if (h > 1 && !(h % 2)) nrri++;
    2838          84 :         elt = gtofp(elt, prec);
    2839          84 :         gel(soli, j) = elt;
    2840             :       }
    2841       17144 :       if (h > 1)
    2842             :       {
    2843             :         GEN r;
    2844        2732 :         if (h == 2)
    2845        2718 :           r = sqrtr(elt);
    2846             :         else
    2847             :         {
    2848          14 :           if (signe(elt) < 0)
    2849           7 :             r = negr(sqrtnr(negr(elt), h));
    2850             :           else
    2851           7 :             r = sqrtnr(elt, h);
    2852             :         }
    2853        2732 :         gel(soli, j) = r;
    2854        2732 :         if (!(h % 2)) gel(soli2, j) = negr(r);
    2855             :       }
    2856             :     }
    2857        4951 :     if (!(h % 2)) soli = shallowconcat(soli, soli2);
    2858        4951 :     if (ex[i] > 1) soli = shallowconcat1( const_vec(ex[i], soli) );
    2859        4951 :     sol = sol? shallowconcat(sol, soli): soli;
    2860        4951 :     nrr += ex[i]*nrri;
    2861             :   }
    2862        6134 :   if (!sol) { avma = av; return cgetg(1,t_COL); }
    2863             : 
    2864        6120 :   if (DEBUGLEVEL > 4)
    2865             :   {
    2866           0 :     err_printf("Number of real roots: %d\n", lg(sol)-1);
    2867           0 :     err_printf(" -- of which 2-integral: %ld\n", nrr);
    2868             :   }
    2869        6120 :   return gerepileupto(av, sort(sol));
    2870             : }
    2871             : 
    2872             : /* P non-constant, squarefree ZX */
    2873             : long
    2874       15589 : ZX_sturm(GEN P)
    2875             : {
    2876       15589 :   pari_sp av = avma;
    2877             :   long h, r;
    2878       15589 :   P = ZX_deflate_max(P, &h);
    2879       15589 :   if (odd(h))
    2880        9667 :     r = itos(ZX_Uspensky(P, NULL, 2, 0));
    2881             :   else
    2882        5922 :     r = 2*itos(ZX_Uspensky(P, gen_0, 2, 0));
    2883       15589 :   avma = av; return r;
    2884             : }
    2885             : /* P non-constant, squarefree ZX */
    2886             : long
    2887        2281 : ZX_sturmpart(GEN P, GEN ab)
    2888             : {
    2889        2281 :   pari_sp av = avma;
    2890             :   long r;
    2891        2281 :   if (!check_ab(ab)) return ZX_sturm(P);
    2892        1588 :   r = itos(ZX_Uspensky(P, ab, 2, 0));
    2893        1588 :   avma = av; return r;
    2894             : }

Generated by: LCOV version 1.11