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 to exceed 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 - kernel/none - mp.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.14.0 lcov report (development 27095-7e0f67a20c) Lines: 1127 1153 97.7 %
Date: 2021-11-29 07:04:58 Functions: 70 70 100.0 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : #line 2 "../src/kernel/none/mp.c"
       2             : /* Copyright (C) 2000-2003 The PARI group.
       3             : 
       4             : This file is part of the PARI/GP package.
       5             : 
       6             : PARI/GP is free software; you can redistribute it and/or modify it under the
       7             : terms of the GNU General Public License as published by the Free Software
       8             : Foundation; either version 2 of the License, or (at your option) any later
       9             : version. It is distributed in the hope that it will be useful, but WITHOUT
      10             : ANY WARRANTY WHATSOEVER.
      11             : 
      12             : Check the License for details. You should have received a copy of it, along
      13             : with the package; see the file 'COPYING'. If not, write to the Free Software
      14             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      15             : 
      16             : /***********************************************************************/
      17             : /**                                                                   **/
      18             : /**                       MULTIPRECISION KERNEL                       **/
      19             : /**                                                                   **/
      20             : /***********************************************************************/
      21             : #include "pari.h"
      22             : #include "paripriv.h"
      23             : #include "../src/kernel/none/tune-gen.h"
      24             : 
      25             : void
      26         722 : pari_kernel_init(void) { }
      27             : void
      28         720 : pari_kernel_close(void) { }
      29             : const char *
      30           2 : pari_kernel_version(void) { return ""; }
      31             : 
      32             : /* NOTE: arguments of "spec" routines (muliispec, addiispec, etc.) aren't
      33             :  * GENs but pairs (long *a, long na) representing a list of digits (in basis
      34             :  * BITS_IN_LONG) : a[0], ..., a[na-1]. [ In ordre to facilitate splitting: no
      35             :  * need to reintroduce codewords ] */
      36             : 
      37             : #define LIMBS(x)  ((x)+2)
      38             : #define NLIMBS(x) (lgefint(x)-2)
      39             : 
      40             : /* Normalize a nonnegative integer */
      41             : GEN
      42   402397899 : int_normalize(GEN x, long known_zero_words)
      43             : {
      44   402397899 :   long i, lx = lgefint(x);
      45             :   GEN x0;
      46   402397899 :   if (lx == 2) { x[1] = evalsigne(0) | evallgefint(2); return x; }
      47   402397899 :   if (!known_zero_words && x[2]) return x;
      48  1680323658 :   for (i = 2+known_zero_words; i < lx; i++)
      49  1641377967 :     if (x[i]) break;
      50   171097716 :   x0 = x; i -= 2; x += i;
      51   171097716 :   if (x0 == (GEN)avma) set_avma((pari_sp)x);
      52    93247959 :   else stackdummy((pari_sp)(x0+i), (pari_sp)x0);
      53   171097716 :   lx -= i;
      54   171097716 :   x[0] = evaltyp(t_INT) | evallg(lx);
      55   171097716 :   if (lx == 2) x[1] = evalsigne(0) | evallgefint(lx);
      56   132152025 :   else         x[1] = evalsigne(1) | evallgefint(lx);
      57   171097716 :   return x;
      58             : }
      59             : 
      60             : /***********************************************************************/
      61             : /**                                                                   **/
      62             : /**                      ADDITION / SUBTRACTION                       **/
      63             : /**                                                                   **/
      64             : /***********************************************************************/
      65             : 
      66             : GEN
      67     2237184 : setloop(GEN a)
      68             : {
      69     2237184 :   pari_sp av = avma;
      70     2237184 :   (void)cgetg(lgefint(a) + 3, t_VECSMALL);
      71     2237184 :   return icopy_avma(a, av); /* two cells of extra space before a */
      72             : }
      73             : 
      74             : /* we had a = setloop(?), then some incloops. Reset a to b */
      75             : GEN
      76      128046 : resetloop(GEN a, GEN b) {
      77      128046 :   long lb = lgefint(b);
      78      128046 :   a += lgefint(a) - lb;
      79      128046 :   a[0] = evaltyp(t_INT) | evallg(lb);
      80      128046 :   affii(b, a); return a;
      81             : }
      82             : 
      83             : /* assume a > 0, initialized by setloop. Do a++ */
      84             : static GEN
      85    20447436 : incpos(GEN a)
      86             : {
      87    20447436 :   long i, l = lgefint(a);
      88    20447439 :   for (i=l-1; i>1; i--)
      89    20447436 :     if (++a[i]) return a;
      90           3 :   l++; a--; /* use extra cell */
      91           3 :   a[0]=evaltyp(t_INT) | _evallg(l);
      92           3 :   a[1]=evalsigne(1) | evallgefint(l);
      93           3 :   a[2]=1; return a;
      94             : }
      95             : 
      96             : /* assume a < 0, initialized by setloop. Do a++ */
      97             : static GEN
      98       24969 : incneg(GEN a)
      99             : {
     100       24969 :   long i, l = lgefint(a)-1;
     101       24969 :   if (uel(a,l)--)
     102             :   {
     103       24966 :     if (l == 2 && !a[2])
     104             :     {
     105        1017 :       a++; /* save one cell */
     106        1017 :       a[0] = evaltyp(t_INT) | _evallg(2);
     107        1017 :       a[1] = evalsigne(0) | evallgefint(2);
     108             :     }
     109       24966 :     return a;
     110             :   }
     111           3 :   for (i = l-1;; i--) /* finishes since a[2] != 0 */
     112           3 :     if (uel(a,i)--) break;
     113           3 :   if (!a[2])
     114             :   {
     115           3 :     a++; /* save one cell */
     116           3 :     a[0] = evaltyp(t_INT) | _evallg(l);
     117           3 :     a[1] = evalsigne(-1) | evallgefint(l);
     118             :   }
     119           3 :   return a;
     120             : }
     121             : 
     122             : /* assume a initialized by setloop. Do a++ */
     123             : GEN
     124    20721141 : incloop(GEN a)
     125             : {
     126    20721141 :   switch(signe(a))
     127             :   {
     128      248736 :     case 0: a--; /* use extra cell */
     129      248736 :       a[0]=evaltyp(t_INT) | _evallg(3);
     130      248736 :       a[1]=evalsigne(1) | evallgefint(3);
     131      248736 :       a[2]=1; return a;
     132       24969 :     case -1: return incneg(a);
     133    20447436 :     default: return incpos(a);
     134             :   }
     135             : }
     136             : 
     137             : INLINE GEN
     138  1592662923 : adduispec(ulong s, GEN x, long nx)
     139             : {
     140  1592662923 :   GEN xd, zd = (GEN)avma;
     141             :   long lz;
     142             : 
     143  1592662923 :   if (nx == 1) return adduu(s, uel(x,0));
     144   177563388 :   lz = nx+3; (void)new_chunk(lz);
     145   177563388 :   xd = x + nx;
     146   177563388 :   *--zd = (ulong)*--xd + s;
     147   177563388 :   if ((ulong)*zd < s)
     148             :     for(;;)
     149             :     {
     150    22899414 :       if (xd == x) { *--zd = 1; break; } /* enlarge z */
     151    19186068 :       *--zd = ((ulong)*--xd) + 1;
     152    19186068 :       if (*zd) { lz--; break; }
     153             :     }
     154   161207982 :   else lz--;
     155   732786738 :   while (xd > x) *--zd = *--xd;
     156   177563388 :   *--zd = evalsigne(1) | evallgefint(lz);
     157   177563388 :   *--zd = evaltyp(t_INT) | evallg(lz);
     158   177563388 :   return gc_const((pari_sp)zd, zd);
     159             : }
     160             : 
     161             : GEN
     162   264995286 : adduispec_offset(ulong s, GEN x, long offset, long nx)
     163             : {
     164   264995286 :   GEN xd = x+lgefint(x)-nx-offset;
     165   336840027 :   while (nx && *xd==0) {xd++; nx--;}
     166   264995286 :   if (!nx) return utoi(s);
     167   236470446 :   return adduispec(s,xd,nx);
     168             : }
     169             : 
     170             : static GEN
     171  3358344078 : addiispec(GEN x, GEN y, long nx, long ny)
     172             : {
     173             :   GEN xd, yd, zd;
     174  3358344078 :   long lz, i = -2;
     175             :   LOCAL_OVERFLOW;
     176             : 
     177  3358344078 :   if (nx < ny) swapspec(x,y, nx,ny);
     178  3358344078 :   if (ny == 1) return adduispec(*y,x,nx);
     179  2061165903 :   zd = (GEN)avma;
     180  2061165903 :   lz = nx+3; (void)new_chunk(lz);
     181  2061165903 :   xd = x + nx;
     182  2061165903 :   yd = y + ny;
     183  2061165903 :   zd[-1] = addll(xd[-1], yd[-1]);
     184             : #ifdef addllx8
     185  1619577532 :   for (  ; i-8 > -ny; i-=8)
     186   932522231 :     addllx8(xd+i, yd+i, zd+i, overflow);
     187             : #endif
     188 25075427987 :   for (  ; i >= -ny; i--) zd[i] = addllx(xd[i], yd[i]);
     189  2061165903 :   if (overflow)
     190             :     for(;;)
     191             :     {
     192   269712342 :       if (i < -nx) { zd[i] = 1; i--; break; } /* enlarge z */
     193   181224450 :       zd[i] = uel(xd,i) + 1;
     194   181224450 :       if (zd[i]) { i--; lz--; break; }
     195    42447714 :       i--;
     196             :     }
     197  1833901275 :   else lz--;
     198 10196758938 :   for (; i >= -nx; i--) zd[i] = xd[i];
     199  2061165903 :   zd += i+1;
     200  2061165903 :   *--zd = evalsigne(1) | evallgefint(lz);
     201  2061165903 :   *--zd = evaltyp(t_INT) | evallg(lz);
     202  2061165903 :   return gc_const((pari_sp)zd, zd);
     203             : }
     204             : 
     205             : /* assume x >= s */
     206             : INLINE GEN
     207  1284190844 : subiuspec(GEN x, ulong s, long nx)
     208             : {
     209  1284190844 :   GEN xd, zd = (GEN)avma;
     210             :   long lz;
     211             :   LOCAL_OVERFLOW;
     212             : 
     213  1284190844 :   if (nx == 1) return utoi(x[0] - s);
     214             : 
     215   172677728 :   lz = nx+2; (void)new_chunk(lz);
     216   172677728 :   xd = x + nx;
     217   172677728 :   *--zd = subll(*--xd, s);
     218   172677728 :   if (overflow)
     219             :     for(;;)
     220             :     {
     221     7696251 :       *--zd = ((ulong)*--xd) - 1;
     222    50191149 :       if (*xd) break;
     223             :     }
     224   172677728 :   if (xd == x)
     225    68259294 :     while (*zd == 0) { zd++; lz--; } /* shorten z */
     226             :   else
     227  2343204524 :     do  *--zd = *--xd; while (xd > x);
     228   172677728 :   *--zd = evalsigne(1) | evallgefint(lz);
     229   172677728 :   *--zd = evaltyp(t_INT) | evallg(lz);
     230   172677728 :   return gc_const((pari_sp)zd, zd);
     231             : }
     232             : 
     233             : /* assume x > y */
     234             : static GEN
     235  2900328430 : subiispec(GEN x, GEN y, long nx, long ny)
     236             : {
     237             :   GEN xd,yd,zd;
     238  2900328430 :   long lz, i = -2;
     239             :   LOCAL_OVERFLOW;
     240             : 
     241  2900328430 :   if (ny==1) return subiuspec(x,*y,nx);
     242  1689736654 :   zd = (GEN)avma;
     243  1689736654 :   lz = nx+2; (void)new_chunk(lz);
     244  1689736654 :   xd = x + nx;
     245  1689736654 :   yd = y + ny;
     246  1689736654 :   zd[-1] = subll(xd[-1], yd[-1]);
     247             : #ifdef subllx8
     248  1657171145 :   for (  ; i-8 > -ny; i-=8)
     249  1093925594 :     subllx8(xd+i, yd+i, zd+i, overflow);
     250             : #endif
     251 25286441210 :   for (  ; i >= -ny; i--) zd[i] = subllx(xd[i], yd[i]);
     252  1689736654 :   if (overflow)
     253             :     for(;;)
     254             :     {
     255   319259097 :       zd[i] = uel(xd,i) - 1;
     256   490417405 :       if (xd[i--]) break;
     257             :     }
     258  1689736654 :   if (i>=-nx)
     259  2408066615 :     for (; i >= -nx; i--) zd[i] = xd[i];
     260             :   else
     261  1861268975 :     while (zd[i+1] == 0) { i++; lz--; } /* shorten z */
     262  1689736654 :   zd += i+1;
     263  1689736654 :   *--zd = evalsigne(1) | evallgefint(lz);
     264  1689736654 :   *--zd = evaltyp(t_INT) | evallg(lz);
     265  1689736654 :   return gc_const((pari_sp)zd, zd);
     266             : }
     267             : 
     268             : static void
     269   227128313 : roundr_up_ip(GEN x, long l)
     270             : {
     271   227128313 :   long i = l;
     272             :   for(;;)
     273             :   {
     274   227868215 :     if (++uel(x,--i)) break;
     275      991041 :     if (i == 2) { x[2] = (long)HIGHBIT; shiftr_inplace(x, 1); break; }
     276             :   }
     277   227128313 : }
     278             : 
     279             : void
     280   296811063 : affir(GEN x, GEN y)
     281             : {
     282   296811063 :   const long s = signe(x), ly = lg(y);
     283             :   long lx, sh, i;
     284             : 
     285   296811063 :   if (!s)
     286             :   {
     287    24898908 :     y[1] = evalexpo(-prec2nbits(ly));
     288    24898908 :     return;
     289             :   }
     290             : 
     291   271912155 :   lx = lgefint(x); sh = bfffo(x[2]);
     292   271912155 :   y[1] = evalsigne(s) | evalexpo(bit_accuracy(lx)-sh-1);
     293   271912155 :   if (sh) {
     294   268790778 :     if (lx <= ly)
     295             :     {
     296   539800185 :       for (i=lx; i<ly; i++) y[i]=0;
     297   168402414 :       shift_left(y,x,2,lx-1, 0,sh);
     298   168402414 :       return;
     299             :     }
     300   100388364 :     shift_left(y,x,2,ly-1, x[ly],sh);
     301             :     /* lx > ly: round properly */
     302   100388364 :     if ((uel(x,ly)<<sh) & HIGHBIT) roundr_up_ip(y, ly);
     303             :   }
     304             :   else {
     305     3121377 :     if (lx <= ly)
     306             :     {
     307     4697775 :       for (i=2; i<lx; i++) y[i]=x[i];
     308     2909874 :       for (   ; i<ly; i++) y[i]=0;
     309     1228788 :       return;
     310             :     }
     311     4360554 :     for (i=2; i<ly; i++) y[i]=x[i];
     312             :     /* lx > ly: round properly */
     313     1892589 :     if (uel(x,ly) & HIGHBIT) roundr_up_ip(y, ly);
     314             :   }
     315             : }
     316             : 
     317             : INLINE GEN
     318   851733769 : shiftispec(GEN x, long nx, long n)
     319             : {
     320             :   long ny, i, m;
     321             :   GEN y, yd;
     322   851733769 :   if (!n)  return icopyspec(x, nx);
     323             : 
     324   800328790 :   if (n > 0)
     325             :   {
     326   474923263 :     GEN z = (GEN)avma;
     327   474923263 :     long d = dvmdsBIL(n, &m);
     328             : 
     329   474923263 :     ny = nx+d; y = new_chunk(ny + 2); yd = y + 2;
     330  3570453598 :     for ( ; d; d--) *--z = 0;
     331   818367703 :     if (!m) for (i=0; i<nx; i++) yd[i]=x[i];
     332             :     else
     333             :     {
     334   467951014 :       const ulong sh = BITS_IN_LONG - m;
     335   467951014 :       shift_left(yd,x, 0,nx-1, 0,m);
     336   467951014 :       i = uel(x,0) >> sh;
     337             :       /* Extend y on the left? */
     338   467951014 :       if (i) { ny++; y = new_chunk(1); y[2] = i; }
     339             :     }
     340             :   }
     341             :   else
     342             :   {
     343   325405527 :     ny = nx - dvmdsBIL(-n, &m);
     344   325405527 :     if (ny<1) return gen_0;
     345   324236652 :     y = new_chunk(ny + 2); yd = y + 2;
     346   324236652 :     if (m) {
     347   240343773 :       shift_right(yd,x, 0,ny, 0,m);
     348   240343773 :       if (yd[0] == 0)
     349             :       {
     350    36803139 :         if (ny==1) return gc_const((pari_sp)(y+3), gen_0);
     351    29931324 :         ny--; set_avma((pari_sp)(++y));
     352             :       }
     353             :     } else {
     354  3616836621 :       for (i=0; i<ny; i++) yd[i]=x[i];
     355             :     }
     356             :   }
     357   792288100 :   y[1] = evalsigne(1)|evallgefint(ny + 2);
     358   792288100 :   y[0] = evaltyp(t_INT)|evallg(ny + 2); return y;
     359             : }
     360             : 
     361             : GEN
     362    38970999 : mantissa2nr(GEN x, long n)
     363             : { /*This is a kludge since x is not an integer*/
     364    38970999 :   long s = signe(x);
     365             :   GEN y;
     366             : 
     367    38970999 :   if(s == 0) return gen_0;
     368    38970084 :   y = shiftispec(x + 2, lg(x) - 2, n);
     369    38970084 :   if (signe(y)) setsigne(y, s);
     370    38970084 :   return y;
     371             : }
     372             : 
     373             : GEN
     374     1460826 : truncr(GEN x)
     375             : {
     376             :   long d,e,i,s,m;
     377             :   GEN y;
     378             : 
     379     1460826 :   if ((s=signe(x)) == 0 || (e=expo(x)) < 0) return gen_0;
     380      804876 :   d = nbits2lg(e+1); m = remsBIL(e);
     381      804876 :   if (d > lg(x)) pari_err_PREC( "truncr (precision loss in truncation)");
     382             : 
     383      804873 :   y=cgeti(d); y[1] = evalsigne(s) | evallgefint(d);
     384      804873 :   if (++m == BITS_IN_LONG)
     385        1155 :     for (i=2; i<d; i++) y[i]=x[i];
     386             :   else
     387      804345 :     shift_right(y,x, 2,d,0, BITS_IN_LONG - m);
     388      804873 :   return y;
     389             : }
     390             : 
     391             : /* integral part */
     392             : GEN
     393     1545102 : floorr(GEN x)
     394             : {
     395             :   long d,e,i,lx,m;
     396             :   GEN y;
     397             : 
     398     1545102 :   if (signe(x) >= 0) return truncr(x);
     399      610851 :   if ((e=expo(x)) < 0) return gen_m1;
     400      228063 :   d = nbits2lg(e+1); m = remsBIL(e);
     401      228063 :   lx=lg(x); if (d>lx) pari_err_PREC( "floorr (precision loss in truncation)");
     402      228063 :   y = new_chunk(d);
     403      228063 :   if (++m == BITS_IN_LONG)
     404             :   {
     405        3318 :     for (i=2; i<d; i++) y[i]=x[i];
     406        8070 :     i=d; while (i<lx && !x[i]) i++;
     407        1425 :     if (i==lx) goto END;
     408             :   }
     409             :   else
     410             :   {
     411      226638 :     shift_right(y,x, 2,d,0, BITS_IN_LONG - m);
     412      226638 :     if (uel(x,d-1)<<m == 0)
     413             :     {
     414      404988 :       i=d; while (i<lx && !x[i]) i++;
     415       78153 :       if (i==lx) goto END;
     416             :     }
     417             :   }
     418             :   /* set y:=y+1 */
     419      181728 :   for (i=d-1; i>=2; i--) { uel(y,i)++; if (y[i]) goto END; }
     420           0 :   y=new_chunk(1); y[2]=1; d++;
     421      228063 : END:
     422      228063 :   y[1] = evalsigne(-1) | evallgefint(d);
     423      228063 :   y[0] = evaltyp(t_INT) | evallg(d); return y;
     424             : }
     425             : 
     426             : INLINE int
     427  3598683537 : cmpiispec(GEN x, GEN y, long lx, long ly)
     428             : {
     429             :   long i;
     430  3598683537 :   if (lx < ly) return -1;
     431  3424784022 :   if (lx > ly) return  1;
     432  3468764093 :   i = 0; while (i<lx && x[i]==y[i]) i++;
     433  3126795572 :   if (i==lx) return 0;
     434  2958518273 :   return (uel(x,i) > uel(y,i))? 1: -1;
     435             : }
     436             : 
     437             : INLINE int
     438   182086833 : equaliispec(GEN x, GEN y, long lx, long ly)
     439             : {
     440             :   long i;
     441   182086833 :   if (lx != ly) return 0;
     442   333481851 :   i = ly-1; while (i>=0 && x[i]==y[i]) i--;
     443   182009322 :   return i < 0;
     444             : }
     445             : 
     446             : /***********************************************************************/
     447             : /**                                                                   **/
     448             : /**                          MULTIPLICATION                           **/
     449             : /**                                                                   **/
     450             : /***********************************************************************/
     451             : /* assume ny > 0 */
     452             : INLINE GEN
     453  3334730934 : muluispec(ulong x, GEN y, long ny)
     454             : {
     455  3334730934 :   GEN yd, z = (GEN)avma;
     456  3334730934 :   long lz = ny+3;
     457             :   LOCAL_HIREMAINDER;
     458             : 
     459  3334730934 :   (void)new_chunk(lz);
     460  3334730934 :   yd = y + ny; *--z = mulll(x, *--yd);
     461 11256159285 :   while (yd > y) *--z = addmul(x,*--yd);
     462  3334730934 :   if (hiremainder) *--z = hiremainder; else lz--;
     463  3334730934 :   *--z = evalsigne(1) | evallgefint(lz);
     464  3334730934 :   *--z = evaltyp(t_INT) | evallg(lz);
     465  3334730934 :   return gc_const((pari_sp)z, z);
     466             : }
     467             : 
     468             : /* a + b*|Y| */
     469             : GEN
     470      419412 : addumului(ulong a, ulong b, GEN Y)
     471             : {
     472             :   GEN yd,y,z;
     473             :   long ny,lz;
     474             :   LOCAL_HIREMAINDER;
     475             :   LOCAL_OVERFLOW;
     476             : 
     477      419412 :   if (!b || !signe(Y)) return utoi(a);
     478             : 
     479      419364 :   y = LIMBS(Y); z = (GEN)avma;
     480      419364 :   ny = NLIMBS(Y);
     481      419364 :   lz = ny+3;
     482             : 
     483      419364 :   (void)new_chunk(lz);
     484      419364 :   yd = y + ny; *--z = addll(a, mulll(b, *--yd));
     485      419364 :   if (overflow) hiremainder++; /* can't overflow */
     486      419364 :   while (yd > y) *--z = addmul(b,*--yd);
     487      419364 :   if (hiremainder) *--z = hiremainder; else lz--;
     488      419364 :   *--z = evalsigne(1) | evallgefint(lz);
     489      419364 :   *--z = evaltyp(t_INT) | evallg(lz);
     490      419364 :   return gc_const((pari_sp)z, z);
     491             : }
     492             : 
     493             : /***********************************************************************/
     494             : /**                                                                   **/
     495             : /**                          DIVISION                                 **/
     496             : /**                                                                   **/
     497             : /***********************************************************************/
     498             : 
     499             : ulong
     500  1287692895 : umodiu(GEN y, ulong x)
     501             : {
     502  1287692895 :   long sy=signe(y),ly,i;
     503             :   ulong xi;
     504             :   LOCAL_HIREMAINDER;
     505             : 
     506  1287692895 :   if (!x) pari_err_INV("umodiu",gen_0);
     507  1287692895 :   if (!sy) return 0;
     508  1040919345 :   ly = lgefint(y);
     509  1040919345 :   if (x <= uel(y,2))
     510             :   {
     511   327548424 :     hiremainder=0;
     512   327548424 :     if (ly==3)
     513             :     {
     514   300645561 :       hiremainder=uel(y,2)%x;
     515   300645561 :       if (!hiremainder) return 0;
     516   252995640 :       return (sy > 0)? hiremainder: x - hiremainder;
     517             :     }
     518             :   }
     519             :   else
     520             :   {
     521   713370921 :     if (ly==3) return (sy > 0)? uel(y,2): x - uel(y,2);
     522    89824542 :     hiremainder=uel(y,2); ly--; y++;
     523             :   }
     524   116727405 :   xi = get_Fl_red(x);
     525   497124360 :   for (i=2; i<ly; i++) (void)divll_pre(y[i],x,xi);
     526   116727405 :   if (!hiremainder) return 0;
     527   110684841 :   return (sy > 0)? hiremainder: x - hiremainder;
     528             : }
     529             : 
     530             : /* return |y| \/ x */
     531             : GEN
     532   278549646 : absdiviu_rem(GEN y, ulong x, ulong *rem)
     533             : {
     534             :   long ly,i;
     535             :   GEN z;
     536             :   ulong xi;
     537             :   LOCAL_HIREMAINDER;
     538             : 
     539   278549646 :   if (!x) pari_err_INV("absdiviu_rem",gen_0);
     540   278549646 :   if (!signe(y)) { *rem = 0; return gen_0; }
     541             : 
     542   269944122 :   ly = lgefint(y);
     543   269944122 :   if (x <= uel(y,2))
     544             :   {
     545   210825354 :     hiremainder=0;
     546   210825354 :     if (ly==3)
     547             :     {
     548   100881822 :       z = cgetipos(3);
     549   100881822 :       z[2] = divll(uel(y,2),x);
     550   100881822 :       *rem = hiremainder; return z;
     551             :     }
     552             :   }
     553             :   else
     554             :   {
     555    59118768 :     if (ly==3) { *rem = uel(y,2); return gen_0; }
     556    44176581 :     hiremainder = uel(y,2); ly--; y++;
     557             :   }
     558   154120113 :   xi = get_Fl_red(x);
     559   154120113 :   z = cgetipos(ly);
     560  1135390446 :   for (i=2; i<ly; i++) z[i]=divll_pre(y[i],x,xi);
     561   154120113 :   *rem = hiremainder; return z;
     562             : }
     563             : 
     564             : GEN
     565    66001749 : divis_rem(GEN y, long x, long *rem)
     566             : {
     567    66001749 :   long sy=signe(y),ly,s,i;
     568             :   GEN z;
     569             :   ulong xi;
     570             :   LOCAL_HIREMAINDER;
     571             : 
     572    66001749 :   if (!x) pari_err_INV("divis_rem",gen_0);
     573    66001749 :   if (!sy) { *rem=0; return gen_0; }
     574    48535083 :   if (x<0) { s = -sy; x = -x; } else s = sy;
     575             : 
     576    48535083 :   ly = lgefint(y);
     577    48535083 :   if ((ulong)x <= uel(y,2))
     578             :   {
     579    34232610 :     hiremainder=0;
     580    34232610 :     if (ly==3)
     581             :     {
     582    33707742 :       z = cgeti(3); z[1] = evallgefint(3) | evalsigne(s);
     583    33707742 :       z[2] = divll(uel(y,2),x);
     584    33707742 :       if (sy<0) hiremainder = - ((long)hiremainder);
     585    33707742 :       *rem = (long)hiremainder; return z;
     586             :     }
     587             :   }
     588             :   else
     589             :   {
     590    14302473 :     if (ly==3) { *rem = itos(y); return gen_0; }
     591      291657 :     hiremainder = uel(y,2); ly--; y++;
     592             :   }
     593      816525 :   xi = get_Fl_red(x);
     594      816525 :   z = cgeti(ly); z[1] = evallgefint(ly) | evalsigne(s);
     595     3161028 :   for (i=2; i<ly; i++) z[i]=divll_pre(y[i],x,xi);
     596      816525 :   if (sy<0) hiremainder = - ((long)hiremainder);
     597      816525 :   *rem = (long)hiremainder; return z;
     598             : }
     599             : 
     600             : GEN
     601     2261853 : divis(GEN y, long x)
     602             : {
     603     2261853 :   long sy=signe(y),ly,s,i;
     604             :   ulong xi;
     605             :   GEN z;
     606             :   LOCAL_HIREMAINDER;
     607             : 
     608     2261853 :   if (!x) pari_err_INV("divis",gen_0);
     609     2261853 :   if (!sy) return gen_0;
     610     2261847 :   if (x<0) { s = -sy; x = -x; } else s = sy;
     611             : 
     612     2261847 :   ly = lgefint(y);
     613     2261847 :   if ((ulong)x <= uel(y,2))
     614             :   {
     615     2242440 :     hiremainder=0;
     616     2242440 :     if (ly==3)
     617             :     {
     618     2064156 :       z = cgeti(3); z[1] = evallgefint(3) | evalsigne(s);
     619     2064156 :       z[2] = divll(y[2],x);
     620     2064156 :       return z;
     621             :     }
     622             :   }
     623             :   else
     624             :   {
     625       19407 :     if (ly==3) return gen_0;
     626       19407 :     hiremainder=y[2]; ly--; y++;
     627             :   }
     628      197691 :   xi = get_Fl_red(x);
     629      197691 :   z = cgeti(ly); z[1] = evallgefint(ly) | evalsigne(s);
     630     3343212 :   for (i=2; i<ly; i++) z[i]=divll_pre(y[i],x, xi);
     631      197691 :   return z;
     632             : }
     633             : 
     634             : GEN
     635    74062125 : divrr(GEN x, GEN y)
     636             : {
     637    74062125 :   long sx=signe(x), sy=signe(y), lx,ly,lr,e,i,j;
     638             :   ulong y0,y1;
     639             :   GEN r, r1;
     640             : 
     641    74062125 :   if (!x) pari_err_INV("divrr",y);
     642    74062125 :   e = expo(x) - expo(y);
     643    74062125 :   if (!sx) return real_0_bit(e);
     644    73838748 :   if (sy<0) sx = -sx;
     645             : 
     646    73838748 :   lx=lg(x); ly=lg(y);
     647    73838748 :   if (ly==3)
     648             :   {
     649    16724457 :     ulong k = x[2], l = (lx>3)? x[3]: 0;
     650             :     LOCAL_HIREMAINDER;
     651    16724457 :     if (k < uel(y,2)) e--;
     652             :     else
     653             :     {
     654     5914650 :       l >>= 1; if (k&1) l |= HIGHBIT;
     655     5914650 :       k >>= 1;
     656             :     }
     657    16724457 :     hiremainder = k; k = divll(l,y[2]);
     658    16724457 :     if (hiremainder > (uel(y,2) >> 1) && !++k) { k = HIGHBIT; e++; }
     659    16724457 :     r = cgetr(3);
     660    16724457 :     r[1] = evalsigne(sx) | evalexpo(e);
     661    16724457 :     r[2] = k; return r;
     662             :   }
     663             : 
     664    57114291 :   lr = minss(lx,ly); r = new_chunk(lr);
     665    57114291 :   r1 = r-1;
     666   364509306 :   r1[1] = 0; for (i=2; i<lr; i++) r1[i]=x[i];
     667    57114291 :   r1[lr] = (lx>ly)? x[lr]: 0;
     668    57114291 :   y0 = y[2]; y1 = y[3];
     669   421623597 :   for (i=0; i<lr-1; i++)
     670             :   { /* r1 = r + (i-1), OK up to r1[2] (accesses at most r[lr]) */
     671             :     ulong k, qp;
     672             :     LOCAL_HIREMAINDER;
     673             :     LOCAL_OVERFLOW;
     674             : 
     675   364509306 :     if (uel(r1,1) == y0) { qp = ULONG_MAX; k = addll(y0,r1[2]); }
     676             :     else
     677             :     {
     678   360524319 :       if (uel(r1,1) > y0) /* can't happen if i=0 */
     679             :       {
     680           0 :         GEN y1 = y+1;
     681           0 :         j = lr-i; r1[j] = subll(r1[j],y1[j]);
     682           0 :         for (j--; j>0; j--) r1[j] = subllx(r1[j],y1[j]);
     683           0 :         j=i; do uel(r,--j)++; while (j && !uel(r,j));
     684             :       }
     685   360524319 :       hiremainder = r1[1]; overflow = 0;
     686   360524319 :       qp = divll(r1[2],y0); k = hiremainder;
     687             :     }
     688   364509306 :     j = lr-i+1;
     689   364509306 :     if (!overflow)
     690             :     {
     691             :       long k3, k4;
     692   361463424 :       k3 = mulll(qp,y1);
     693   361463424 :       if (j == 3) /* i = lr - 2 maximal, r1[3] undefined -> 0 */
     694    57054813 :         k4 = subll(hiremainder,k);
     695             :       else
     696             :       {
     697   304408611 :         k3 = subll(k3, r1[3]);
     698   304408611 :         k4 = subllx(hiremainder,k);
     699             :       }
     700   462982858 :       while (!overflow && k4) { qp--; k3 = subll(k3,y1); k4 = subllx(k4,y0); }
     701             :     }
     702   364509306 :     if (j<ly) (void)mulll(qp,y[j]); else { hiremainder = 0 ; j = ly; }
     703  2878183569 :     for (j--; j>1; j--)
     704             :     {
     705  2513674263 :       r1[j] = subll(r1[j], addmul(qp,y[j]));
     706  2513674263 :       hiremainder += overflow;
     707             :     }
     708   364509306 :     if (uel(r1,1) != hiremainder)
     709             :     {
     710      658332 :       if (uel(r1,1) < hiremainder)
     711             :       {
     712      658332 :         qp--;
     713      658332 :         j = lr-i-(lr-i>=ly); r1[j] = addll(r1[j], y[j]);
     714     5787633 :         for (j--; j>1; j--) r1[j] = addllx(r1[j], y[j]);
     715             :       }
     716             :       else
     717             :       {
     718           0 :         r1[1] -= hiremainder;
     719           0 :         while (r1[1])
     720             :         {
     721           0 :           qp++; if (!qp) { j=i; do uel(r,--j)++; while (j && !r[j]); }
     722           0 :           j = lr-i-(lr-i>=ly); r1[j] = subll(r1[j],y[j]);
     723           0 :           for (j--; j>1; j--) r1[j] = subllx(r1[j],y[j]);
     724           0 :           r1[1] -= overflow;
     725             :         }
     726             :       }
     727             :     }
     728   364509306 :     *++r1 = qp;
     729             :   }
     730             :   /* i = lr-1 */
     731             :   /* round correctly */
     732    57114291 :   if (uel(r1,1) > (y0>>1))
     733             :   {
     734    27005891 :     j=i; do uel(r,--j)++; while (j && !r[j]);
     735             :   }
     736   364509306 :   r1 = r-1; for (j=i; j>=2; j--) r[j]=r1[j];
     737    57114291 :   if (r[0] == 0) e--;
     738    28130379 :   else if (r[0] == 1) { shift_right(r,r, 2,lr, 1,1); }
     739             :   else { /* possible only when rounding up to 0x2 0x0 ... */
     740           0 :     r[2] = (long)HIGHBIT; e++;
     741             :   }
     742    57114291 :   r[0] = evaltyp(t_REAL)|evallg(lr);
     743    57114291 :   r[1] = evalsigne(sx) | evalexpo(e);
     744    57114291 :   return r;
     745             : }
     746             : 
     747             : GEN
     748    35651256 : divri(GEN x, GEN y)
     749             : {
     750    35651256 :   long lx, s = signe(y);
     751             :   pari_sp av;
     752             :   GEN z;
     753             : 
     754    35651256 :   if (!s) pari_err_INV("divri",y);
     755    35651256 :   if (!signe(x)) return real_0_bit(expo(x) - expi(y));
     756    35491254 :   if (!is_bigint(y)) {
     757    23872353 :     GEN z = divru(x, y[2]);
     758    23872353 :     if (s < 0) togglesign(z);
     759    23872353 :     return z;
     760             :   }
     761    11618901 :   lx = lg(x); z = cgetr(lx); av = avma;
     762    11618901 :   affrr(divrr(x, itor(y, lx+1)), z);
     763    11618901 :   return gc_const(av, z);
     764             : }
     765             : 
     766             : /* Integer division x / y: such that sign(r) = sign(x)
     767             :  *   if z = ONLY_REM return remainder, otherwise return quotient
     768             :  *   if z != NULL set *z to remainder
     769             :  *   *z is the last object on stack (and thus can be disposed of with cgiv
     770             :  *   instead of gerepile)
     771             :  * If *z is zero, we put gen_0 here and no copy.
     772             :  * space needed: lx + ly */
     773             : GEN
     774  1052226300 : dvmdii(GEN x, GEN y, GEN *z)
     775             : {
     776  1052226300 :   long sx = signe(x), sy = signe(y);
     777  1052226300 :   long lx, ly = lgefint(y), lz, i, j, sh, lq, lr;
     778             :   pari_sp av;
     779             :   ulong y0,y0i,y1, *xd,*rd,*qd;
     780             :   GEN q, r, r1;
     781             : 
     782  1052226300 :   if (!sx)
     783             :   {
     784    41253855 :     if (ly < 3) pari_err_INV("dvmdii",gen_0);
     785    41253852 :     if (!z || z == ONLY_REM) return gen_0;
     786    14181060 :     *z=gen_0; return gen_0;
     787             :   }
     788  1010972445 :   if (ly <= 3)
     789             :   {
     790             :     ulong rem;
     791   538795434 :     if (ly < 3) pari_err_INV("dvmdii",gen_0);
     792   538795434 :     if (z == ONLY_REM)
     793             :     {
     794   438380802 :       rem = umodiu(x,uel(y,2));
     795   438380802 :       if (!rem) return gen_0;
     796   395276490 :       return (sx < 0)? utoineg(uel(y,2) - rem): utoipos(rem);
     797             :     }
     798   100414632 :     q = absdiviu_rem(x, uel(y,2), &rem);
     799   100414632 :     if (sx != sy) togglesign(q);
     800   100414632 :     if (!z) return q;
     801    97537731 :     if (!rem) *z = gen_0;
     802    39008910 :     else *z = sx < 0? utoineg(rem): utoipos(rem);
     803    97537731 :     return q;
     804             :   }
     805   472177011 :   lx=lgefint(x);
     806   472177011 :   lz=lx-ly;
     807   472177011 :   if (lz <= 0)
     808             :   {
     809   234755295 :     if (lz == 0)
     810             :     {
     811   211323756 :       for (i=2; i<lx; i++)
     812   210577257 :         if (x[i] != y[i])
     813             :         {
     814   196456653 :           if (uel(x,i) > uel(y,i)) goto DIVIDE;
     815    52619064 :           goto TRIVIAL;
     816             :         }
     817      746499 :       if (z == ONLY_REM) return gen_0;
     818       45180 :       if (z) *z = gen_0;
     819       45180 :       if (sx < 0) sy = -sy;
     820       45180 :       return stoi(sy);
     821             :     }
     822    37552143 : TRIVIAL:
     823    90171207 :     if (z == ONLY_REM) return icopy(x);
     824     1496994 :     if (z) *z = icopy(x);
     825     1496994 :     return gen_0;
     826             :   }
     827   237421716 : DIVIDE: /* quotient is nonzero */
     828   381259305 :   av=avma; if (sx<0) sy = -sy;
     829   381259305 :   r1 = new_chunk(lx); sh = bfffo(y[2]);
     830   381259305 :   if (sh)
     831             :   { /* normalize so that highbit(y) = 1 (shift left x and y by sh bits)*/
     832   372442830 :     const ulong m = BITS_IN_LONG - sh;
     833   372442830 :     r = new_chunk(ly);
     834   372442830 :     shift_left(r, y,2,ly-1, 0,sh); y = r;
     835   372442830 :     shift_left(r1,x,2,lx-1, 0,sh);
     836   372442830 :     r1[1] = uel(x,2) >> m;
     837             :   }
     838             :   else
     839             :   {
     840    83192583 :     r1[1] = 0; for (j=2; j<lx; j++) r1[j] = x[j];
     841             :   }
     842   381259305 :   x = r1;
     843   381259305 :   y0 = y[2]; y0i = get_Fl_red(y0);
     844   381259305 :   y1 = y[3];
     845  1775308242 :   for (i=0; i<=lz; i++)
     846             :   { /* r1 = x + i */
     847             :     ulong k, qp;
     848             :     LOCAL_HIREMAINDER;
     849             :     LOCAL_OVERFLOW;
     850             : 
     851  1394048937 :     if (uel(r1,1) == y0)
     852             :     {
     853       22371 :       qp = ULONG_MAX; k = addll(y0,r1[2]);
     854             :     }
     855             :     else
     856             :     {
     857  1394026566 :       hiremainder = r1[1]; overflow = 0;
     858  1394026566 :       qp = divll_pre(r1[2],y0,y0i); k = hiremainder;
     859             :     }
     860  1394048937 :     if (!overflow)
     861             :     {
     862  1394027661 :       long k3 = subll(mulll(qp,y1), r1[3]);
     863  1394027661 :       long k4 = subllx(hiremainder,k);
     864  1624643378 :       while (!overflow && k4) { qp--; k3 = subll(k3,y1); k4 = subllx(k4,y0); }
     865             :     }
     866  1394048937 :     hiremainder = 0; j = ly;
     867 35390497311 :     for (j--; j>1; j--)
     868             :     {
     869 33996448374 :       r1[j] = subll(r1[j], addmul(qp,y[j]));
     870 33996448374 :       hiremainder += overflow;
     871             :     }
     872  1394048937 :     if (uel(r1,1) < hiremainder)
     873             :     {
     874     5390754 :       qp--;
     875     5390754 :       j = ly-1; r1[j] = addll(r1[j],y[j]);
     876    24061851 :       for (j--; j>1; j--) r1[j] = addllx(r1[j],y[j]);
     877             :     }
     878  1394048937 :     *++r1 = qp;
     879             :   }
     880             : 
     881   381259305 :   lq = lz+2;
     882   381259305 :   if (!z)
     883             :   {
     884     2277432 :     qd = (ulong*)av;
     885     2277432 :     xd = (ulong*)(x + lq);
     886     2277432 :     if (x[1]) { lz++; lq++; }
     887    31174122 :     while (lz--) *--qd = *--xd;
     888     2277432 :     *--qd = evalsigne(sy) | evallgefint(lq);
     889     2277432 :     *--qd = evaltyp(t_INT) | evallg(lq);
     890     2277432 :     return gc_const((pari_sp)qd, (GEN)qd);
     891             :   }
     892             : 
     893   458926416 :   j=lq; while (j<lx && !x[j]) j++;
     894   378981873 :   lz = lx-j;
     895   378981873 :   if (z == ONLY_REM)
     896             :   {
     897   231102477 :     if (lz==0) return gc_const(av, gen_0);
     898   222215994 :     rd = (ulong*)av; lr = lz+2;
     899   222215994 :     xd = (ulong*)(x + lx);
     900   254821110 :     if (!sh) while (lz--) *--rd = *--xd;
     901             :     else
     902             :     { /* shift remainder right by sh bits */
     903   214340283 :       const ulong shl = BITS_IN_LONG - sh;
     904             :       ulong l;
     905   214340283 :       xd--;
     906   900149895 :       while (--lz) /* fill r[3..] */
     907             :       {
     908   685809612 :         l = *xd >> sh;
     909   685809612 :         *--rd = l | (*--xd << shl);
     910             :       }
     911   214340283 :       l = *xd >> sh;
     912   214340283 :       if (l) *--rd = l; else lr--;
     913             :     }
     914   222215994 :     *--rd = evalsigne(sx) | evallgefint(lr);
     915   222215994 :     *--rd = evaltyp(t_INT) | evallg(lr);
     916   222215994 :     return gc_const((pari_sp)rd, (GEN)rd);
     917             :   }
     918             : 
     919   147879396 :   lr = lz+2;
     920   147879396 :   rd = NULL; /* gcc -Wall */
     921   147879396 :   if (lz)
     922             :   { /* non zero remainder: initialize rd */
     923   145082907 :     xd = (ulong*)(x + lx);
     924   145082907 :     if (!sh)
     925             :     {
     926      761871 :       rd = (ulong*)avma; (void)new_chunk(lr);
     927     3136728 :       while (lz--) *--rd = *--xd;
     928             :     }
     929             :     else
     930             :     { /* shift remainder right by sh bits */
     931   144321036 :       const ulong shl = BITS_IN_LONG - sh;
     932             :       ulong l;
     933   144321036 :       rd = (ulong*)x; /* overwrite shifted y */
     934   144321036 :       xd--;
     935   503197791 :       while (--lz)
     936             :       {
     937   358876755 :         l = *xd >> sh;
     938   358876755 :         *--rd = l | (*--xd << shl);
     939             :       }
     940   144321036 :       l = *xd >> sh;
     941   144321036 :       if (l) *--rd = l; else lr--;
     942             :     }
     943   145082907 :     *--rd = evalsigne(sx) | evallgefint(lr);
     944   145082907 :     *--rd = evaltyp(t_INT) | evallg(lr);
     945   145082907 :     rd += lr;
     946             :   }
     947   147879396 :   qd = (ulong*)av;
     948   147879396 :   xd = (ulong*)(x + lq);
     949   147879396 :   if (x[1]) lq++;
     950   415650822 :   j = lq-2; while (j--) *--qd = *--xd;
     951   147879396 :   *--qd = evalsigne(sy) | evallgefint(lq);
     952   147879396 :   *--qd = evaltyp(t_INT) | evallg(lq);
     953   147879396 :   q = (GEN)qd;
     954   147879396 :   if (lr==2) *z = gen_0;
     955             :   else
     956             :   { /* rd has been properly initialized: we had lz > 0 */
     957   865658197 :     while (lr--) *--qd = *--rd;
     958   145082907 :     *z = (GEN)qd;
     959             :   }
     960   147879396 :   return gc_const((pari_sp)qd, q);
     961             : }
     962             : 
     963             : /* Montgomery reduction.
     964             :  * N has k words, assume T >= 0 has less than 2k.
     965             :  * Return res := T / B^k mod N, where B = 2^BIL
     966             :  * such that 0 <= res < T/B^k + N  and  res has less than k words */
     967             : GEN
     968    25047678 : red_montgomery(GEN T, GEN N, ulong inv)
     969             : {
     970             :   pari_sp av;
     971             :   GEN Te, Td, Ne, Nd, scratch;
     972    25047678 :   ulong i, j, m, t, d, k = NLIMBS(N);
     973             :   int carry;
     974             :   LOCAL_HIREMAINDER;
     975             :   LOCAL_OVERFLOW;
     976             : 
     977    25047678 :   if (k == 0) return gen_0;
     978    25047678 :   d = NLIMBS(T); /* <= 2*k */
     979    25047678 :   if (d == 0) return gen_0;
     980             : #ifdef DEBUG
     981             :   if (d > 2*k) pari_err_BUG("red_montgomery");
     982             : #endif
     983    25047669 :   if (k == 1)
     984             :   { /* as below, special cased for efficiency */
     985       94269 :     ulong n = uel(N,2);
     986       94269 :     if (d == 1) {
     987       94173 :       hiremainder = uel(T,2);
     988       94173 :       m = hiremainder * inv;
     989       94173 :       (void)addmul(m, n); /* t + m*n = 0 */
     990       94173 :       return utoi(hiremainder);
     991             :     } else { /* d = 2 */
     992          96 :       hiremainder = uel(T,3);
     993          96 :       m = hiremainder * inv;
     994          96 :       (void)addmul(m, n); /* t + m*n = 0 */
     995          96 :       t = addll(hiremainder, uel(T,2));
     996          96 :       if (overflow) t -= n; /* t > n doesn't fit in 1 word */
     997          96 :       return utoi(t);
     998             :     }
     999             :   }
    1000             :   /* assume k >= 2 */
    1001    24953400 :   av = avma; scratch = new_chunk(k<<1); /* >= k + 2: result fits */
    1002             : 
    1003             :   /* copy T to scratch space (pad with zeroes to 2k words) */
    1004    24953400 :   Td = (GEN)av;
    1005    24953400 :   Te = T + (d+2);
    1006   665458377 :   for (i=0; i < d     ; i++) *--Td = *--Te;
    1007    44902131 :   for (   ; i < (k<<1); i++) *--Td = 0;
    1008             : 
    1009    24953400 :   Te = (GEN)av; /* 1 beyond end of current T mantissa (in scratch) */
    1010    24953400 :   Ne = N + k+2; /* 1 beyond end of N mantissa */
    1011             : 
    1012    24953400 :   carry = 0;
    1013   355180254 :   for (i=0; i<k; i++) /* set T := T/B nod N, k times */
    1014             :   {
    1015   330226854 :     Td = Te; /* one beyond end of (new) T mantissa */
    1016   330226854 :     Nd = Ne;
    1017   330226854 :     hiremainder = *--Td;
    1018   330226854 :     m = hiremainder * inv; /* solve T + m N = O(B) */
    1019             : 
    1020             :     /* set T := (T + mN) / B */
    1021   330226854 :     Te = Td;
    1022   330226854 :     (void)addmul(m, *--Nd); /* = 0 */
    1023  6024637956 :     for (j=1; j<k; j++)
    1024             :     {
    1025  5694411102 :       t = addll(addmul(m, *--Nd), *--Td);
    1026  5694411102 :       *Td = t;
    1027  5694411102 :       hiremainder += overflow;
    1028             :     }
    1029   330226854 :     t = addll(hiremainder, *--Td); *Td = t + carry;
    1030   330226854 :     carry = (overflow || (carry && *Td == 0));
    1031             :   }
    1032    24953400 :   if (carry)
    1033             :   { /* Td > N overflows (k+1 words), set Td := Td - N */
    1034      353751 :     Td = Te;
    1035      353751 :     Nd = Ne;
    1036      353751 :     t = subll(*--Td, *--Nd); *Td = t;
    1037     6742518 :     while (Td > scratch) { t = subllx(*--Td, *--Nd); *Td = t; }
    1038             :   }
    1039             : 
    1040             :   /* copy result */
    1041    24953400 :   Td = (GEN)av;
    1042    28856910 :   while (*scratch == 0 && Te > scratch) scratch++; /* strip leading 0s */
    1043   351276744 :   while (Te > scratch) *--Td = *--Te;
    1044    24953400 :   k = (GEN)av - Td; if (!k) return gc_const(av, gen_0);
    1045    24953400 :   k += 2;
    1046    24953400 :   *--Td = evalsigne(1) | evallgefint(k);
    1047    24953400 :   *--Td = evaltyp(t_INT) | evallg(k);
    1048             : #ifdef DEBUG
    1049             : {
    1050             :   long l = NLIMBS(N), s = BITS_IN_LONG*l;
    1051             :   GEN R = int2n(s);
    1052             :   GEN res = remii(mulii(T, Fp_inv(R, N)), N);
    1053             :   if (k > lgefint(N)
    1054             :     || !equalii(remii(Td,N),res)
    1055             :     || cmpii(Td, addii(shifti(T, -s), N)) >= 0) pari_err_BUG("red_montgomery");
    1056             : }
    1057             : #endif
    1058    24953400 :   return gc_const((pari_sp)Td, Td);
    1059             : }
    1060             : 
    1061             : /* EXACT INTEGER DIVISION */
    1062             : 
    1063             : /* assume xy>0, the division is exact and y is odd. Destroy x */
    1064             : static GEN
    1065    30928593 : diviuexact_i(GEN x, ulong y)
    1066             : {
    1067             :   long i, lz, lx;
    1068             :   ulong q, yinv;
    1069             :   GEN z, z0, x0, x0min;
    1070             : 
    1071    30928593 :   if (y == 1) return icopy(x);
    1072    25542876 :   lx = lgefint(x);
    1073    25542876 :   if (lx == 3)
    1074             :   {
    1075      850542 :     q = uel(x,2) / y;
    1076      850542 :     if (!q) pari_err_OP("exact division", x, utoi(y));
    1077      850542 :     return utoipos(q);
    1078             :   }
    1079    24692334 :   yinv = invmod2BIL(y);
    1080    24692334 :   lz = (y <= uel(x,2)) ? lx : lx-1;
    1081    24692334 :   z = new_chunk(lz);
    1082    24692334 :   z0 = z + lz;
    1083    24692334 :   x0 = x + lx; x0min = x + lx-lz+2;
    1084             : 
    1085    79302948 :   while (x0 > x0min)
    1086             :   {
    1087    54610614 :     *--z0 = q = yinv*uel(--x0,0); /* i-th quotient */
    1088    54610614 :     if (!q) continue;
    1089             :     /* x := x - q * y */
    1090             :     { /* update neither lowest word (could set it to 0) nor highest ones */
    1091    54093609 :       GEN x1 = x0 - 1;
    1092             :       LOCAL_HIREMAINDER;
    1093    54093609 :       (void)mulll(q,y);
    1094    54093609 :       if (hiremainder)
    1095             :       {
    1096    44254983 :         if (uel(x1,0) < hiremainder)
    1097             :         {
    1098       81339 :           uel(x1,0) -= hiremainder;
    1099       83295 :           do uel(--x1,0)--; while (uel(x1,0) == ULONG_MAX);
    1100             :         }
    1101             :         else
    1102    44173644 :           uel(x1,0) -= hiremainder;
    1103             :       }
    1104             :     }
    1105             :   }
    1106    24692334 :   i=2; while(!z[i]) i++;
    1107    24692334 :   z += i-2; lz -= i-2;
    1108    24692334 :   z[0] = evaltyp(t_INT)|evallg(lz);
    1109    24692334 :   z[1] = evalsigne(1)|evallg(lz);
    1110    24692334 :   if (lz == 2) pari_err_OP("exact division", x, utoi(y));
    1111    24692334 :   return gc_const((pari_sp)z, z);
    1112             : }
    1113             : 
    1114             : /* assume y != 0 and the division is exact */
    1115             : GEN
    1116    18666213 : diviuexact(GEN x, ulong y)
    1117             : {
    1118             :   pari_sp av;
    1119    18666213 :   long lx, vy, s = signe(x);
    1120             :   GEN z;
    1121             : 
    1122    18666213 :   if (!s) return gen_0;
    1123    17829000 :   if (y == 1) return icopy(x);
    1124    14852229 :   lx = lgefint(x);
    1125    14852229 :   if (lx == 3) {
    1126    11576184 :     ulong q = uel(x,2) / y;
    1127    11576184 :     if (!q) pari_err_OP("exact division", x, utoi(y));
    1128    11576184 :     return (s > 0)? utoipos(q): utoineg(q);
    1129             :   }
    1130     3276045 :   av = avma; (void)new_chunk(lx); vy = vals(y);
    1131     3276045 :   if (vy) {
    1132     1429542 :     y >>= vy;
    1133     1429542 :     if (y == 1) { set_avma(av); return shifti(x, -vy); }
    1134      741351 :     x = shifti(x, -vy);
    1135      741351 :     if (lx == 3) {
    1136           0 :       ulong q = uel(x,2) / y;
    1137           0 :       set_avma(av);
    1138           0 :       if (!q) pari_err_OP("exact division", x, utoi(y));
    1139           0 :       return (s > 0)? utoipos(q): utoineg(q);
    1140             :     }
    1141     1846503 :   } else x = icopy(x);
    1142     2587854 :   set_avma(av);
    1143     2587854 :   z = diviuexact_i(x, y);
    1144     2587854 :   setsigne(z, s); return z;
    1145             : }
    1146             : 
    1147             : /* Find z such that x=y*z, knowing that y | x (unchecked)
    1148             :  * Method: y0 z0 = x0 mod B = 2^BITS_IN_LONG ==> z0 = 1/y0 mod B.
    1149             :  *    Set x := (x - z0 y) / B, updating only relevant words, and repeat */
    1150             : GEN
    1151   299604600 : diviiexact(GEN x, GEN y)
    1152             : {
    1153   299604600 :   long lx, ly, lz, vy, i, ii, sx = signe(x), sy = signe(y);
    1154             :   pari_sp av;
    1155             :   ulong y0inv,q;
    1156             :   GEN z;
    1157             : 
    1158   299604600 :   if (!sy) pari_err_INV("diviiexact",gen_0);
    1159   299604600 :   if (!sx) return gen_0;
    1160   238008387 :   lx = lgefint(x);
    1161   238008387 :   if (lx == 3) {
    1162   190071684 :     q = uel(x,2) / uel(y,2);
    1163   190071684 :     if (!q) pari_err_OP("exact division", x, y);
    1164   190071684 :     return (sx+sy) ? utoipos(q): utoineg(q);
    1165             :   }
    1166    47936703 :   vy = vali(y); av = avma;
    1167    47936703 :   (void)new_chunk(lx); /* enough room for z */
    1168    47936703 :   if (vy)
    1169             :   { /* make y odd */
    1170    22314102 :     y = shifti(y,-vy);
    1171    22314102 :     x = shifti(x,-vy); lx = lgefint(x);
    1172             :   }
    1173    25622601 :   else x = icopy(x); /* necessary because we destroy x */
    1174    47936703 :   set_avma(av); /* will erase our x,y when exiting */
    1175             :   /* now y is odd */
    1176    47936703 :   ly = lgefint(y);
    1177    47936703 :   if (ly == 3)
    1178             :   {
    1179    28340739 :     z = diviuexact_i(x,uel(y,2)); /* x != 0 */
    1180    28340739 :     setsigne(z, (sx+sy)? 1: -1); return z;
    1181             :   }
    1182    19595964 :   y0inv = invmod2BIL(y[ly-1]);
    1183    34077882 :   i=2; while (i<ly && y[i]==x[i]) i++;
    1184    19595964 :   lz = (i==ly || uel(y,i) < uel(x,i)) ? lx-ly+3 : lx-ly+2;
    1185    19595964 :   z = new_chunk(lz);
    1186             : 
    1187    19595964 :   y += ly - 1; /* now y[-i] = i-th word of y */
    1188    85749870 :   for (ii=lx-1,i=lz-1; i>=2; i--,ii--)
    1189             :   {
    1190             :     long limj;
    1191             :     LOCAL_HIREMAINDER;
    1192             :     LOCAL_OVERFLOW;
    1193             : 
    1194    66153906 :     z[i] = q = y0inv*uel(x,ii); /* i-th quotient */
    1195    66153906 :     if (!q) continue;
    1196             : 
    1197             :     /* x := x - q * y */
    1198    66080898 :     (void)mulll(q,y[0]); limj = maxss(lx - lz, ii+3-ly);
    1199             :     { /* update neither lowest word (could set it to 0) nor highest ones */
    1200    66080898 :       GEN x0 = x + (ii - 1), y0 = y - 1, xlim = x + limj;
    1201   762398934 :       for (; x0 >= xlim; x0--, y0--)
    1202             :       {
    1203   696318036 :         *x0 = subll(*x0, addmul(q,*y0));
    1204   696318036 :         hiremainder += overflow;
    1205             :       }
    1206    66080898 :       if (hiremainder && limj != lx - lz)
    1207             :       {
    1208    34680147 :         if ((ulong)*x0 < hiremainder)
    1209             :         {
    1210      402024 :           *x0 -= hiremainder;
    1211      402024 :           do (*--x0)--; while ((ulong)*x0 == ULONG_MAX);
    1212             :         }
    1213             :         else
    1214    34278123 :           *x0 -= hiremainder;
    1215             :       }
    1216             :     }
    1217             :   }
    1218    19595964 :   i=2; while(!z[i]) i++;
    1219    19595964 :   z += i-2; lz -= (i-2);
    1220    19595964 :   z[0] = evaltyp(t_INT)|evallg(lz);
    1221    19595964 :   z[1] = evalsigne((sx+sy)? 1: -1) | evallg(lz);
    1222    19595964 :   if (lz == 2) pari_err_OP("exact division", x, y);
    1223    19595964 :   return gc_const((pari_sp)z, z);
    1224             : }
    1225             : 
    1226             : /* assume yz != and yz | x */
    1227             : GEN
    1228      141813 : diviuuexact(GEN x, ulong y, ulong z)
    1229             : {
    1230             :   long tmp[4];
    1231             :   ulong t;
    1232             :   LOCAL_HIREMAINDER;
    1233      141813 :   t = mulll(y, z);
    1234      141813 :   if (!hiremainder) return diviuexact(x, t);
    1235           0 :   tmp[0] = evaltyp(t_INT)|_evallg(4);
    1236           0 :   tmp[1] = evalsigne(1)|evallgefint(4);
    1237           0 :   tmp[2] = hiremainder;
    1238           0 :   tmp[3] = t;
    1239           0 :   return diviiexact(x, tmp);
    1240             : }
    1241             : 
    1242             : /********************************************************************/
    1243             : /**                                                                **/
    1244             : /**               INTEGER MULTIPLICATION (BASECASE)                **/
    1245             : /**                                                                **/
    1246             : /********************************************************************/
    1247             : /* nx >= ny = num. of digits of x, y (not GEN, see mulii) */
    1248             : INLINE GEN
    1249  3566275872 : muliispec_basecase(GEN x, GEN y, long nx, long ny)
    1250             : {
    1251             :   GEN z2e,z2d,yd,xd,ye,zd;
    1252             :   long p1,lz;
    1253             :   LOCAL_HIREMAINDER;
    1254             : 
    1255  3566275872 :   if (ny == 1) return muluispec((ulong)*y, x, nx);
    1256   864863124 :   if (ny == 0) return gen_0;
    1257   864030174 :   zd = (GEN)avma; lz = nx+ny+2;
    1258   864030174 :   (void)new_chunk(lz);
    1259   864030174 :   xd = x + nx;
    1260   864030174 :   yd = y + ny;
    1261   864030174 :   ye = yd; p1 = *--xd;
    1262             : 
    1263   864030174 :   *--zd = mulll(p1, *--yd); z2e = zd;
    1264  5531819241 :   while (yd > y) *--zd = addmul(p1, *--yd);
    1265   864030174 :   *--zd = hiremainder;
    1266             : 
    1267  7082716248 :   while (xd > x)
    1268             :   {
    1269             :     LOCAL_OVERFLOW;
    1270  6218686074 :     yd = ye; p1 = *--xd;
    1271             : 
    1272  6218686074 :     z2d = --z2e;
    1273  6218686074 :     *z2d = addll(mulll(p1, *--yd), *z2d); z2d--;
    1274 68701792446 :     while (yd > y)
    1275             :     {
    1276 62483106372 :       hiremainder += overflow;
    1277 62483106372 :       *z2d = addll(addmul(p1, *--yd), *z2d); z2d--;
    1278             :     }
    1279  6218686074 :     *--zd = hiremainder + overflow;
    1280             :   }
    1281   864030174 :   if (*zd == 0) { zd++; lz--; } /* normalize */
    1282   864030174 :   *--zd = evalsigne(1) | evallgefint(lz);
    1283   864030174 :   *--zd = evaltyp(t_INT) | evallg(lz);
    1284   864030174 :   return gc_const((pari_sp)zd, zd);
    1285             : }
    1286             : 
    1287             : INLINE GEN
    1288   709155135 : sqrispec_basecase(GEN x, long nx)
    1289             : {
    1290             :   GEN z2e,z2d,yd,xd,zd,x0,z0;
    1291             :   long p1,lz;
    1292             :   LOCAL_HIREMAINDER;
    1293             :   LOCAL_OVERFLOW;
    1294             : 
    1295   709155135 :   if (nx == 1) return sqru((ulong)*x);
    1296   492939840 :   if (nx == 0) return gen_0;
    1297   114369291 :   zd = (GEN)avma; lz = (nx+1) << 1;
    1298   114369291 :   z0 = new_chunk(lz);
    1299   114369291 :   if (nx == 1)
    1300             :   {
    1301           0 :     *--zd = mulll(*x, *x);
    1302           0 :     *--zd = hiremainder; goto END;
    1303             :   }
    1304   114369291 :   xd = x + nx;
    1305             : 
    1306             :   /* compute double products --> zd */
    1307   114369291 :   p1 = *--xd; yd = xd; --zd;
    1308   114369291 :   *--zd = mulll(p1, *--yd); z2e = zd;
    1309   825283608 :   while (yd > x) *--zd = addmul(p1, *--yd);
    1310   114369291 :   *--zd = hiremainder;
    1311             : 
    1312   114369291 :   x0 = x+1;
    1313   825283608 :   while (xd > x0)
    1314             :   {
    1315             :     LOCAL_OVERFLOW;
    1316   710914317 :     p1 = *--xd; yd = xd;
    1317             : 
    1318   710914317 :     z2e -= 2; z2d = z2e;
    1319   710914317 :     *z2d = addll(mulll(p1, *--yd), *z2d); z2d--;
    1320  6558528300 :     while (yd > x)
    1321             :     {
    1322  5847613983 :       hiremainder += overflow;
    1323  5847613983 :       *z2d = addll(addmul(p1, *--yd), *z2d); z2d--;
    1324             :     }
    1325   710914317 :     *--zd = hiremainder + overflow;
    1326             :   }
    1327             :   /* multiply zd by 2 (put result in zd - 1) */
    1328   114369291 :   zd[-1] = ((*zd & HIGHBIT) != 0);
    1329   114369291 :   shift_left(zd, zd, 0, (nx<<1)-3, 0, 1);
    1330             : 
    1331             :   /* add the squares */
    1332   114369291 :   xd = x + nx; zd = z0 + lz;
    1333   114369291 :   p1 = *--xd;
    1334   114369291 :   zd--; *zd = mulll(p1,p1);
    1335   114369291 :   zd--; *zd = addll(hiremainder, *zd);
    1336   939652899 :   while (xd > x)
    1337             :   {
    1338   825283608 :     p1 = *--xd;
    1339   825283608 :     zd--; *zd = addll(mulll(p1,p1)+ overflow, *zd);
    1340   825283608 :     zd--; *zd = addll(hiremainder + overflow, *zd);
    1341             :   }
    1342             : 
    1343   114369291 : END:
    1344   114369291 :   if (*zd == 0) { zd++; lz--; } /* normalize */
    1345   114369291 :   *--zd = evalsigne(1) | evallgefint(lz);
    1346   114369291 :   *--zd = evaltyp(t_INT) | evallg(lz);
    1347   114369291 :   return gc_const((pari_sp)zd, zd);
    1348             : }
    1349             : 
    1350             : /********************************************************************/
    1351             : /**                                                                **/
    1352             : /**               INTEGER MULTIPLICATION (FFT)                     **/
    1353             : /**                                                                **/
    1354             : /********************************************************************/
    1355             : 
    1356             : /*
    1357             :  Compute parameters for FFT:
    1358             :    len: result length
    1359             :    k: FFT depth.
    1360             :    n: number of blocks (2^k)
    1361             :    bs: block size
    1362             :    mod: Modulus is M=2^(BIL*mod)+1
    1363             :    ord: order of 2 in Z/MZ.
    1364             :  We must have:
    1365             :    bs*n >= l
    1366             :    2^(BIL*mod) > nb*2^(2*BIL*bs)
    1367             :    2^k | 2*BIL*mod
    1368             : */
    1369             : static void
    1370       47313 : mulliifft_params(long len, long *k, long *mod, long *bs, long *n, ulong *ord)
    1371             : {
    1372             :   long r;
    1373       47313 :   *k = expu((3*len)>>2)-3;
    1374             :   do {
    1375       47313 :     (*k)--;
    1376       47313 :     r  = *k-(TWOPOTBITS_IN_LONG+2);
    1377       47313 :     *n = 1L<<*k;
    1378       47313 :     *bs = (len+*n-1)>>*k;
    1379       47313 :     *mod= 2**bs+1;
    1380       47313 :     if (r>0)
    1381        1407 :       *mod=((*mod+(1L<<r)-1)>>r)<<r;
    1382       47313 :   } while(*mod>=3**bs);
    1383       47313 :   *ord= 4**mod*BITS_IN_LONG;
    1384       47313 : }
    1385             : 
    1386             : /* Zf_: arithmetic in ring Z/MZ where M= 2^(BITS_IN_LONG*mod)+1
    1387             :  * for some mod.
    1388             :  * Do not garbage collect.
    1389             :  */
    1390             : 
    1391             : static GEN
    1392    88338048 : Zf_add(GEN a, GEN b, GEN M)
    1393             : {
    1394    88338048 :   GEN y, z = addii(a,b);
    1395    88338048 :   long mod = lgefint(M)-3;
    1396    88338048 :   long l = NLIMBS(z);
    1397    88338048 :   if (l<=mod) return z;
    1398    33879606 :   y = subiu(z, 1);
    1399    33879606 :   if (NLIMBS(y)<=mod) return z;
    1400    33879606 :   return int_normalize(y,1);
    1401             : }
    1402             : 
    1403             : static GEN
    1404    89664429 : Zf_sub(GEN a, GEN b, GEN M)
    1405             : {
    1406    89664429 :   GEN z = subii(a,b);
    1407    89664429 :   return signe(z)>=0? z: addii(M,z);
    1408             : }
    1409             : 
    1410             : /* destroy z */
    1411             : static GEN
    1412   188017812 : Zf_red_destroy(GEN z, GEN M)
    1413             : {
    1414   188017812 :   long mod = lgefint(M)-3;
    1415   188017812 :   long l = NLIMBS(z);
    1416             :   GEN y;
    1417   188017812 :   if (l<=mod) return z;
    1418    83011740 :   y = shifti(z, -mod*BITS_IN_LONG);
    1419    83011740 :   z = int_normalize(z, NLIMBS(y));
    1420    83011740 :   y = Zf_red_destroy(y, M);
    1421    83011740 :   z = subii(z, y);
    1422    83011740 :   if (signe(z)<0) z = addii(z, M);
    1423    83011740 :   return z;
    1424             : }
    1425             : 
    1426             : INLINE GEN
    1427    97181688 : Zf_shift(GEN a, ulong s, GEN M) { return Zf_red_destroy(shifti(a, s), M); }
    1428             : 
    1429             : /*
    1430             :  Multiply by sqrt(2)^s
    1431             :  We use the formula sqrt(2)=z_8*(1-z_4)) && z_8=2^(ord/16) [2^(ord/4)+1]
    1432             : */
    1433             : 
    1434             : static GEN
    1435    88338048 : Zf_mulsqrt2(GEN a, ulong s, ulong ord, GEN M)
    1436             : {
    1437    88338048 :   ulong hord = ord>>1;
    1438    88338048 :   if (!signe(a)) return gen_0;
    1439    86704542 :   if (odd(s)) /* Multiply by 2^(s/2) */
    1440             :   {
    1441     1326381 :     GEN az8  = Zf_shift(a,   ord>>4, M);
    1442     1326381 :     GEN az83 = Zf_shift(az8, ord>>3, M);
    1443     1326381 :     a = Zf_sub(az8, az83, M);
    1444     1326381 :     s--;
    1445             :   }
    1446    86704542 :   if (s < hord)
    1447    64729071 :     return Zf_shift(a, s>>1, M);
    1448             :   else
    1449    21975471 :     return subii(M,Zf_shift(a, (s-hord)>>1, M));
    1450             : }
    1451             : 
    1452             : INLINE GEN
    1453      354048 : Zf_sqr(GEN a, GEN M) { return Zf_red_destroy(sqri(a), M); }
    1454             : 
    1455             : INLINE GEN
    1456     7470336 : Zf_mul(GEN a, GEN b, GEN M) { return Zf_red_destroy(mulii(a,b), M); }
    1457             : 
    1458             : /* In place, bit reversing FFT */
    1459             : static void
    1460    15201579 : muliifft_dit(ulong o, ulong ord, GEN M, GEN FFT, long d, long step)
    1461             : {
    1462    15201579 :   pari_sp av = avma;
    1463             :   long i;
    1464    15201579 :   ulong j, no = (o<<1)%ord;
    1465    15201579 :   long hstep=step>>1;
    1466    73599147 :   for (i = d+1, j = 0; i <= d+hstep; ++i, j =(j+o)%ord)
    1467             :   {
    1468    58397568 :     GEN a = Zf_add(gel(FFT,i), gel(FFT,i+hstep), M);
    1469    58397568 :     GEN b = Zf_mulsqrt2(Zf_sub(gel(FFT,i), gel(FFT,i+hstep), M), j, ord, M);
    1470    58397568 :     affii(a,gel(FFT,i));
    1471    58397568 :     affii(b,gel(FFT,i+hstep));
    1472    58397568 :     set_avma(av);
    1473             :   }
    1474    15201579 :   if (hstep>1)
    1475             :   {
    1476     7554219 :     muliifft_dit(no, ord, M, FFT, d, hstep);
    1477     7554219 :     muliifft_dit(no, ord, M, FFT, d+hstep, hstep);
    1478             :   }
    1479    15201579 : }
    1480             : 
    1481             : /* In place, bit reversed FFT, inverse of muliifft_dit */
    1482             : static void
    1483     7777071 : muliifft_dis(ulong o, ulong ord, GEN M, GEN FFT, long d, long step)
    1484             : {
    1485     7777071 :   pari_sp av = avma;
    1486             :   long i;
    1487     7777071 :   ulong j, no = (o<<1)%ord;
    1488     7777071 :   long hstep=step>>1;
    1489     7777071 :   if (hstep>1)
    1490             :   {
    1491     3864879 :     muliifft_dis(no, ord, M, FFT, d, hstep);
    1492     3864879 :     muliifft_dis(no, ord, M, FFT, d+hstep, hstep);
    1493             :   }
    1494    37717551 :   for (i = d+1, j = 0; i <= d+hstep; ++i, j =(j+o)%ord)
    1495             :   {
    1496    29940480 :     GEN z = Zf_mulsqrt2(gel(FFT,i+hstep), j, ord, M);
    1497    29940480 :     GEN a = Zf_add(gel(FFT,i), z, M);
    1498    29940480 :     GEN b = Zf_sub(gel(FFT,i), z, M);
    1499    29940480 :     affii(a,gel(FFT,i));
    1500    29940480 :     affii(b,gel(FFT,i+hstep));
    1501    29940480 :     set_avma(av);
    1502             :   }
    1503     7777071 : }
    1504             : 
    1505             : static GEN
    1506       93141 : muliifft_spliti(GEN a, long na, long bs, long n, long mod)
    1507             : {
    1508       93141 :   GEN ap = a+na-1;
    1509       93141 :   GEN c  = cgetg(n+1, t_VEC);
    1510             :   long i,j;
    1511    15387861 :   for(i=1;i<=n;i++)
    1512             :   {
    1513    15294720 :     GEN z = cgeti(mod+3);
    1514    15294720 :     if (na)
    1515             :     {
    1516     7560588 :       long m = minss(bs, na), v=0;
    1517     7560588 :       GEN zp, aa=ap-m+1;
    1518    31340796 :       while (!*aa && v<m) {aa++; v++;}
    1519     7560588 :       zp = z+m-v+1;
    1520   203127159 :       for (j=v; j < m; j++)
    1521   195566571 :         *zp-- = *ap--;
    1522     7560588 :       ap -= v; na -= m;
    1523     7560588 :       z[1] = evalsigne(m!=v) | evallgefint(m-v+2);
    1524             :     } else
    1525     7734132 :       z[1] = evalsigne(0) | evallgefint(2);
    1526    15294720 :     gel(c, i) = z;
    1527             :   }
    1528       93141 :   return c;
    1529             : }
    1530             : 
    1531             : static GEN
    1532       47313 : muliifft_unspliti(GEN V, long bs, long len)
    1533             : {
    1534       47313 :   long s, i, j, l = lg(V);
    1535       47313 :   GEN a = cgeti(len);
    1536       47313 :   a[1] = evalsigne(1)|evallgefint(len);
    1537   224698593 :   for(i=2;i<len;i++)
    1538   224651280 :     a[i] = 0;
    1539     7871697 :   for(i=1, s=0; i<l; i++, s+=bs)
    1540             :   {
    1541     7824384 :     GEN  u = gel(V,i);
    1542     7824384 :     if (signe(u))
    1543             :     {
    1544     7564203 :       GEN ap = int_W(a,s);
    1545     7564203 :       GEN up = int_LSW(u);
    1546     7564203 :       long lu = NLIMBS(u);
    1547             :       LOCAL_OVERFLOW;
    1548     7564203 :       *ap = addll(*ap, *up--); ap--;
    1549   442210293 :       for (j=1; j<lu; j++)
    1550   434646090 :        { *ap = addllx(*ap, *up--); ap--; }
    1551     7564344 :       while (overflow)
    1552         141 :        { *ap = addll(*ap, 1); ap--; }
    1553             :     }
    1554             :   }
    1555       47313 :   return int_normalize(a,0);
    1556             : }
    1557             : 
    1558             : static GEN
    1559        1485 : sqrispec_fft(GEN a, long na)
    1560             : {
    1561        1485 :   pari_sp av, ltop = avma;
    1562        1485 :   long len = 2*na;
    1563             :   long k, mod, bs, n;
    1564             :   GEN  FFT, M;
    1565             :   long i;
    1566             :   ulong o, ord;
    1567             : 
    1568        1485 :   mulliifft_params(len,&k,&mod,&bs,&n,&ord);
    1569        1485 :   o = ord>>k;
    1570        1485 :   M = int2n(mod*BITS_IN_LONG);
    1571        1485 :   M[2+mod] = 1;
    1572        1485 :   FFT = muliifft_spliti(a, na, bs, n, mod);
    1573        1485 :   muliifft_dit(o, ord, M, FFT, 0, n);
    1574        1485 :   av = avma;
    1575      355533 :   for(i=1; i<=n; i++)
    1576             :   {
    1577      354048 :     affii(Zf_sqr(gel(FFT,i), M), gel(FFT,i));
    1578      354048 :     set_avma(av);
    1579             :   }
    1580        1485 :   muliifft_dis(ord-o, ord, M, FFT, 0, n);
    1581      355533 :   for(i=1; i<=n; i++)
    1582             :   {
    1583      354048 :     affii(Zf_shift(gel(FFT,i), (ord>>1)-k, M), gel(FFT,i));
    1584      354048 :     set_avma(av);
    1585             :   }
    1586        1485 :   return gerepileuptoint(ltop, muliifft_unspliti(FFT,bs,2+len));
    1587             : }
    1588             : 
    1589             : static GEN
    1590       45828 : muliispec_fft(GEN a, GEN b, long na, long nb)
    1591             : {
    1592       45828 :   pari_sp av, av2, ltop = avma;
    1593       45828 :   long len = na+nb;
    1594             :   long k, mod, bs, n;
    1595             :   GEN FFT, FFTb, M;
    1596             :   long i;
    1597             :   ulong o, ord;
    1598             : 
    1599       45828 :   mulliifft_params(len,&k,&mod,&bs,&n,&ord);
    1600       45828 :   o = ord>>k;
    1601       45828 :   M = int2n(mod*BITS_IN_LONG);
    1602       45828 :   M[2+mod] = 1;
    1603       45828 :   FFT = muliifft_spliti(a, na, bs, n, mod);
    1604       45828 :   av=avma;
    1605       45828 :   muliifft_dit(o, ord, M, FFT, 0, n);
    1606       45828 :   FFTb = muliifft_spliti(b, nb, bs, n, mod);
    1607       45828 :   av2 = avma;
    1608       45828 :   muliifft_dit(o, ord, M, FFTb, 0, n);
    1609     7516164 :   for(i=1; i<=n; i++)
    1610             :   {
    1611     7470336 :     affii(Zf_mul(gel(FFT,i), gel(FFTb,i), M), gel(FFT,i));
    1612     7470336 :     set_avma(av2);
    1613             :   }
    1614       45828 :   set_avma(av);
    1615       45828 :   muliifft_dis(ord-o, ord, M, FFT, 0, n);
    1616     7516164 :   for(i=1; i<=n; i++)
    1617             :   {
    1618     7470336 :     affii(Zf_shift(gel(FFT,i),(ord>>1)-k,M), gel(FFT,i));
    1619     7470336 :     set_avma(av);
    1620             :   }
    1621       45828 :   return gerepileuptoint(ltop, muliifft_unspliti(FFT,bs,2+len));
    1622             : }
    1623             : 
    1624             : /********************************************************************/
    1625             : /**                                                                **/
    1626             : /**               INTEGER MULTIPLICATION (KARATSUBA)               **/
    1627             : /**                                                                **/
    1628             : /********************************************************************/
    1629             : 
    1630             : /* return (x shifted left d words) + y. Assume d > 0, x > 0 and y >= 0 */
    1631             : static GEN
    1632   300178500 : addshiftw(GEN x, GEN y, long d)
    1633             : {
    1634   300178500 :   GEN z,z0,y0,yd, zd = (GEN)avma;
    1635   300178500 :   long a,lz,ly = lgefint(y);
    1636             : 
    1637   300178500 :   z0 = new_chunk(d);
    1638   300178500 :   a = ly-2; yd = y+ly;
    1639   300178500 :   if (a >= d)
    1640             :   {
    1641  5498000319 :     y0 = yd-d; while (yd > y0) *--zd = *--yd; /* copy last d words of y */
    1642   297414342 :     a -= d;
    1643   297414342 :     if (a)
    1644   203545143 :       z = addiispec(LIMBS(x), LIMBS(y), NLIMBS(x), a);
    1645             :     else
    1646    93869199 :       z = icopy(x);
    1647             :   }
    1648             :   else
    1649             :   {
    1650    13076199 :     y0 = yd-a; while (yd > y0) *--zd = *--yd; /* copy last a words of y */
    1651    43216407 :     while (zd > z0) *--zd = 0;    /* complete with 0s */
    1652     2764158 :     z = icopy(x);
    1653             :   }
    1654   300178500 :   lz = lgefint(z)+d;
    1655   300178500 :   z[1] = evalsigne(1) | evallgefint(lz);
    1656   300178500 :   z[0] = evaltyp(t_INT) | evallg(lz); return z;
    1657             : }
    1658             : 
    1659             : /* Fast product (Karatsuba) of integers. a and b are "special" GENs
    1660             :  * c,c0,c1,c2 are genuine GENs.
    1661             :  */
    1662             : GEN
    1663  3664842978 : muliispec(GEN a, GEN b, long na, long nb)
    1664             : {
    1665             :   GEN a0,c,c0;
    1666             :   long n0, n0a, i;
    1667             :   pari_sp av;
    1668             : 
    1669  3664842978 :   if (na < nb) swapspec(a,b, na,nb);
    1670  3664842978 :   if (nb < MULII_KARATSUBA_LIMIT) return muliispec_basecase(a,b,na,nb);
    1671    98567106 :   if (nb >= MULII_FFT_LIMIT)      return muliispec_fft(a,b,na,nb);
    1672    98521278 :   i=(na>>1); n0=na-i; na=i;
    1673    98521278 :   av=avma; a0=a+na; n0a=n0;
    1674   165943356 :   while (n0a && !*a0) { a0++; n0a--; }
    1675             : 
    1676    98521278 :   if (n0a && nb > n0)
    1677    96233607 :   { /* nb <= na <= n0 */
    1678             :     GEN b0,c1,c2;
    1679             :     long n0b;
    1680             : 
    1681    96233607 :     nb -= n0;
    1682    96233607 :     c = muliispec(a,b,na,nb);
    1683    96233607 :     b0 = b+nb; n0b = n0;
    1684   160891947 :     while (n0b && !*b0) { b0++; n0b--; }
    1685    96233607 :     if (n0b)
    1686             :     {
    1687    95619840 :       c0 = muliispec(a0,b0, n0a,n0b);
    1688             : 
    1689    95619840 :       c2 = addiispec(a0,a, n0a,na);
    1690    95619840 :       c1 = addiispec(b0,b, n0b,nb);
    1691    95619840 :       c1 = muliispec(LIMBS(c1),LIMBS(c2), NLIMBS(c1),NLIMBS(c2));
    1692    95619840 :       c2 = addiispec(LIMBS(c0),LIMBS(c),  NLIMBS(c0),NLIMBS(c));
    1693             : 
    1694    95619840 :       c1 = subiispec(LIMBS(c1),LIMBS(c2), NLIMBS(c1),NLIMBS(c2));
    1695             :     }
    1696             :     else
    1697             :     {
    1698      613767 :       c0 = gen_0;
    1699      613767 :       c1 = muliispec(a0,b, n0a,nb);
    1700             :     }
    1701    96233607 :     c = addshiftw(c,c1, n0);
    1702             :   }
    1703             :   else
    1704             :   {
    1705     2287671 :     c = muliispec(a,b,na,nb);
    1706     2287671 :     c0 = muliispec(a0,b,n0a,nb);
    1707             :   }
    1708    98521278 :   return gerepileuptoint(av, addshiftw(c,c0, n0));
    1709             : }
    1710             : GEN
    1711      158481 : muluui(ulong x, ulong y, GEN z)
    1712             : {
    1713      158481 :   long t, s = signe(z);
    1714             :   GEN r;
    1715             :   LOCAL_HIREMAINDER;
    1716             : 
    1717      158481 :   if (!x || !y || !signe(z)) return gen_0;
    1718      158184 :   t = mulll(x,y);
    1719      158184 :   if (!hiremainder)
    1720      158184 :     r = muluispec(t, z+2, lgefint(z)-2);
    1721             :   else
    1722             :   {
    1723             :     long tmp[2];
    1724           0 :     tmp[0] = hiremainder;
    1725           0 :     tmp[1] = t;
    1726           0 :     r = muliispec(z+2,tmp,lgefint(z)-2,2);
    1727             :   }
    1728      158184 :   setsigne(r,s); return r;
    1729             : }
    1730             : 
    1731             : #define sqrispec_mirror sqrispec
    1732             : #define muliispec_mirror muliispec
    1733             : 
    1734             : /* x % (2^n), assuming n >= 0 */
    1735             : GEN
    1736    19549413 : remi2n(GEN x, long n)
    1737             : {
    1738    19549413 :   long hi,l,k,lx,ly, sx = signe(x);
    1739             :   GEN z, xd, zd;
    1740             : 
    1741    19549413 :   if (!sx || !n) return gen_0;
    1742             : 
    1743    19520337 :   k = dvmdsBIL(n, &l);
    1744    19520337 :   lx = lgefint(x);
    1745    19520337 :   if (lx < k+3) return icopy(x);
    1746             : 
    1747    19422504 :   xd = x + (lx-k-1);
    1748             :   /* x = |_|...|#|1|...|k| : copy the last l bits of # and the last k words
    1749             :    *            ^--- initial xd  */
    1750    19422504 :   hi = ((ulong)*xd) & ((1UL<<l)-1); /* last l bits of # = top bits of result */
    1751    19422504 :   if (!hi)
    1752             :   { /* strip leading zeroes from result */
    1753      520161 :     xd++; while (k && !*xd) { k--; xd++; }
    1754      511872 :     if (!k) return gen_0;
    1755      314022 :     ly = k+2; xd--;
    1756             :   }
    1757             :   else
    1758    18910632 :     ly = k+3;
    1759             : 
    1760    19224654 :   zd = z = cgeti(ly);
    1761    19224654 :   *++zd = evalsigne(sx) | evallgefint(ly);
    1762    19224654 :   if (hi) *++zd = hi;
    1763    81904368 :   for ( ;k; k--) *++zd = *++xd;
    1764    19224654 :   return z;
    1765             : }
    1766             : 
    1767             : GEN
    1768   714612096 : sqrispec(GEN a, long na)
    1769             : {
    1770             :   GEN a0,c;
    1771             :   long n0, n0a, i;
    1772             :   pari_sp av;
    1773             : 
    1774   714612096 :   if (na < SQRI_KARATSUBA_LIMIT) return sqrispec_basecase(a,na);
    1775     5456961 :   if (na >= SQRI_FFT_LIMIT) return sqrispec_fft(a,na);
    1776     5455476 :   i=(na>>1); n0=na-i; na=i;
    1777     5455476 :   av=avma; a0=a+na; n0a=n0;
    1778     7667889 :   while (n0a && !*a0) { a0++; n0a--; }
    1779     5455476 :   c = sqrispec(a,na);
    1780     5455476 :   if (n0a)
    1781             :   {
    1782     5449380 :     GEN t, c1, c0 = sqrispec(a0,n0a);
    1783             : #if 0
    1784             :     c1 = shifti(muliispec(a0,a, n0a,na),1);
    1785             : #else /* faster */
    1786     5449380 :     t = addiispec(a0,a,n0a,na);
    1787     5449380 :     t = sqrispec(LIMBS(t),NLIMBS(t));
    1788     5449380 :     c1= addiispec(LIMBS(c0),LIMBS(c), NLIMBS(c0), NLIMBS(c));
    1789     5449380 :     c1= subiispec(LIMBS(t),LIMBS(c1), NLIMBS(t), NLIMBS(c1));
    1790             : #endif
    1791     5449380 :     c = addshiftw(c,c1, n0);
    1792     5449380 :     c = addshiftw(c,c0, n0);
    1793             :   }
    1794             :   else
    1795        6096 :     c = addshiftw(c,gen_0,n0<<1);
    1796     5455476 :   return gerepileuptoint(av, c);
    1797             : }
    1798             : 
    1799             : /********************************************************************/
    1800             : /**                                                                **/
    1801             : /**                    KARATSUBA SQUARE ROOT                       **/
    1802             : /**      adapted from Paul Zimmermann's implementation of          **/
    1803             : /**      his algorithm in GMP (mpn_sqrtrem)                        **/
    1804             : /**                                                                **/
    1805             : /********************************************************************/
    1806             : 
    1807             : /* Square roots table */
    1808             : static const unsigned char approx_tab[192] = {
    1809             :   128,128,129,130,131,132,133,134,135,136,137,138,139,140,141,142,
    1810             :   143,144,144,145,146,147,148,149,150,150,151,152,153,154,155,155,
    1811             :   156,157,158,159,160,160,161,162,163,163,164,165,166,167,167,168,
    1812             :   169,170,170,171,172,173,173,174,175,176,176,177,178,178,179,180,
    1813             :   181,181,182,183,183,184,185,185,186,187,187,188,189,189,190,191,
    1814             :   192,192,193,193,194,195,195,196,197,197,198,199,199,200,201,201,
    1815             :   202,203,203,204,204,205,206,206,207,208,208,209,209,210,211,211,
    1816             :   212,212,213,214,214,215,215,216,217,217,218,218,219,219,220,221,
    1817             :   221,222,222,223,224,224,225,225,226,226,227,227,228,229,229,230,
    1818             :   230,231,231,232,232,233,234,234,235,235,236,236,237,237,238,238,
    1819             :   239,240,240,241,241,242,242,243,243,244,244,245,245,246,246,247,
    1820             :   247,248,248,249,249,250,250,251,251,252,252,253,253,254,254,255
    1821             : };
    1822             : 
    1823             : /* N[0], assume N[0] >= 2^(BIL-2).
    1824             :  * Return r,s such that s^2 + r = N, 0 <= r <= 2s */
    1825             : static void
    1826    49669809 : p_sqrtu1(ulong *N, ulong *ps, ulong *pr)
    1827             : {
    1828    49669809 :   ulong prec, r, s, q, u, n0 = N[0];
    1829             : 
    1830    49669809 :   q = n0 >> (BITS_IN_LONG - 8);
    1831             :   /* 2^6 = 64 <= q < 256 = 2^8 */
    1832    49669809 :   s = approx_tab[q - 64];                                /* 128 <= s < 255 */
    1833    49669809 :   r = (n0 >> (BITS_IN_LONG - 16)) - s * s;                /* r <= 2*s */
    1834    49669809 :   if (r > (s << 1)) { r -= (s << 1) | 1; s++; }
    1835             : 
    1836             :   /* 8-bit approximation from the high 8-bits of N[0] */
    1837    49669809 :   prec = 8;
    1838    49669809 :   n0 <<= 2 * prec;
    1839   149009427 :   while (2 * prec < BITS_IN_LONG)
    1840             :   { /* invariant: s has prec bits, and r <= 2*s */
    1841    99339618 :     r = (r << prec) + (n0 >> (BITS_IN_LONG - prec));
    1842    99339618 :     n0 <<= prec;
    1843    99339618 :     u = 2 * s;
    1844    99339618 :     q = r / u; u = r - q * u;
    1845    99339618 :     s = (s << prec) + q;
    1846    99339618 :     u = (u << prec) + (n0 >> (BITS_IN_LONG - prec));
    1847    99339618 :     q = q * q;
    1848    99339618 :     r = u - q;
    1849    99339618 :     if (u < q) { s--; r += (s << 1) | 1; }
    1850    99339618 :     n0 <<= prec;
    1851    99339618 :     prec = 2 * prec;
    1852             :   }
    1853    49669809 :   *ps = s;
    1854    49669809 :   *pr = r;
    1855    49669809 : }
    1856             : 
    1857             : /* N[0..1], assume N[0] >= 2^(BIL-2).
    1858             :  * Return 1 if remainder overflows, 0 otherwise */
    1859             : static int
    1860    48053004 : p_sqrtu2(ulong *N, ulong *ps, ulong *pr)
    1861             : {
    1862    48053004 :   ulong cc, qhl, r, s, q, u, n1 = N[1];
    1863             :   LOCAL_OVERFLOW;
    1864             : 
    1865    48053004 :   p_sqrtu1(N, &s, &r); /* r <= 2s */
    1866    72419295 :   qhl = 0; while (r >= s) { qhl++; r -= s; }
    1867             :   /* now r < s < 2^(BIL/2) */
    1868    48053004 :   r = (r << BITS_IN_HALFULONG) | (n1 >> BITS_IN_HALFULONG);
    1869    48053004 :   u = s << 1;
    1870    48053004 :   q = r / u; u = r - q * u;
    1871    48053004 :   q += (qhl & 1) << (BITS_IN_HALFULONG - 1);
    1872    48053004 :   qhl >>= 1;
    1873             :   /* (initial r)<<(BIL/2) + n1>>(BIL/2) = (qhl<<(BIL/2) + q) * 2s + u */
    1874    48053004 :   s = ((s + qhl) << BITS_IN_HALFULONG) + q;
    1875    48053004 :   cc = u >> BITS_IN_HALFULONG;
    1876    48053004 :   r = (u << BITS_IN_HALFULONG) | (n1 & LOWMASK);
    1877    48053004 :   r = subll(r, q * q);
    1878    48053004 :   cc -= overflow + qhl;
    1879             :   /* now subtract 2*q*2^(BIL/2) + 2^BIL if qhl is set */
    1880    48053004 :   if ((long)cc < 0)
    1881             :   {
    1882    12846135 :     if (s) {
    1883    12825423 :       r = addll(r, s);
    1884    12825423 :       cc += overflow;
    1885    12825423 :       s--;
    1886             :     } else {
    1887       20712 :       cc++;
    1888       20712 :       s = ~0UL;
    1889             :     }
    1890    12846135 :     r = addll(r, s);
    1891    12846135 :     cc += overflow;
    1892             :   }
    1893    48053004 :   *ps = s;
    1894    48053004 :   *pr = r; return cc;
    1895             : }
    1896             : 
    1897             : static void
    1898    44001660 : xmpn_zero(GEN x, long n)
    1899             : {
    1900   274415106 :   while (--n >= 0) x[n]=0;
    1901    44001660 : }
    1902             : static void
    1903   441339198 : xmpn_copy(GEN z, GEN x, long n)
    1904             : {
    1905   441339198 :   long k = n;
    1906  1601322595 :   while (--k >= 0) z[k] = x[k];
    1907   441339198 : }
    1908             : /* a[0..la-1] * 2^(lb BIL) | b[0..lb-1] */
    1909             : static GEN
    1910   189037518 : catii(GEN a, long la, GEN b, long lb)
    1911             : {
    1912   189037518 :   long l = la + lb + 2;
    1913   189037518 :   GEN z = cgetipos(l);
    1914   189037518 :   xmpn_copy(LIMBS(z), a, la);
    1915   189037518 :   xmpn_copy(LIMBS(z) + la, b, lb);
    1916   189037518 :   return int_normalize(z, 0);
    1917             : }
    1918             : 
    1919             : /* sqrt n[0..1], assume n normalized */
    1920             : static GEN
    1921    46478655 : sqrtispec2(GEN n, GEN *pr)
    1922             : {
    1923             :   ulong s, r;
    1924    46478655 :   int hi = p_sqrtu2((ulong*)n, &s, &r);
    1925    46478655 :   GEN S = utoi(s);
    1926    46478655 :   *pr = hi? uutoi(1,r): utoi(r);
    1927    46478655 :   return S;
    1928             : }
    1929             : 
    1930             : /* sqrt n[0], _dont_ assume n normalized */
    1931             : static GEN
    1932     1616805 : sqrtispec1_sh(GEN n, GEN *pr)
    1933             : {
    1934             :   GEN S;
    1935     1616805 :   ulong r, s, u0 = uel(n,0);
    1936     1616805 :   int sh = bfffo(u0) & ~1UL;
    1937     1616805 :   if (sh) u0 <<= sh;
    1938     1616805 :   p_sqrtu1(&u0, &s, &r);
    1939             :   /* s^2 + r = u0, s < 2^(BIL/2). Rescale back:
    1940             :    * 2^(2k) n = S^2 + R
    1941             :    * so 2^(2k) n = (S - s0)^2 + (2*S*s0 - s0^2 + R), s0 = S mod 2^k. */
    1942     1616805 :   if (sh) {
    1943     1355205 :     int k = sh >> 1;
    1944     1355205 :     ulong s0 = s & ((1L<<k) - 1);
    1945     1355205 :     r += s * (s0<<1);
    1946     1355205 :     s >>= k;
    1947     1355205 :     r >>= sh;
    1948             :   }
    1949     1616805 :   S = utoi(s);
    1950     1616805 :   if (pr) *pr = utoi(r);
    1951     1616805 :   return S;
    1952             : }
    1953             : 
    1954             : /* sqrt n[0..1], _dont_ assume n normalized */
    1955             : static GEN
    1956     1574349 : sqrtispec2_sh(GEN n, GEN *pr)
    1957             : {
    1958             :   GEN S;
    1959     1574349 :   ulong U[2], r, s, u0 = uel(n,0), u1 = uel(n,1);
    1960     1574349 :   int hi, sh = bfffo(u0) & ~1UL;
    1961     1574349 :   if (sh) {
    1962     1332993 :     u0 = (u0 << sh) | (u1 >> (BITS_IN_LONG-sh));
    1963     1332993 :     u1 <<= sh;
    1964             :   }
    1965     1574349 :   U[0] = u0;
    1966     1574349 :   U[1] = u1; hi = p_sqrtu2(U, &s, &r);
    1967             :   /* s^2 + R = u0|u1. Rescale back:
    1968             :    * 2^(2k) n = S^2 + R
    1969             :    * so 2^(2k) n = (S - s0)^2 + (2*S*s0 - s0^2 + R), s0 = S mod 2^k. */
    1970     1574349 :   if (sh) {
    1971     1332993 :     int k = sh >> 1;
    1972     1332993 :     ulong s0 = s & ((1L<<k) - 1);
    1973             :     LOCAL_HIREMAINDER;
    1974             :     LOCAL_OVERFLOW;
    1975     1332993 :     r = addll(r, mulll(s, (s0<<1)));
    1976     1332993 :     if (overflow) hiremainder++;
    1977     1332993 :     hiremainder += hi; /* + 0 or 1 */
    1978     1332993 :     s >>= k;
    1979     1332993 :     r = (r>>sh) | (hiremainder << (BITS_IN_LONG-sh));
    1980     1332993 :     hi = (hiremainder & (1L<<sh));
    1981             :   }
    1982     1574349 :   S = utoi(s);
    1983     1574349 :   if (pr) *pr = hi? uutoi(1,r): utoi(r);
    1984     1574349 :   return S;
    1985             : }
    1986             : 
    1987             : /* Let N = N[0..2n-1]. Return S (and set R) s.t S^2 + R = N, 0 <= R <= 2S
    1988             :  * Assume N normalized */
    1989             : static GEN
    1990   140997414 : sqrtispec(GEN N, long n, GEN *r)
    1991             : {
    1992             :   GEN S, R, q, z, u;
    1993             :   long l, h;
    1994             : 
    1995   140997414 :   if (n == 1) return sqrtispec2(N, r);
    1996    94518759 :   l = n >> 1;
    1997    94518759 :   h = n - l; /* N = a3(h) | a2(h) | a1(l) | a0(l words) */
    1998    94518759 :   S = sqrtispec(N, h, &R); /* S^2 + R = a3|a2 */
    1999             : 
    2000    94518759 :   z = catii(LIMBS(R), NLIMBS(R), N + 2*h, l); /* = R | a1(l) */
    2001    94518759 :   q = dvmdii(z, shifti(S,1), &u);
    2002    94518759 :   z = catii(LIMBS(u), NLIMBS(u), N + n + h, l); /* = u | a0(l) */
    2003             : 
    2004    94518759 :   S = addshiftw(S, q, l);
    2005    94518759 :   R = subii(z, sqri(q));
    2006    94518759 :   if (signe(R) < 0)
    2007             :   {
    2008    16423654 :     GEN S2 = shifti(S,1);
    2009    16423654 :     R = addis(subiispec(LIMBS(S2),LIMBS(R), NLIMBS(S2),NLIMBS(R)), -1);
    2010    16423654 :     S = addis(S, -1);
    2011             :   }
    2012    94518759 :   *r = R; return S;
    2013             : }
    2014             : 
    2015             : /* Return S (and set R) s.t S^2 + R = N, 0 <= R <= 2S.
    2016             :  * As for dvmdii, R is last on stack and guaranteed to be gen_0 in case the
    2017             :  * remainder is 0. R = NULL is allowed. */
    2018             : GEN
    2019     5668692 : sqrtremi(GEN N, GEN *r)
    2020             : {
    2021             :   pari_sp av;
    2022     5668692 :   GEN S, R, n = N+2;
    2023     5668692 :   long k, l2, ln = NLIMBS(N);
    2024             :   int sh;
    2025             : 
    2026     5668692 :   if (ln <= 2)
    2027             :   {
    2028     3191697 :     if (ln == 2) return sqrtispec2_sh(n, r);
    2029     1617348 :     if (ln == 1) return sqrtispec1_sh(n, r);
    2030         543 :     if (r) *r = gen_0;
    2031         543 :     return gen_0;
    2032             :   }
    2033     2476995 :   av = avma;
    2034     2476995 :   sh = bfffo(n[0]) >> 1;
    2035     2476995 :   l2 = (ln + 1) >> 1;
    2036     2476995 :   if (sh || (ln & 1)) { /* normalize n, so that n[0] >= 2^BIL / 4 */
    2037     2476605 :     GEN s0, t = new_chunk(ln + 1);
    2038     2476605 :     t[ln] = 0;
    2039     2476605 :     if (sh)
    2040     2475975 :       shift_left(t, n, 0,ln-1, 0, sh << 1);
    2041             :     else
    2042         630 :       xmpn_copy(t, n, ln);
    2043     2476605 :     S = sqrtispec(t, l2, &R); /* t normalized, 2 * l2 words */
    2044             :     /* Rescale back:
    2045             :      * 2^(2k) n = S^2 + R, k = sh + (ln & 1)*BIL/2
    2046             :      * so 2^(2k) n = (S - s0)^2 + (2*S*s0 - s0^2 + R), s0 = S mod 2^k. */
    2047     2476605 :     k = sh + (ln & 1) * (BITS_IN_LONG/2);
    2048     2476605 :     s0 = remi2n(S, k);
    2049     2476605 :     R = addii(shifti(R,-1), mulii(s0, S));
    2050     2476605 :     R = shifti(R, 1 - (k<<1));
    2051     2476605 :     S = shifti(S, -k);
    2052             :   }
    2053             :   else
    2054         390 :     S = sqrtispec(n, l2, &R);
    2055             : 
    2056     2476995 :   if (!r) { set_avma((pari_sp)S); return gerepileuptoint(av, S); }
    2057     2430195 :   gerepileall(av, 2, &S, &R); *r = R; return S;
    2058             : }
    2059             : 
    2060             : /* compute sqrt(|a|), assuming a != 0 */
    2061             : 
    2062             : #if 1
    2063             : GEN
    2064    44001660 : sqrtr_abs(GEN x)
    2065             : {
    2066    44001660 :   long l = realprec(x) - 2, e = expo(x), er = e>>1;
    2067    44001660 :   GEN b, c, res = cgetr(2 + l);
    2068    44001660 :   res[1] = evalsigne(1) | evalexpo(er);
    2069    44001660 :   if (e&1) {
    2070    19261872 :     b = new_chunk(l << 1);
    2071    19261872 :     xmpn_copy(b, x+2, l);
    2072    19261872 :     xmpn_zero(b + l,l);
    2073    19261872 :     b = sqrtispec(b, l, &c);
    2074    19261872 :     xmpn_copy(res+2, b+2, l);
    2075    19261872 :     if (cmpii(c, b) > 0) roundr_up_ip(res, l+2);
    2076             :   } else {
    2077             :     ulong u;
    2078    24739788 :     b = new_chunk(2 + (l << 1));
    2079    24739788 :     shift_left(b+1, x+2, 0,l-1, 0, BITS_IN_LONG-1);
    2080    24739788 :     b[0] = uel(x,2)>>1;
    2081    24739788 :     xmpn_zero(b + l+1,l+1);
    2082    24739788 :     b = sqrtispec(b, l+1, &c);
    2083    24739788 :     xmpn_copy(res+2, b+2, l);
    2084    24739788 :     u = uel(b,l+2);
    2085    24739788 :     if ( u&HIGHBIT || (u == ~HIGHBIT && cmpii(c,b) > 0))
    2086    12155447 :       roundr_up_ip(res, l+2);
    2087             :   }
    2088    44001660 :   return gc_const((pari_sp)res, res);
    2089             : }
    2090             : 
    2091             : #else /* use t_REAL: currently much slower (quadratic division) */
    2092             : 
    2093             : #ifdef LONG_IS_64BIT
    2094             : /* 64 bits of b = sqrt(a[0] * 2^64 + a[1])  [ up to 1ulp ] */
    2095             : static ulong
    2096             : sqrtu2(ulong *a)
    2097             : {
    2098             :   ulong c, b = dblmantissa( sqrt((double)a[0]) );
    2099             :   LOCAL_HIREMAINDER;
    2100             :   LOCAL_OVERFLOW;
    2101             : 
    2102             :   /* > 32 correct bits, 1 Newton iteration to reach 64 */
    2103             :   if (b <= a[0]) return HIGHBIT | (a[0] >> 1);
    2104             :   hiremainder = a[0]; c = divll(a[1], b);
    2105             :   return (addll(c, b) >> 1) | HIGHBIT;
    2106             : }
    2107             : /* 64 bits of sqrt(a[0] * 2^63) */
    2108             : static ulong
    2109             : sqrtu2_1(ulong *a)
    2110             : {
    2111             :   ulong t[2];
    2112             :   t[0] = (a[0] >> 1);
    2113             :   t[1] = (a[0] << (BITS_IN_LONG-1)) | (a[1] >> 1);
    2114             :   return sqrtu2(t);
    2115             : }
    2116             : #else
    2117             : /* 32 bits of sqrt(a[0] * 2^32) */
    2118             : static ulong
    2119             : sqrtu2(ulong *a)   { return dblmantissa( sqrt((double)a[0]) ); }
    2120             : /* 32 bits of sqrt(a[0] * 2^31) */
    2121             : static ulong
    2122             : sqrtu2_1(ulong *a) { return dblmantissa( sqrt(2. * a[0]) ); }
    2123             : #endif
    2124             : 
    2125             : GEN
    2126             : sqrtr_abs(GEN x)
    2127             : {
    2128             :   long l1, i, l = lg(x), ex = expo(x);
    2129             :   GEN a, t, y = cgetr(l);
    2130             :   pari_sp av, av0 = avma;
    2131             : 
    2132             :   a = rtor(x, l+1);
    2133             :   t = cgetr(l+1);
    2134             :   if (ex & 1) { /* odd exponent */
    2135             :     a[1] = evalsigne(1) | _evalexpo(1);
    2136             :     t[2] = (long)sqrtu2((ulong*)a + 2);
    2137             :   } else { /* even exponent */
    2138             :     a[1] = evalsigne(1) | _evalexpo(0);
    2139             :     t[2] = (long)sqrtu2_1((ulong*)a + 2);
    2140             :   }
    2141             :   t[1] = evalsigne(1) | _evalexpo(0);
    2142             :   for (i = 3; i <= l; i++) t[i] = 0;
    2143             : 
    2144             :   /* |x| = 2^(ex/2) a, t ~ sqrt(a) */
    2145             :   l--; l1 = 1; av = avma;
    2146             :   while (l1 < l) { /* let t := (t + a/t)/2 */
    2147             :     l1 <<= 1; if (l1 > l) l1 = l;
    2148             :     setlg(a, l1 + 2);
    2149             :     setlg(t, l1 + 2);
    2150             :     affrr(addrr(t, divrr(a,t)), t); shiftr_inplace(t, -1);
    2151             :     set_avma(av);
    2152             :   }
    2153             :   affrr(t,y); shiftr_inplace(y, (ex>>1));
    2154             :   return gc_const(av0, y);
    2155             : }
    2156             : 
    2157             : #endif
    2158             : 
    2159             : /*******************************************************************
    2160             :  *                                                                 *
    2161             :  *                           Base Conversion                       *
    2162             :  *                                                                 *
    2163             :  *******************************************************************/
    2164             : 
    2165             : static void
    2166      700431 : convi_dac(GEN x, ulong l, ulong *res)
    2167             : {
    2168      700431 :   pari_sp ltop=avma;
    2169             :   ulong m;
    2170             :   GEN x1,x2;
    2171      700431 :   if (l==1) { *res=itou(x); return; }
    2172      332517 :   m=l>>1;
    2173      332517 :   x1=dvmdii(x,powuu(1000000000UL,m),&x2);
    2174      332517 :   convi_dac(x1,l-m,res+m);
    2175      332517 :   convi_dac(x2,m,res);
    2176      332517 :   set_avma(ltop);
    2177             : }
    2178             : 
    2179             : /* convert integer --> base 10^9 [not memory clean] */
    2180             : ulong *
    2181      297354 : convi(GEN x, long *l)
    2182             : {
    2183      297354 :   long lz, lx = lgefint(x);
    2184             :   ulong *z;
    2185      297354 :   if (lx == 3 && uel(x,2) < 1000000000UL) {
    2186      261957 :     z = (ulong*)new_chunk(1);
    2187      261957 :     *z = x[2];
    2188      261957 :     *l = 1; return z+1;
    2189             :   }
    2190       35397 :   lz = 1 + (long)bit_accuracy_mul(lx, LOG10_2/9);
    2191       35397 :   z = (ulong*)new_chunk(lz);
    2192       35397 :   convi_dac(x,(ulong)lz,z);
    2193       64338 :   while (z[lz-1]==0) lz--;
    2194       35397 :   *l=lz; return z+lz;
    2195             : }
    2196             : 

Generated by: LCOV version 1.13