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/gmp - mp.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.12.1 lcov report (development 24988-2584e74448) Lines: 646 714 90.5 %
Date: 2020-01-26 05:57:03 Functions: 54 55 98.2 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : #line 2 "../src/kernel/gmp/mp.c"
       2             : /* Copyright (C) 2002-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. It is distributed in the hope that it will be useful, but WITHOUT
       9             : ANY WARRANTY WHATSOEVER.
      10             : 
      11             : Check the License for details. You should have received a copy of it, along
      12             : with the package; see the file 'COPYING'. If not, write to the Free Software
      13             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      14             : 
      15             : /***********************************************************************/
      16             : /**                                                                   **/
      17             : /**                               GMP KERNEL                          **/
      18             : /** BA2002Sep24                                                       **/
      19             : /***********************************************************************/
      20             : /* GMP t_INT as just like normal t_INT, just the mantissa is the other way
      21             :  * round
      22             :  *
      23             :  *   `How would you like to live in Looking-glass House, Kitty?  I
      24             :  *   wonder if they'd give you milk in there?  Perhaps Looking-glass
      25             :  *   milk isn't good to drink--But oh, Kitty! now we come to the
      26             :  *   passage.  You can just see a little PEEP of the passage in
      27             :  *   Looking-glass House, if you leave the door of our drawing-room
      28             :  *   wide open:  and it's very like our passage as far as you can see,
      29             :  *   only you know it may be quite different on beyond.  Oh, Kitty!
      30             :  *   how nice it would be if we could only get through into Looking-
      31             :  *   glass House!  I'm sure it's got, oh! such beautiful things in it!
      32             :  *
      33             :  *  Through the Looking Glass,  Lewis Carrol
      34             :  *
      35             :  *  (pityful attempt to beat GN code/comments rate)
      36             :  *  */
      37             : 
      38             : #include <gmp.h>
      39             : #include "pari.h"
      40             : #include "paripriv.h"
      41             : #include "../src/kernel/none/tune-gen.h"
      42             : 
      43             : /*We need PARI invmod renamed to invmod_pari*/
      44             : #define INVMOD_PARI
      45             : 
      46           0 : static void *pari_gmp_realloc(void *ptr, size_t old_size, size_t new_size) {
      47           0 :   (void)old_size; return (void *) pari_realloc(ptr,new_size);
      48             : }
      49             : 
      50      310807 : static void pari_gmp_free(void *ptr, size_t old_size){
      51      310807 :   (void)old_size; pari_free(ptr);
      52      310807 : }
      53             : 
      54             : static void *(*old_gmp_malloc)(size_t new_size);
      55             : static void *(*old_gmp_realloc)(void *ptr, size_t old_size, size_t new_size);
      56             : static void (*old_gmp_free)(void *ptr, size_t old_size);
      57             : 
      58             : void
      59         940 : pari_kernel_init(void)
      60             : {
      61         940 :   mp_get_memory_functions (&old_gmp_malloc, &old_gmp_realloc, &old_gmp_free);
      62         940 :   mp_set_memory_functions((void *(*)(size_t)) pari_malloc, pari_gmp_realloc, pari_gmp_free);
      63         940 : }
      64             : 
      65             : void
      66         932 : pari_kernel_close(void)
      67             : {
      68             :   void *(*new_gmp_malloc)(size_t new_size);
      69             :   void *(*new_gmp_realloc)(void *ptr, size_t old_size, size_t new_size);
      70             :   void (*new_gmp_free)(void *ptr, size_t old_size);
      71         932 :   mp_get_memory_functions (&new_gmp_malloc, &new_gmp_realloc, &new_gmp_free);
      72         932 :   if (new_gmp_malloc==pari_malloc) new_gmp_malloc = old_gmp_malloc;
      73         932 :   if (new_gmp_realloc==pari_gmp_realloc) new_gmp_realloc = old_gmp_realloc;
      74         932 :   if (new_gmp_free==pari_gmp_free) new_gmp_free = old_gmp_free;
      75         932 :   mp_set_memory_functions(new_gmp_malloc, new_gmp_realloc, new_gmp_free);
      76         932 : }
      77             : 
      78             : #define LIMBS(x)  ((mp_limb_t *)((x)+2))
      79             : #define NLIMBS(x) (lgefint(x)-2)
      80             : /*This one is for t_REALs to emphasize they are not t_INTs*/
      81             : #define RLIMBS(x)  ((mp_limb_t *)((x)+2))
      82             : #define RNLIMBS(x) (lg(x)-2)
      83             : 
      84             : INLINE void
      85     2577704 : xmpn_copy(mp_limb_t *x, mp_limb_t *y, long n)
      86             : {
      87     2577704 :   while (--n >= 0) x[n]=y[n];
      88     2577704 : }
      89             : 
      90             : INLINE void
      91   249940390 : xmpn_mirror(mp_limb_t *x, long n)
      92             : {
      93             :   long i;
      94  1484017602 :   for(i=0;i<(n>>1);i++)
      95             :   {
      96  1234077212 :     ulong m=x[i];
      97  1234077212 :     x[i]=x[n-1-i];
      98  1234077212 :     x[n-1-i]=m;
      99             :   }
     100   249940390 : }
     101             : 
     102             : INLINE void
     103   318900419 : xmpn_mirrorcopy(mp_limb_t *z, mp_limb_t *x, long n)
     104             : {
     105             :   long i;
     106  3653217757 :   for(i=0;i<n;i++)
     107  3334317338 :     z[i]=x[n-1-i];
     108   318900419 : }
     109             : 
     110             : INLINE void
     111    95548396 : xmpn_zero(mp_ptr x, mp_size_t n)
     112             : {
     113    95548396 :   while (--n >= 0) x[n]=0;
     114    95548396 : }
     115             : 
     116             : INLINE GEN
     117    24346167 : icopy_ef(GEN x, long l)
     118             : {
     119    24346167 :   register long lx = lgefint(x);
     120    24346167 :   const GEN y = cgeti(l);
     121             : 
     122    24345775 :   while (--lx > 0) y[lx]=x[lx];
     123    24345775 :   return y;
     124             : }
     125             : 
     126             : /* NOTE: arguments of "spec" routines (muliispec, addiispec, etc.) aren't
     127             :  * GENs but pairs (long *a, long na) representing a list of digits (in basis
     128             :  * BITS_IN_LONG) : a[0], ..., a[na-1]. [ In ordre to facilitate splitting: no
     129             :  * need to reintroduce codewords ]
     130             :  * Use speci(a,na) to visualize the corresponding GEN.
     131             :  */
     132             : 
     133             : /***********************************************************************/
     134             : /**                                                                   **/
     135             : /**                     ADDITION / SUBTRACTION                        **/
     136             : /**                                                                   **/
     137             : /***********************************************************************/
     138             : 
     139             : GEN
     140     2980154 : setloop(GEN a)
     141             : {
     142     2980154 :   pari_sp av = avma - 2 * sizeof(long);
     143     2980154 :   (void)cgetg(lgefint(a) + 3, t_VECSMALL);
     144     2980154 :   return icopy_avma(a, av); /* two cells of extra space after a */
     145             : }
     146             : 
     147             : /* we had a = setloop(?), then some incloops. Reset a to b */
     148             : GEN
     149      170696 : resetloop(GEN a, GEN b) {
     150      170696 :   a[0] = evaltyp(t_INT) | evallg(lgefint(b));
     151      170696 :   affii(b, a); return a;
     152             : }
     153             : 
     154             : /* assume a > 0, initialized by setloop. Do a++ */
     155             : static GEN
     156    31807598 : incpos(GEN a)
     157             : {
     158    31807598 :   long i, l = lgefint(a);
     159    31807603 :   for (i=2; i<l; i++)
     160    31807599 :     if (++uel(a,i)) return a;
     161           4 :   a[l] = 1; l++;
     162           4 :   a[0]=evaltyp(t_INT) | _evallg(l);
     163           4 :   a[1]=evalsigne(1) | evallgefint(l);
     164           4 :   return a;
     165             : }
     166             : 
     167             : /* assume a < 0, initialized by setloop. Do a++ */
     168             : static GEN
     169       10952 : incneg(GEN a)
     170             : {
     171       10952 :   long i, l = lgefint(a)-1;
     172       10952 :   if (uel(a,2)--)
     173             :   {
     174       10948 :     if (!a[l]) /* implies l = 2 */
     175             :     {
     176         604 :       a[0] = evaltyp(t_INT) | _evallg(2);
     177         604 :       a[1] = evalsigne(0) | evallgefint(2);
     178             :     }
     179       10948 :     return a;
     180             :   }
     181           5 :   for (i=3; i<=l; i++)
     182           5 :     if (uel(a,i)--) break;
     183           4 :   if (!a[l])
     184             :   {
     185           4 :     a[0] = evaltyp(t_INT) | _evallg(l);
     186           4 :     a[1] = evalsigne(-1) | evallgefint(l);
     187             :   }
     188           4 :   return a;
     189             : }
     190             : 
     191             : /* assume a initialized by setloop. Do a++ */
     192             : GEN
     193    32135100 : incloop(GEN a)
     194             : {
     195    32135100 :   switch(signe(a))
     196             :   {
     197             :     case 0:
     198      330856 :       a[0]=evaltyp(t_INT) | _evallg(3);
     199      330856 :       a[1]=evalsigne(1) | evallgefint(3);
     200      330856 :       a[2]=1; return a;
     201       10952 :     case -1: return incneg(a);
     202    31793292 :     default: return incpos(a);
     203             :   }
     204             : }
     205             : 
     206             : INLINE GEN
     207  1393871903 : adduispec(ulong s, GEN x, long nx)
     208             : {
     209             :   GEN  zd;
     210             :   long lz;
     211             : 
     212  1393871903 :   if (nx == 1) return adduu(uel(x,0), s);
     213   237319989 :   lz = nx+3; zd = cgeti(lz);
     214   235845853 :   if (mpn_add_1(LIMBS(zd),(mp_limb_t *)x,nx,s))
     215     4974894 :     zd[lz-1]=1;
     216             :   else
     217   230902345 :     lz--;
     218   235877239 :   zd[1] = evalsigne(1) | evallgefint(lz);
     219   235877239 :   return zd;
     220             : }
     221             : 
     222             : GEN
     223   278396801 : adduispec_offset(ulong s, GEN x, long offset, long nx)
     224             : {
     225   278396801 :   GEN xd=x+2+offset;
     226   278396801 :   while (nx && *(xd+nx-1)==0) nx--;
     227   278396801 :   if (!nx) return utoi(s);
     228   258039422 :   return adduispec(s,xd,nx);
     229             : }
     230             : 
     231             : INLINE GEN
     232  2012818544 : addiispec(GEN x, GEN y, long nx, long ny)
     233             : {
     234             :   GEN zd;
     235             :   long lz;
     236             : 
     237  2012818544 :   if (nx < ny) swapspec(x,y, nx,ny);
     238  2012818544 :   if (ny == 1) return adduispec(*y,x,nx);
     239   942263308 :   lz = nx+3; zd = cgeti(lz);
     240             : 
     241   937932225 :   if (mpn_add(LIMBS(zd),(mp_limb_t *)x,nx,(mp_limb_t *)y,ny))
     242    12285343 :     zd[lz-1]=1;
     243             :   else
     244   926461186 :     lz--;
     245             : 
     246   938746529 :   zd[1] = evalsigne(1) | evallgefint(lz);
     247   938746529 :   return zd;
     248             : }
     249             : 
     250             : /* assume x >= y */
     251             : INLINE GEN
     252  1032478075 : subiuspec(GEN x, ulong s, long nx)
     253             : {
     254             :   GEN zd;
     255             :   long lz;
     256             : 
     257  1032478075 :   if (nx == 1) return utoi(x[0] - s);
     258             : 
     259    98126401 :   lz = nx + 2; zd = cgeti(lz);
     260    98070042 :   mpn_sub_1 (LIMBS(zd), (mp_limb_t *)x, nx, s);
     261    98070123 :   if (! zd[lz - 1]) { --lz; }
     262             : 
     263    98070123 :   zd[1] = evalsigne(1) | evallgefint(lz);
     264    98070123 :   return zd;
     265             : }
     266             : 
     267             : /* assume x > y */
     268             : INLINE GEN
     269  2006696575 : subiispec(GEN x, GEN y, long nx, long ny)
     270             : {
     271             :   GEN zd;
     272             :   long lz;
     273  2006696575 :   if (ny==1) return subiuspec(x,*y,nx);
     274   984850391 :   lz = nx+2; zd = cgeti(lz);
     275             : 
     276   974897748 :   mpn_sub (LIMBS(zd), (mp_limb_t *)x, nx, (mp_limb_t *)y, ny);
     277   975837049 :   while (lz >= 3 && zd[lz - 1] == 0) { lz--; }
     278             : 
     279   975837049 :   zd[1] = evalsigne(1) | evallgefint(lz);
     280   975837049 :   return zd;
     281             : }
     282             : 
     283             : static void
     284   354808014 : roundr_up_ip(GEN x, long l)
     285             : {
     286   354808014 :   long i = l;
     287             :   for(;;)
     288             :   {
     289   359113590 :     if (++((ulong*)x)[--i]) break;
     290     2334418 :     if (i == 2) { x[2] = HIGHBIT; shiftr_inplace(x, 1); break; }
     291             :   }
     292   354808014 : }
     293             : 
     294             : void
     295   204369592 : affir(GEN x, GEN y)
     296             : {
     297   204369592 :   const long s = signe(x), ly = lg(y);
     298             :   long lx, sh, i;
     299             : 
     300   204369592 :   if (!s)
     301             :   {
     302     5215187 :     y[1] = evalexpo(-bit_accuracy(ly));
     303     5214977 :     return;
     304             :   }
     305   199154405 :   lx = lgefint(x); sh = bfffo(*int_MSW(x));
     306   199154405 :   y[1] = evalsigne(s) | evalexpo(bit_accuracy(lx)-sh-1);
     307   199154347 :   if (sh) {
     308   195963450 :     if (lx <= ly)
     309             :     {
     310   107164028 :       for (i=lx; i<ly; i++) y[i]=0;
     311   107164028 :       mpn_lshift(LIMBS(y),LIMBS(x),lx-2,sh);
     312   107164061 :       xmpn_mirror(LIMBS(y),lx-2);
     313   107164061 :       return;
     314             :     }
     315    88799422 :     mpn_lshift(LIMBS(y),LIMBS(x)+lx-ly,ly-2,sh);
     316    88799425 :     uel(y,2) |= uel(x,lx-ly+1) >> (BITS_IN_LONG-sh);
     317    88799425 :     xmpn_mirror(LIMBS(y),ly-2);
     318             :     /* lx > ly: round properly */
     319    88799436 :     if ((uel(x,lx-ly+1)<<sh) & HIGHBIT) roundr_up_ip(y, ly);
     320             :   }
     321             :   else {
     322     3190897 :     GEN xd=int_MSW(x);
     323     3190897 :     if (lx <= ly)
     324             :     {
     325     1650311 :       for (i=2; i<lx; i++,xd=int_precW(xd)) y[i]=*xd;
     326     1650311 :       for (   ; i<ly; i++) y[i]=0;
     327     1650311 :       return;
     328             :     }
     329     1540586 :     for (i=2; i<ly; i++,xd=int_precW(xd)) y[i]=*xd;
     330             :     /* lx > ly: round properly */
     331     1540586 :     if (uel(x,lx-ly+1) & HIGHBIT) roundr_up_ip(y, ly);
     332             :   }
     333             : }
     334             : 
     335             : INLINE GEN
     336   468361733 : shiftispec(GEN x, long nx, long n)
     337             : {
     338             :   long ny,m;
     339             :   GEN yd, y;
     340             : 
     341   468361733 :   if (!n) return icopyspec(x, nx);
     342             : 
     343   454102838 :   if (n > 0)
     344             :   {
     345   292884134 :     long d = dvmdsBIL(n, &m);
     346             :     long i;
     347             : 
     348   292810678 :     ny = nx + d + (m!=0);
     349   292810678 :     y = cgeti(ny + 2); yd = y + 2;
     350   292207899 :     for (i=0; i<d; i++) yd[i] = 0;
     351             : 
     352   292207899 :     if (!m) xmpn_copy((mp_limb_t *) (yd + d), (mp_limb_t *) x, nx);
     353             :     else
     354             :     {
     355   290860487 :       ulong carryd = mpn_lshift((mp_limb_t *) (yd + d), (mp_limb_t *) x, nx, m);
     356   291210734 :       if (carryd) yd[ny - 1] = carryd;
     357   272330959 :       else ny--;
     358             :     }
     359             :   }
     360             :   else
     361             :   {
     362   161218704 :     long d = dvmdsBIL(-n, &m);
     363             : 
     364   161217547 :     ny = nx - d;
     365   161217547 :     if (ny < 1) return gen_0;
     366   160440268 :     y = cgeti(ny + 2); yd = y + 2;
     367             : 
     368   160433355 :     if (!m) xmpn_copy((mp_limb_t *) yd, (mp_limb_t *) (x + d), nx - d);
     369             :     else
     370             :     {
     371   159357614 :       mpn_rshift((mp_limb_t *) yd, (mp_limb_t *) (x + d), nx - d, m);
     372   159357707 :       if (yd[ny - 1] == 0)
     373             :       {
     374    18265759 :         if (ny == 1) { set_avma((pari_sp)(yd + 1)); return gen_0; }
     375    16611694 :         ny--;
     376             :       }
     377             :     }
     378             :   }
     379   451216558 :   y[1] = evalsigne(1)|evallgefint(ny + 2);
     380   451216558 :   return y;
     381             : }
     382             : 
     383             : GEN
     384    41075238 : mantissa2nr(GEN x, long n)
     385             : {
     386    41075238 :   long ly, i, m, s = signe(x), lx = lg(x);
     387             :   GEN y;
     388    41075238 :   if (!s) return gen_0;
     389    41073936 :   if (!n)
     390             :   {
     391    20760047 :     y = cgeti(lx);
     392    20760047 :     y[1] = evalsigne(s) | evallgefint(lx);
     393    20760047 :     xmpn_mirrorcopy(LIMBS(y),RLIMBS(x),lx-2);
     394    20760047 :     return y;
     395             :   }
     396    20313889 :   if (n > 0)
     397             :   {
     398      458982 :     GEN z = (GEN)avma;
     399      458982 :     long d = dvmdsBIL(n, &m);
     400             : 
     401      458982 :     ly = lx+d; y = new_chunk(ly);
     402      458982 :     for ( ; d; d--) *--z = 0;
     403      458982 :     if (!m) for (i=2; i<lx; i++) y[i]=x[i];
     404             :     else
     405             :     {
     406      458426 :       register const ulong sh = BITS_IN_LONG - m;
     407      458426 :       shift_left(y,x, 2,lx-1, 0,m);
     408      458426 :       i = uel(x,2) >> sh;
     409             :       /* Extend y on the left? */
     410      458426 :       if (i) { ly++; y = new_chunk(1); y[2] = i; }
     411             :     }
     412             :   }
     413             :   else
     414             :   {
     415    19854907 :     ly = lx - dvmdsBIL(-n, &m);
     416    19854747 :     if (ly<3) return gen_0;
     417    19854747 :     y = new_chunk(ly);
     418    19854585 :     if (m) {
     419    19797514 :       shift_right(y,x, 2,ly, 0,m);
     420    19797910 :       if (y[2] == 0)
     421             :       {
     422           0 :         if (ly==3) { set_avma((pari_sp)(y+3)); return gen_0; }
     423           0 :         ly--; set_avma((pari_sp)(++y));
     424             :       }
     425             :     } else {
     426       57071 :       for (i=2; i<ly; i++) y[i]=x[i];
     427             :     }
     428             :   }
     429    20313963 :   xmpn_mirror(LIMBS(y),ly-2);
     430    20314052 :   y[1] = evalsigne(s)|evallgefint(ly);
     431    20314052 :   y[0] = evaltyp(t_INT)|evallg(ly); return y;
     432             : }
     433             : 
     434             : GEN
     435     1908056 : truncr(GEN x)
     436             : {
     437             :   long s, e, d, m, i;
     438             :   GEN y;
     439     1908056 :   if ((s=signe(x)) == 0 || (e=expo(x)) < 0) return gen_0;
     440     1515658 :   d = nbits2lg(e+1); m = remsBIL(e);
     441     1515658 :   if (d > lg(x)) pari_err_PREC( "truncr (precision loss in truncation)");
     442             : 
     443     1515658 :   y=cgeti(d); y[1] = evalsigne(s) | evallgefint(d);
     444     1515658 :   if (++m == BITS_IN_LONG)
     445         275 :     for (i=2; i<d; i++) y[d-i+1]=x[i];
     446             :   else
     447             :   {
     448     1515383 :     GEN z=cgeti(d);
     449     1515383 :     for (i=2; i<d; i++) z[d-i+1]=x[i];
     450     1515383 :     mpn_rshift(LIMBS(y),LIMBS(z),d-2,BITS_IN_LONG-m);
     451     1515383 :     set_avma((pari_sp)y);
     452             :   }
     453     1515658 :   return y;
     454             : }
     455             : 
     456             : /* integral part */
     457             : GEN
     458      767944 : floorr(GEN x)
     459             : {
     460             :   long e, d, m, i, lx;
     461             :   GEN y;
     462      767944 :   if (signe(x) >= 0) return truncr(x);
     463      257205 :   if ((e=expo(x)) < 0) return gen_m1;
     464       92367 :   d = nbits2lg(e+1); m = remsBIL(e);
     465       92367 :   lx=lg(x); if (d>lx) pari_err_PREC( "floorr (precision loss in truncation)");
     466       92367 :   y = cgeti(d+1);
     467       92367 :   if (++m == BITS_IN_LONG)
     468             :   {
     469         569 :     for (i=2; i<d; i++) y[d-i+1]=x[i];
     470         569 :     i=d; while (i<lx && !x[i]) i++;
     471         569 :     if (i==lx) goto END;
     472             :   }
     473             :   else
     474             :   {
     475       91798 :     GEN z=cgeti(d);
     476       91798 :     for (i=2; i<d; i++) z[d-i+1]=x[i];
     477       91798 :     mpn_rshift(LIMBS(y),LIMBS(z),d-2,BITS_IN_LONG-m);
     478       91798 :     if (uel(x,d-1)<<m == 0)
     479             :     {
     480       25865 :       i=d; while (i<lx && !x[i]) i++;
     481       25865 :       if (i==lx) goto END;
     482             :     }
     483             :   }
     484       80495 :   if (mpn_add_1(LIMBS(y),LIMBS(y),d-2,1))
     485           0 :     y[d++]=1;
     486             : END:
     487       92367 :   y[1] = evalsigne(-1) | evallgefint(d);
     488       92367 :   return y;
     489             : }
     490             : 
     491             : INLINE int
     492  2452492950 : cmpiispec(GEN x, GEN y, long lx, long ly)
     493             : {
     494  2452492950 :   if (lx < ly) return -1;
     495  2332843699 :   if (lx > ly) return  1;
     496  2224528315 :   return mpn_cmp((mp_limb_t*)x,(mp_limb_t*)y, lx);
     497             : }
     498             : 
     499             : INLINE int
     500   176777153 : equaliispec(GEN x, GEN y, long lx, long ly)
     501             : {
     502   176777153 :   if (lx != ly) return 0;
     503   176760437 :   return !mpn_cmp((mp_limb_t*)x,(mp_limb_t*)y, lx);
     504             : }
     505             : 
     506             : /***********************************************************************/
     507             : /**                                                                   **/
     508             : /**                          MULTIPLICATION                           **/
     509             : /**                                                                   **/
     510             : /***********************************************************************/
     511             : /* assume ny > 0 */
     512             : INLINE GEN
     513  2661417169 : muluispec(ulong x, GEN y, long ny)
     514             : {
     515  2661417169 :   if (ny == 1)
     516  2031396109 :     return muluu(x, *y);
     517             :   else
     518             :   {
     519   630021060 :     long lz = ny+3;
     520   630021060 :     GEN z = cgeti(lz);
     521   626932512 :     ulong hi = mpn_mul_1 (LIMBS(z), (mp_limb_t *)y, ny, x);
     522   627033171 :     if (hi) { z[lz - 1] = hi; } else lz--;
     523   627033171 :     z[1] = evalsigne(1) | evallgefint(lz);
     524   627033171 :     return z;
     525             :   }
     526             : }
     527             : 
     528             : /* a + b*|y| */
     529             : GEN
     530      700026 : addumului(ulong a, ulong b, GEN y)
     531             : {
     532             :   GEN z;
     533             :   long i, lz;
     534             :   ulong hi;
     535      700026 :   if (!b || !signe(y)) return utoi(a);
     536      699958 :   lz = lgefint(y)+1;
     537      699958 :   z = cgeti(lz);
     538      699958 :   z[2]=a;
     539      699958 :   for(i=3;i<lz;i++) z[i]=0;
     540      699958 :   hi=mpn_addmul_1(LIMBS(z), LIMBS(y), NLIMBS(y), b);
     541      699958 :   if (hi) z[lz-1]=hi; else lz--;
     542      699958 :   z[1] = evalsigne(1) | evallgefint(lz);
     543      699958 :   set_avma((pari_sp)z); return z;
     544             : }
     545             : 
     546             : /***********************************************************************/
     547             : /**                                                                   **/
     548             : /**                          DIVISION                                 **/
     549             : /**                                                                   **/
     550             : /***********************************************************************/
     551             : 
     552             : ulong
     553   935718160 : umodiu(GEN y, ulong x)
     554             : {
     555   935718160 :   long sy=signe(y);
     556             :   ulong hi;
     557   935718160 :   if (!x) pari_err_INV("umodiu",gen_0);
     558   937288528 :   if (!sy) return 0;
     559   707685271 :   hi = mpn_mod_1(LIMBS(y),NLIMBS(y),x);
     560   708192267 :   if (!hi) return 0;
     561   697728055 :   return (sy > 0)? hi: x - hi;
     562             : }
     563             : 
     564             : /* return |y| \/ x */
     565             : GEN
     566   155324521 : absdiviu_rem(GEN y, ulong x, ulong *rem)
     567             : {
     568             :   long ly;
     569             :   GEN z;
     570             : 
     571   155324521 :   if (!x) pari_err_INV("absdiviu_rem",gen_0);
     572   155325081 :   if (!signe(y)) { *rem = 0; return gen_0; }
     573             : 
     574   154736908 :   ly = lgefint(y);
     575   154736908 :   if (ly == 3 && (ulong)x > uel(y,2)) { *rem = uel(y,2); return gen_0; }
     576             : 
     577   154575554 :   z = cgeti(ly);
     578   154500294 :   *rem = mpn_divrem_1(LIMBS(z), 0, LIMBS(y), NLIMBS(y), x);
     579   154525887 :   if (z [ly - 1] == 0) ly--;
     580   154525887 :   z[1] = evallgefint(ly) | evalsigne(1);
     581   154525887 :   return z;
     582             : }
     583             : 
     584             : GEN
     585    59550043 : divis_rem(GEN y, long x, long *rem)
     586             : {
     587    59550043 :   long sy=signe(y),ly,s;
     588             :   GEN z;
     589             : 
     590    59550043 :   if (!x) pari_err_INV("divis_rem",gen_0);
     591    59550074 :   if (!sy) { *rem = 0; return gen_0; }
     592    47126104 :   if (x<0) { s = -sy; x = -x; } else s = sy;
     593             : 
     594    47126104 :   ly = lgefint(y);
     595    47126104 :   if (ly == 3 && (ulong)x > uel(y,2)) { *rem = itos(y); return gen_0; }
     596             : 
     597    33767983 :   z = cgeti(ly);
     598    33767743 :   *rem = mpn_divrem_1(LIMBS(z), 0, LIMBS(y), NLIMBS(y), x);
     599    33767760 :   if (sy<0) *rem = -  *rem;
     600    33767760 :   if (z[ly - 1] == 0) ly--;
     601    33767760 :   z[1] = evallgefint(ly) | evalsigne(s);
     602    33767760 :   return z;
     603             : }
     604             : 
     605             : GEN
     606     1189057 : divis(GEN y, long x)
     607             : {
     608     1189057 :   long sy=signe(y),ly,s;
     609             :   GEN z;
     610             : 
     611     1189057 :   if (!x) pari_err_INV("divis",gen_0);
     612     1189057 :   if (!sy) return gen_0;
     613     1189053 :   if (x<0) { s = -sy; x = -x; } else s=sy;
     614             : 
     615     1189053 :   ly = lgefint(y);
     616     1189053 :   if (ly == 3 && (ulong)x > uel(y,2)) return gen_0;
     617             : 
     618     1189053 :   z = cgeti(ly);
     619     1189053 :   (void)mpn_divrem_1(LIMBS(z), 0, LIMBS(y), NLIMBS(y), x);
     620     1189053 :   if (z[ly - 1] == 0) ly--;
     621     1189053 :   z[1] = evallgefint(ly) | evalsigne(s);
     622     1189053 :   return z;
     623             : }
     624             : 
     625             : /* We keep llx bits of x and lly bits of y*/
     626             : static GEN
     627    61658105 : divrr_with_gmp(GEN x, GEN y)
     628             : {
     629    61658105 :   long lx=RNLIMBS(x),ly=RNLIMBS(y);
     630    61658105 :   long lw=minss(lx,ly);
     631    61655602 :   long lly=minss(lw+1,ly);
     632    61652634 :   GEN  w=cgetr(lw+2);
     633    61636066 :   long lu=lw+lly;
     634    61636066 :   long llx=minss(lu,lx);
     635    61627858 :   mp_limb_t *u=(mp_limb_t *)new_chunk(lu);
     636    61593150 :   mp_limb_t *z=(mp_limb_t *)new_chunk(lly);
     637             :   mp_limb_t *q,*r;
     638    61556601 :   long e=expo(x)-expo(y);
     639    61556601 :   long sx=signe(x),sy=signe(y);
     640    61556601 :   xmpn_mirrorcopy(z,RLIMBS(y),lly);
     641    61591671 :   xmpn_mirrorcopy(u+lu-llx,RLIMBS(x),llx);
     642    61611878 :   xmpn_zero(u,lu-llx);
     643    61662584 :   q = (mp_limb_t *)new_chunk(lw+1);
     644    61659419 :   r = (mp_limb_t *)new_chunk(lly);
     645             : 
     646    61615813 :   mpn_tdiv_qr(q,r,0,u,lu,z,lly);
     647             : 
     648             :   /*Round up: This is not exactly correct we should test 2*r>z*/
     649    61729344 :   if (uel(r,lly-1) > (uel(z,lly-1)>>1))
     650    30589820 :     mpn_add_1(q,q,lw+1,1);
     651             : 
     652    61729392 :   xmpn_mirrorcopy(RLIMBS(w),q,lw);
     653             : 
     654    61729116 :   if (q[lw] == 0) e--;
     655    26264869 :   else if (q[lw] == 1) { shift_right(w,w, 2,lw+2, 1,1); }
     656          14 :   else { w[2] = HIGHBIT; e++; }
     657    61716199 :   if (sy < 0) sx = -sx;
     658    61716199 :   w[1] = evalsigne(sx) | evalexpo(e);
     659    61701187 :   set_avma((pari_sp)w); return w;
     660             : }
     661             : 
     662             : /* We keep llx bits of x and lly bits of y*/
     663             : static GEN
     664     5759800 : divri_with_gmp(GEN x, GEN y)
     665             : {
     666     5759800 :   long llx=RNLIMBS(x),ly=NLIMBS(y);
     667     5759800 :   long lly=minss(llx+1,ly);
     668     5759794 :   GEN  w=cgetr(llx+2);
     669     5759730 :   long lu=llx+lly, ld=ly-lly;
     670     5759730 :   mp_limb_t *u=(mp_limb_t *)new_chunk(lu);
     671     5759690 :   mp_limb_t *z=(mp_limb_t *)new_chunk(lly);
     672             :   mp_limb_t *q,*r;
     673     5759635 :   long sh=bfffo(y[ly+1]);
     674     5759635 :   long e=expo(x)-expi(y);
     675     5759696 :   long sx=signe(x),sy=signe(y);
     676     5759696 :   if (sh) mpn_lshift(z,LIMBS(y)+ld,lly,sh);
     677      154552 :   else xmpn_copy(z,LIMBS(y)+ld,lly);
     678     5759719 :   xmpn_mirrorcopy(u+lu-llx,RLIMBS(x),llx);
     679     5759791 :   xmpn_zero(u,lu-llx);
     680     5759807 :   q = (mp_limb_t *)new_chunk(llx+1);
     681     5759774 :   r = (mp_limb_t *)new_chunk(lly);
     682             : 
     683     5759718 :   mpn_tdiv_qr(q,r,0,u,lu,z,lly);
     684             : 
     685             :   /*Round up: This is not exactly correct we should test 2*r>z*/
     686     5759824 :   if (uel(r,lly-1) > (uel(z,lly-1)>>1))
     687     2325616 :     mpn_add_1(q,q,llx+1,1);
     688             : 
     689     5759825 :   xmpn_mirrorcopy(RLIMBS(w),q,llx);
     690             : 
     691     5759826 :   if (q[llx] == 0) e--;
     692     3723108 :   else if (q[llx] == 1) { shift_right(w,w, 2,llx+2, 1,1); }
     693           0 :   else { w[2] = HIGHBIT; e++; }
     694     5759807 :   if (sy < 0) sx = -sx;
     695     5759807 :   w[1] = evalsigne(sx) | evalexpo(e);
     696     5759799 :   set_avma((pari_sp)w); return w;
     697             : }
     698             : 
     699             : GEN
     700    30865247 : divri(GEN x, GEN y)
     701             : {
     702    30865247 :   long  s = signe(y);
     703             : 
     704    30865247 :   if (!s) pari_err_INV("divri",gen_0);
     705    30865252 :   if (!signe(x)) return real_0_bit(expo(x) - expi(y));
     706    30725700 :   if (!is_bigint(y)) {
     707    24965934 :     GEN z = divru(x, y[2]);
     708    24965815 :     if (s < 0) togglesign(z);
     709    24965815 :     return z;
     710             :   }
     711     5759794 :   return divri_with_gmp(x,y);
     712             : }
     713             : 
     714             : GEN
     715   140507521 : divrr(GEN x, GEN y)
     716             : {
     717   140507521 :   long sx=signe(x), sy=signe(y), lx,ly,lr,e,i,j;
     718             :   ulong y0,y1;
     719             :   GEN r, r1;
     720             : 
     721   140507521 :   if (!sy) pari_err_INV("divrr",y);
     722   140523845 :   e = expo(x) - expo(y);
     723   140523845 :   if (!sx) return real_0_bit(e);
     724   138794062 :   if (sy<0) sx = -sx;
     725             : 
     726   138794062 :   lx=lg(x); ly=lg(y);
     727   138794062 :   if (ly==3)
     728             :   {
     729    77135875 :     ulong k = x[2], l = (lx>3)? x[3]: 0;
     730             :     LOCAL_HIREMAINDER;
     731    77135875 :     if (k < uel(y,2)) e--;
     732             :     else
     733             :     {
     734    35880390 :       l >>= 1; if (k&1) l |= HIGHBIT;
     735    35880390 :       k >>= 1;
     736             :     }
     737    77135875 :     hiremainder = k; k = divll(l,y[2]);
     738    77135875 :     if (hiremainder > (uel(y,2) >> 1) && !++k) { k = HIGHBIT; e++; }
     739    77135875 :     r = cgetr(3);
     740    77135660 :     r[1] = evalsigne(sx) | evalexpo(e);
     741    77135615 :     r[2] = k; return r;
     742             :   }
     743             : 
     744    61658187 :   if (ly>=DIVRR_GMP_LIMIT)
     745    61658187 :     return divrr_with_gmp(x,y);
     746             : 
     747           0 :   lr = minss(lx,ly); r = new_chunk(lr);
     748           0 :   r1 = r-1;
     749           0 :   r1[1] = 0; for (i=2; i<lr; i++) r1[i]=x[i];
     750           0 :   r1[lr] = (lx>ly)? x[lr]: 0;
     751           0 :   y0 = y[2]; y1 = y[3];
     752           0 :   for (i=0; i<lr-1; i++)
     753             :   { /* r1 = r + (i-1), OK up to r1[2] (accesses at most r[lr]) */
     754             :     ulong k, qp;
     755             :     LOCAL_HIREMAINDER;
     756             :     LOCAL_OVERFLOW;
     757             : 
     758           0 :     if (uel(r1,1) == y0) { qp = ULONG_MAX; k = addll(y0,r1[2]); }
     759             :     else
     760             :     {
     761           0 :       if (uel(r1,1) > y0) /* can't happen if i=0 */
     762             :       {
     763           0 :         GEN y1 = y+1;
     764           0 :         j = lr-i; r1[j] = subll(r1[j],y1[j]);
     765           0 :         for (j--; j>0; j--) r1[j] = subllx(r1[j],y1[j]);
     766           0 :         j=i; do uel(r,--j)++; while (j && !r[j]);
     767             :       }
     768           0 :       hiremainder = r1[1]; overflow = 0;
     769           0 :       qp = divll(r1[2],y0); k = hiremainder;
     770             :     }
     771           0 :     j = lr-i+1;
     772           0 :     if (!overflow)
     773             :     {
     774             :       long k3, k4;
     775           0 :       k3 = mulll(qp,y1);
     776           0 :       if (j == 3) /* i = lr - 2 maximal, r1[3] undefined -> 0 */
     777           0 :         k4 = subll(hiremainder,k);
     778             :       else
     779             :       {
     780           0 :         k3 = subll(k3, r1[3]);
     781           0 :         k4 = subllx(hiremainder,k);
     782             :       }
     783           0 :       while (!overflow && k4) { qp--; k3=subll(k3,y1); k4=subllx(k4,y0); }
     784             :     }
     785           0 :     if (j<ly) (void)mulll(qp,y[j]); else { hiremainder = 0 ; j = ly; }
     786           0 :     for (j--; j>1; j--)
     787             :     {
     788           0 :       r1[j] = subll(r1[j], addmul(qp,y[j]));
     789           0 :       hiremainder += overflow;
     790             :     }
     791           0 :     if (uel(r1,1) != hiremainder)
     792             :     {
     793           0 :       if (uel(r1,1) < hiremainder)
     794             :       {
     795           0 :         qp--;
     796           0 :         j = lr-i-(lr-i>=ly); r1[j] = addll(r1[j], y[j]);
     797           0 :         for (j--; j>1; j--) r1[j] = addllx(r1[j], y[j]);
     798             :       }
     799             :       else
     800             :       {
     801           0 :         uel(r1,1) -= hiremainder;
     802           0 :         while (r1[1])
     803             :         {
     804           0 :           qp++; if (!qp) { j=i; do uel(r,--j)++; while (j && !r[j]); }
     805           0 :           j = lr-i-(lr-i>=ly); r1[j] = subll(r1[j],y[j]);
     806           0 :           for (j--; j>1; j--) r1[j] = subllx(r1[j],y[j]);
     807           0 :           uel(r1,1) -= overflow;
     808             :         }
     809             :       }
     810             :     }
     811           0 :     *++r1 = qp;
     812             :   }
     813             :   /* i = lr-1 */
     814             :   /* round correctly */
     815           0 :   if (uel(r1,1) > (y0>>1))
     816             :   {
     817           0 :     j=i; do uel(r,--j)++; while (j && !r[j]);
     818             :   }
     819           0 :   r1 = r-1; for (j=i; j>=2; j--) r[j]=r1[j];
     820           0 :   if (r[0] == 0) e--;
     821           0 :   else if (r[0] == 1) { shift_right(r,r, 2,lr, 1,1); }
     822             :   else { /* possible only when rounding up to 0x2 0x0 ... */
     823           0 :     r[2] = (long)HIGHBIT; e++;
     824             :   }
     825           0 :   r[0] = evaltyp(t_REAL)|evallg(lr);
     826           0 :   r[1] = evalsigne(sx) | evalexpo(e);
     827           0 :   return r;
     828             : }
     829             : 
     830             : /* Integer division x / y: such that sign(r) = sign(x)
     831             :  *   if z = ONLY_REM return remainder, otherwise return quotient
     832             :  *   if z != NULL set *z to remainder
     833             :  *   *z is the last object on stack (and thus can be disposed of with cgiv
     834             :  *   instead of gerepile)
     835             :  * If *z is zero, we put gen_0 here and no copy.
     836             :  * space needed: lx + ly
     837             :  */
     838             : GEN
     839   990577348 : dvmdii(GEN x, GEN y, GEN *z)
     840             : {
     841   990577348 :   long sx=signe(x),sy=signe(y);
     842             :   long lx, ly, lq;
     843             :   pari_sp av;
     844             :   GEN r,q;
     845             : 
     846   990577348 :   if (!sy) pari_err_INV("dvmdii",y);
     847   995666990 :   if (!sx)
     848             :   {
     849    29418406 :     if (!z || z == ONLY_REM) return gen_0;
     850    14374986 :     *z=gen_0; return gen_0;
     851             :   }
     852   966248584 :   lx=lgefint(x);
     853   966248584 :   ly=lgefint(y); lq=lx-ly;
     854   966248584 :   if (lq <= 0)
     855             :   {
     856   596502907 :     if (lq == 0)
     857             :     {
     858   563398283 :       long s=mpn_cmp(LIMBS(x),LIMBS(y),NLIMBS(x));
     859   565783351 :       if (s>0) goto DIVIDE;
     860   223811527 :       if (s==0)
     861             :       {
     862    12266985 :         if (z == ONLY_REM) return gen_0;
     863     6434300 :         if (z) *z = gen_0;
     864     6434300 :         if (sx < 0) sy = -sy;
     865     6434300 :         return stoi(sy);
     866             :       }
     867             :     }
     868   244649166 :     if (z == ONLY_REM) return icopy(x);
     869    11198218 :     if (z) *z = icopy(x);
     870    11196323 :     return gen_0;
     871             :   }
     872             : DIVIDE: /* quotient is non-zero */
     873   711717501 :   av=avma; if (sx<0) sy = -sy;
     874   711717501 :   if (ly==3)
     875             :   {
     876   306753751 :     ulong lq = lx;
     877             :     ulong si;
     878   306753751 :     q  = cgeti(lq);
     879   306349020 :     si = mpn_divrem_1(LIMBS(q), 0, LIMBS(x), NLIMBS(x), y[2]);
     880   306696916 :     if (q[lq - 1] == 0) lq--;
     881   306696916 :     if (z == ONLY_REM)
     882             :     {
     883   251702102 :       set_avma(av); if (!si) return gen_0;
     884   224280860 :       r=cgeti(3);
     885   223920125 :       r[1] = evalsigne(sx) | evallgefint(3);
     886   223920125 :       r[2] = si; return r;
     887             :     }
     888    54994814 :     q[1] = evalsigne(sy) | evallgefint(lq);
     889    54994814 :     if (!z) return q;
     890    53305347 :     if (!si) { *z=gen_0; return q; }
     891    23238469 :     r=cgeti(3);
     892    23238400 :     r[1] = evalsigne(sx) | evallgefint(3);
     893    23238400 :     r[2] = si; *z=r; return q;
     894             :   }
     895   404963750 :   if (z == ONLY_REM)
     896             :   {
     897   395100863 :     ulong lr = lgefint(y);
     898   395100863 :     ulong lq = lgefint(x)-lgefint(y)+3;
     899   395100863 :     GEN r = cgeti(lr);
     900   385284365 :     GEN q = cgeti(lq);
     901   381141054 :     mpn_tdiv_qr(LIMBS(q), LIMBS(r),0, LIMBS(x), NLIMBS(x), LIMBS(y), NLIMBS(y));
     902   404271146 :     if (!r[lr - 1])
     903             :     {
     904    94027858 :       while(lr>2 && !r[lr - 1]) lr--;
     905    94027858 :       if (lr == 2) {set_avma(av); return gen_0;} /* exact division */
     906             :     }
     907   395519154 :     r[1] = evalsigne(sx) | evallgefint(lr);
     908   395519154 :     set_avma((pari_sp)r); return r;
     909             :   }
     910             :   else
     911             :   {
     912     9862887 :     ulong lq = lgefint(x)-lgefint(y)+3;
     913     9862887 :     ulong lr = lgefint(y);
     914     9862887 :     GEN q = cgeti(lq);
     915     9862546 :     GEN r = cgeti(lr);
     916     9862315 :     mpn_tdiv_qr(LIMBS(q), LIMBS(r),0, LIMBS(x), NLIMBS(x), LIMBS(y), NLIMBS(y));
     917     9862919 :     if (q[lq - 1] == 0) lq--;
     918     9862919 :     q[1] = evalsigne(sy) | evallgefint(lq);
     919     9862919 :     if (!z) { set_avma((pari_sp)q); return q; }
     920     7859725 :     if (!r[lr - 1])
     921             :     {
     922     2778221 :       while(lr>2 && !r[lr - 1]) lr--;
     923     2778221 :       if (lr == 2) { set_avma((pari_sp)q); *z=gen_0; return q; } /* exact division */
     924             :     }
     925     5679287 :     r[1] = evalsigne(sx) | evallgefint(lr);
     926     5679287 :     set_avma((pari_sp)r); *z = r; return q;
     927             :   }
     928             : }
     929             : 
     930             : /* Montgomery reduction.
     931             :  * N has k words, assume T >= 0 has less than 2k.
     932             :  * Return res := T / B^k mod N, where B = 2^BIL
     933             :  * such that 0 <= res < T/B^k + N  and  res has less than k words */
     934             : GEN
     935    14302307 : red_montgomery(GEN T, GEN N, ulong inv)
     936             : {
     937             :   pari_sp av;
     938             :   GEN Te, Td, Ne, Nd, scratch;
     939    14302307 :   ulong i, j, m, t, d, k = NLIMBS(N);
     940             :   int carry;
     941             :   LOCAL_HIREMAINDER;
     942             :   LOCAL_OVERFLOW;
     943             : 
     944    14302307 :   if (k == 0) return gen_0;
     945    14302307 :   d = NLIMBS(T); /* <= 2*k */
     946    14302307 :   if (d == 0) return gen_0;
     947             : #ifdef DEBUG
     948             :   if (d > 2*k) pari_err_BUG("red_montgomery");
     949             : #endif
     950    14302298 :   if (k == 1)
     951             :   { /* as below, special cased for efficiency */
     952         492 :     ulong n = uel(N,2);
     953         492 :     if (d == 1) {
     954         492 :       hiremainder = uel(T,2);
     955         492 :       m = hiremainder * inv;
     956         492 :       (void)addmul(m, n); /* t + m*n = 0 */
     957         492 :       return utoi(hiremainder);
     958             :     } else { /* d = 2 */
     959           0 :       hiremainder = uel(T,2);
     960           0 :       m = hiremainder * inv;
     961           0 :       (void)addmul(m, n); /* t + m*n = 0 */
     962           0 :       t = addll(hiremainder, uel(T,3));
     963           0 :       if (overflow) t -= n; /* t > n doesn't fit in 1 word */
     964           0 :       return utoi(t);
     965             :     }
     966             :   }
     967             :   /* assume k >= 2 */
     968    14301806 :   av = avma; scratch = new_chunk(k<<1); /* >= k + 2: result fits */
     969             : 
     970             :   /* copy T to scratch space (pad with zeroes to 2k words) */
     971    13956269 :   Td = scratch;
     972    13956269 :   Te = T + 2;
     973    13956269 :   for (i=0; i < d     ; i++) *Td++ = *Te++;
     974    13956269 :   for (   ; i < (k<<1); i++) *Td++ = 0;
     975             : 
     976    13956269 :   Te = scratch - 1; /* 1 beyond end of current T mantissa (in scratch) */
     977    13956269 :   Ne = N + 1;       /* 1 beyond end of N mantissa */
     978             : 
     979    13956269 :   carry = 0;
     980   100453285 :   for (i=0; i<k; i++) /* set T := T/B nod N, k times */
     981             :   {
     982    86497016 :     Td = Te; /* one beyond end of (new) T mantissa */
     983    86497016 :     Nd = Ne;
     984    86497016 :     hiremainder = *++Td;
     985    86497016 :     m = hiremainder * inv; /* solve T + m N = O(B) */
     986             : 
     987             :     /* set T := (T + mN) / B */
     988    86497016 :     Te = Td;
     989    86497016 :     (void)addmul(m, *++Nd); /* = 0 */
     990   776978447 :     for (j=1; j<k; j++)
     991             :     {
     992   690481431 :       t = addll(addmul(m, *++Nd), *++Td);
     993   690481431 :       *Td = t;
     994   690481431 :       hiremainder += overflow;
     995             :     }
     996    86497016 :     t = addll(hiremainder, *++Td); *Td = t + carry;
     997    86497016 :     carry = (overflow || (carry && *Td == 0));
     998             :   }
     999    13956269 :   if (carry)
    1000             :   { /* Td > N overflows (k+1 words), set Td := Td - N */
    1001       35455 :     GEN NE = N + k+1;
    1002       35455 :     Td = Te;
    1003       35455 :     Nd = Ne;
    1004       35455 :     t = subll(*++Td, *++Nd); *Td = t;
    1005       35455 :     while (Nd < NE) { t = subllx(*++Td, *++Nd); *Td = t; }
    1006             :   }
    1007             : 
    1008             :   /* copy result */
    1009    13956269 :   Td = (GEN)av - 1; /* *Td = high word of final result */
    1010    13956269 :   while (*Td == 0 && Te < Td) Td--; /* strip leading 0s */
    1011    13956269 :   k = Td - Te; if (!k) { set_avma(av); return gen_0; }
    1012    13956269 :   Td = (GEN)av - k; /* will write mantissa there */
    1013    13956269 :   (void)memmove(Td, Te+1, k*sizeof(long));
    1014    13956269 :   Td -= 2;
    1015    13956269 :   Td[0] = evaltyp(t_INT) | evallg(k+2);
    1016    13931036 :   Td[1] = evalsigne(1) | evallgefint(k+2);
    1017             : #ifdef DEBUG
    1018             : {
    1019             :   long l = NLIMBS(N), s = BITS_IN_LONG*l;
    1020             :   GEN R = int2n(s);
    1021             :   GEN res = remii(mulii(T, Fp_inv(R, N)), N);
    1022             :   if (k > lgefint(N)
    1023             :     || !equalii(remii(Td,N),res)
    1024             :     || cmpii(Td, addii(shifti(T, -s), N)) >= 0) pari_err_BUG("red_montgomery");
    1025             : }
    1026             : #endif
    1027    13931036 :   set_avma((pari_sp)Td); return Td;
    1028             : }
    1029             : 
    1030             : /* EXACT INTEGER DIVISION */
    1031             : 
    1032             : /* use undocumented GMP interface */
    1033             : static void
    1034    58587219 : GEN2mpz(mpz_t X, GEN x)
    1035             : {
    1036    58587219 :   long l = lgefint(x)-2;
    1037    58587219 :   X->_mp_alloc = l;
    1038    58587219 :   X->_mp_size = signe(x) > 0? l: -l;
    1039    58587219 :   X->_mp_d = LIMBS(x);
    1040    58587219 : }
    1041             : static void
    1042    29293618 : mpz2GEN(GEN z, mpz_t Z)
    1043             : {
    1044    29293618 :   long l = Z->_mp_size;
    1045    29293618 :   z[1] = evalsigne(l > 0? 1: -1) | evallgefint(labs(l)+2);
    1046    29293618 : }
    1047             : 
    1048             : #ifdef mpn_divexact_1
    1049             : static GEN
    1050   293485792 : diviuexact_i(GEN x, ulong y)
    1051             : {
    1052   293485792 :   long l = lgefint(x);
    1053   293485792 :   GEN z = cgeti(l);
    1054   293463050 :   mpn_divexact_1(LIMBS(z), LIMBS(x), NLIMBS(x), y);
    1055   293463458 :   if (z[l-1] == 0) l--;
    1056   293463458 :   z[1] = evallgefint(l) | evalsigne(signe(x));
    1057   293463458 :   return z;
    1058             : }
    1059             : #elif 1 && !defined(_WIN64) /* mpz_divexact_ui is not LLP64 friendly   */
    1060             :                             /* assume y != 0 and the division is exact */
    1061             : static GEN
    1062             : diviuexact_i(GEN x, ulong y)
    1063             : {
    1064             :   long l = lgefint(x);
    1065             :   GEN z = cgeti(l);
    1066             :   mpz_t X, Z;
    1067             :   GEN2mpz(X, x);
    1068             :   Z->_mp_alloc = l-2;
    1069             :   Z->_mp_size  = l-2;
    1070             :   Z->_mp_d = LIMBS(z);
    1071             :   mpz_divexact_ui(Z, X, y);
    1072             :   mpz2GEN(z, Z); return z;
    1073             : }
    1074             : #else
    1075             : /* assume y != 0 and the division is exact */
    1076             : static GEN
    1077             : diviuexact_i(GEN x, ulong y)
    1078             : {
    1079             :   /*TODO: implement true exact division.*/
    1080             :   return divii(x,utoi(y));
    1081             : }
    1082             : #endif
    1083             : 
    1084             : GEN
    1085    40256438 : diviuexact(GEN x, ulong y)
    1086             : {
    1087             :   GEN z;
    1088    40256438 :   if (!signe(x)) return gen_0;
    1089    39983928 :   z = diviuexact_i(x, y);
    1090    39981177 :   if (lgefint(z) == 2) pari_err_OP("exact division", x, utoi(y));
    1091    39981061 :   return z;
    1092             : }
    1093             : 
    1094             : /* Find z such that x=y*z, knowing that y | x (unchecked) */
    1095             : GEN
    1096   337028752 : diviiexact(GEN x, GEN y)
    1097             : {
    1098             :   GEN z;
    1099   337028752 :   if (!signe(y)) pari_err_INV("diviiexact",y);
    1100   337049394 :   if (!signe(x)) return gen_0;
    1101   282786277 :   if (lgefint(y) == 3)
    1102             :   {
    1103   253492663 :     z = diviuexact_i(x, y[2]);
    1104   253482912 :     if (signe(y) < 0) togglesign(z);
    1105             :   }
    1106             :   else
    1107             :   {
    1108    29293614 :     long l = lgefint(x);
    1109             :     mpz_t X, Y, Z;
    1110    29293614 :     z = cgeti(l);
    1111    29293611 :     GEN2mpz(X, x);
    1112    29293611 :     GEN2mpz(Y, y);
    1113    29293611 :     Z->_mp_alloc = l-2;
    1114    29293611 :     Z->_mp_size  = l-2;
    1115    29293611 :     Z->_mp_d = LIMBS(z);
    1116    29293611 :     mpz_divexact(Z, X, Y);
    1117    29293619 :     mpz2GEN(z, Z);
    1118             :   }
    1119   282776955 :   if (lgefint(z) == 2) pari_err_OP("exact division", x, y);
    1120   282776419 :   return z;
    1121             : }
    1122             : 
    1123             : /* assume yz != and yz | x */
    1124             : GEN
    1125      170930 : diviuuexact(GEN x, ulong y, ulong z)
    1126             : {
    1127             :   long tmp[4];
    1128             :   ulong t;
    1129             :   LOCAL_HIREMAINDER;
    1130      170930 :   t = mulll(y, z);
    1131      170930 :   if (!hiremainder) return diviuexact(x, t);
    1132           0 :   tmp[0] = evaltyp(t_INT)|_evallg(4);
    1133           0 :   tmp[1] = evalsigne(1)|evallgefint(4);
    1134           0 :   tmp[2] = t;
    1135           0 :   tmp[3] = hiremainder;
    1136           0 :   return diviiexact(x, tmp);
    1137             : }
    1138             : 
    1139             : 
    1140             : /********************************************************************/
    1141             : /**                                                                **/
    1142             : /**               INTEGER MULTIPLICATION                           **/
    1143             : /**                                                                **/
    1144             : /********************************************************************/
    1145             : 
    1146             : /* nx >= ny = num. of digits of x, y (not GEN, see mulii) */
    1147             : GEN
    1148  2751964277 : muliispec(GEN x, GEN y, long nx, long ny)
    1149             : {
    1150             :   GEN zd;
    1151             :   long lz;
    1152             :   ulong hi;
    1153             : 
    1154  2751964277 :   if (nx < ny) swapspec(x,y, nx,ny);
    1155  2751964277 :   if (!ny) return gen_0;
    1156  2751964277 :   if (ny == 1) return muluispec((ulong)*y, x, nx);
    1157             : 
    1158   663539596 :   lz = nx+ny+2;
    1159   663539596 :   zd = cgeti(lz);
    1160   659001825 :   hi = mpn_mul(LIMBS(zd), (mp_limb_t *)x, nx, (mp_limb_t *)y, ny);
    1161   666944689 :   if (!hi) lz--;
    1162             :   /*else zd[lz-1]=hi; GH tell me it is not necessary.*/
    1163             : 
    1164   666944689 :   zd[1] = evalsigne(1) | evallgefint(lz);
    1165   666944689 :   return zd;
    1166             : }
    1167             : GEN
    1168      193154 : muluui(ulong x, ulong y, GEN z)
    1169             : {
    1170      193154 :   long t, s = signe(z);
    1171             :   GEN r;
    1172             :   LOCAL_HIREMAINDER;
    1173             : 
    1174      193154 :   if (!x || !y || !signe(z)) return gen_0;
    1175      192758 :   t = mulll(x,y);
    1176      192758 :   if (!hiremainder)
    1177      192758 :     r = muluispec(t, z+2, lgefint(z)-2);
    1178             :   else
    1179             :   {
    1180             :     long tmp[2];
    1181           0 :     tmp[1] = hiremainder;
    1182           0 :     tmp[0] = t;
    1183           0 :     r = muliispec(z+2,tmp, lgefint(z)-2, 2);
    1184             :   }
    1185      192758 :   setsigne(r,s); return r;
    1186             : }
    1187             : 
    1188             : GEN
    1189   668598314 : sqrispec(GEN x, long nx)
    1190             : {
    1191             :   GEN zd;
    1192             :   long lz;
    1193             : 
    1194   668598314 :   if (!nx) return gen_0;
    1195   184597381 :   if (nx==1) return sqru(*x);
    1196             : 
    1197   107004570 :   lz = (nx<<1)+2;
    1198   107004570 :   zd = cgeti(lz);
    1199             : #ifdef mpn_sqr
    1200   102551588 :   mpn_sqr(LIMBS(zd), (mp_limb_t *)x, nx);
    1201             : #else
    1202             :   mpn_mul_n(LIMBS(zd), (mp_limb_t *)x, (mp_limb_t *)x, nx);
    1203             : #endif
    1204   103923819 :   if (zd[lz-1]==0) lz--;
    1205             : 
    1206   103923819 :   zd[1] = evalsigne(1) | evallgefint(lz);
    1207   103923819 :   return zd;
    1208             : }
    1209             : 
    1210             : INLINE GEN
    1211     7475900 : sqrispec_mirror(GEN x, long nx)
    1212             : {
    1213     7475900 :   GEN cx=new_chunk(nx);
    1214             :   GEN z;
    1215     7475900 :   xmpn_mirrorcopy((mp_limb_t *)cx,(mp_limb_t *)x,nx);
    1216     7475900 :   z=sqrispec(cx, nx);
    1217     7475900 :   xmpn_mirror(LIMBS(z), NLIMBS(z));
    1218     7475900 :   return z;
    1219             : }
    1220             : 
    1221             : /* leaves garbage on the stack. */
    1222             : INLINE GEN
    1223    26188310 : muliispec_mirror(GEN x, GEN y, long nx, long ny)
    1224             : {
    1225             :   GEN cx, cy, z;
    1226    26188310 :   long s = 0;
    1227    26188310 :   while (nx && x[nx-1]==0) { nx--; s++; }
    1228    26188310 :   while (ny && y[ny-1]==0) { ny--; s++; }
    1229    26188310 :   cx=new_chunk(nx); cy=new_chunk(ny);
    1230    26188310 :   xmpn_mirrorcopy((mp_limb_t *)cx,(mp_limb_t *)x,nx);
    1231    26188310 :   xmpn_mirrorcopy((mp_limb_t *)cy,(mp_limb_t *)y,ny);
    1232    26188310 :   z =  nx>=ny ? muliispec(cx, cy, nx, ny): muliispec(cy, cx, ny, nx);
    1233    26188310 :   if (s)
    1234             :   {
    1235     2296506 :     long i, lz = lgefint(z) + s;
    1236     2296506 :     (void)new_chunk(s);
    1237     2296506 :     z -= s;
    1238     2296506 :     for (i=0; i<s; i++) z[2+i]=0;
    1239     2296506 :     z[1] = evalsigne(1) | evallgefint(lz);
    1240     2296506 :     z[0] = evaltyp(t_INT) | evallg(lz);
    1241             :   }
    1242    26188310 :   xmpn_mirror(LIMBS(z), NLIMBS(z));
    1243    26188310 :   return z;
    1244             : }
    1245             : 
    1246             : /* x % (2^n), assuming n >= 0 */
    1247             : GEN
    1248    26081894 : remi2n(GEN x, long n)
    1249             : {
    1250             :   ulong hi;
    1251    26081894 :   long l, k, lx, ly, sx = signe(x);
    1252             :   GEN z, xd, zd;
    1253             : 
    1254    26081894 :   if (!sx || !n) return gen_0;
    1255             : 
    1256    25928511 :   k = dvmdsBIL(n, &l);
    1257    26566878 :   lx = lgefint(x);
    1258    26566878 :   if (lx < k+3) return icopy(x);
    1259             : 
    1260    26255047 :   xd = x + (2 + k);
    1261             :   /* x = |k|...|1|#|... : copy the last l bits of # and the first k words
    1262             :    *              ^--- initial xd  */
    1263    26255047 :   hi = ((ulong)*xd) & ((1UL<<l)-1); /* last l bits of # = top bits of result */
    1264    26255047 :   if (!hi)
    1265             :   { /* strip leading zeroes from result */
    1266     1799929 :     xd--; while (k && !*xd) { k--; xd--; }
    1267     1799929 :     if (!k) return gen_0;
    1268     1031430 :     ly = k+2;
    1269             :   }
    1270             :   else
    1271    24455118 :     ly = k+3;
    1272             : 
    1273    25486548 :   zd = z = cgeti(ly);
    1274    25005858 :   *++zd = evalsigne(sx) | evallgefint(ly);
    1275    25005858 :   xd = x+1; for ( ;k; k--) *++zd = *++xd;
    1276    25005858 :   if (hi) *++zd = hi;
    1277    25005858 :   return z;
    1278             : }
    1279             : 
    1280             : /********************************************************************/
    1281             : /**                                                                **/
    1282             : /**                      INTEGER SQUARE ROOT                       **/
    1283             : /**                                                                **/
    1284             : /********************************************************************/
    1285             : 
    1286             : /* Return S (and set R) s.t S^2 + R = N, 0 <= R <= 2S.
    1287             :  * As for dvmdii, R is last on stack and guaranteed to be gen_0 in case the
    1288             :  * remainder is 0. R = NULL is allowed. */
    1289             : GEN
    1290     3756157 : sqrtremi(GEN a, GEN *r)
    1291             : {
    1292     3756157 :   long l, na = NLIMBS(a);
    1293             :   mp_size_t nr;
    1294             :   GEN S;
    1295     3756157 :   if (!na) {
    1296          60 :     if (r) *r = gen_0;
    1297          60 :     return gen_0;
    1298             :   }
    1299     3756097 :   l = (na + 5) >> 1; /* 2 + ceil(na/2) */
    1300     3756097 :   S = cgetipos(l);
    1301     3756050 :   if (r) {
    1302      349161 :     GEN R = cgeti(2 + na);
    1303      349161 :     nr = mpn_sqrtrem(LIMBS(S), LIMBS(R), LIMBS(a), na);
    1304      349161 :     if (nr) R[1] = evalsigne(1) | evallgefint(nr+2);
    1305       10948 :     else    { set_avma((pari_sp)S); R = gen_0; }
    1306      349161 :     *r = R;
    1307             :   }
    1308             :   else
    1309     3406889 :     (void)mpn_sqrtrem(LIMBS(S), NULL, LIMBS(a), na);
    1310     3756051 :   return S;
    1311             : }
    1312             : 
    1313             : /* compute sqrt(|a|), assuming a != 0 */
    1314             : GEN
    1315    28148229 : sqrtr_abs(GEN a)
    1316             : {
    1317             :   GEN res;
    1318             :   mp_limb_t *b, *c;
    1319    28148229 :   long l = RNLIMBS(a), e = expo(a), er = e>>1;
    1320             :   long n;
    1321    28148229 :   res = cgetr(2 + l);
    1322    28147207 :   res[1] = evalsigne(1) | evalexpo(er);
    1323    28148540 :   if (e&1)
    1324             :   {
    1325    13877509 :     b = (mp_limb_t *) new_chunk(l<<1);
    1326    13877392 :     xmpn_zero(b,l);
    1327    13877442 :     xmpn_mirrorcopy(b+l, RLIMBS(a), l);
    1328    13877596 :     c = (mp_limb_t *) new_chunk(l);
    1329    13877517 :     n = mpn_sqrtrem(c,b,b,l<<1); /* c <- sqrt; b <- rem */
    1330    13878025 :     if (n>l || (n==l && mpn_cmp(b,c,l) > 0)) mpn_add_1(c,c,l,1);
    1331             :   }
    1332             :   else
    1333             :   {
    1334             :     ulong u;
    1335    14271031 :     b = (mp_limb_t *) mantissa2nr(a,-1);
    1336    14271141 :     b[1] = uel(a,l+1)<<(BITS_IN_LONG-1);
    1337    14271141 :     b = (mp_limb_t *) new_chunk(l);
    1338    14271090 :     xmpn_zero(b,l+1); /* overwrites the former b[0] */
    1339    14271073 :     c = (mp_limb_t *) new_chunk(l + 1);
    1340    14270659 :     n = mpn_sqrtrem(c,b,b,(l<<1)+2); /* c <- sqrt; b <- rem */
    1341    14272306 :     u = (ulong)*c++;
    1342    14272306 :     if ( u&HIGHBIT || (u == ~HIGHBIT &&
    1343           0 :              (n>l || (n==l && mpn_cmp(b,c,l)>0))))
    1344     6981207 :       mpn_add_1(c,c,l,1);
    1345             :   }
    1346    28150342 :   xmpn_mirrorcopy(RLIMBS(res),c,l);
    1347    28149998 :   set_avma((pari_sp)res); return res;
    1348             : }
    1349             : 
    1350             : /* Normalize a non-negative integer */
    1351             : GEN
    1352   171370172 : int_normalize(GEN x, long known_zero_words)
    1353             : {
    1354   171370172 :   long i =  lgefint(x) - 1 - known_zero_words;
    1355  1126154471 :   for ( ; i > 1; i--)
    1356  1097374551 :     if (x[i]) { setlgefint(x, i+1); return x; }
    1357    28779920 :   x[1] = evalsigne(0) | evallgefint(2); return x;
    1358             : }
    1359             : 
    1360             : /********************************************************************
    1361             :  **                                                                **
    1362             :  **                           Base Conversion                      **
    1363             :  **                                                                **
    1364             :  ********************************************************************/
    1365             : 
    1366             : ulong *
    1367      343287 : convi(GEN x, long *l)
    1368             : {
    1369      343287 :   long n = nchar2nlong(2 + (long)(NLIMBS(x) * (BITS_IN_LONG * LOG10_2)));
    1370      343287 :   GEN str = cgetg(n+1, t_VECSMALL);
    1371      343287 :   unsigned char *res = (unsigned char*) GSTR(str);
    1372      343287 :   long llz = mpn_get_str(res, 10, LIMBS(icopy(x)), NLIMBS(x));
    1373             :   long lz;
    1374             :   ulong *z;
    1375             :   long i, j;
    1376             :   unsigned char *t;
    1377      343287 :   while (!*res) {res++; llz--;} /*Strip leading zeros*/
    1378      343287 :   lz  = (8+llz)/9;
    1379      343287 :   z = (ulong*)new_chunk(1+lz);
    1380      343287 :   t=res+llz+9;
    1381      738843 :   for(i=0;i<llz-8;i+=9)
    1382             :   {
    1383             :     ulong s;
    1384      395556 :     t-=18;
    1385      395556 :     s=*t++;
    1386     3560004 :     for (j=1; j<9;j++)
    1387     3164448 :       s=10*s+*t++;
    1388      395556 :     *z++=s;
    1389             :   }
    1390      343287 :   if (i<llz)
    1391             :   {
    1392      339616 :     unsigned char *t = res;
    1393      339616 :     ulong s=*t++;
    1394      876305 :     for (j=i+1; j<llz;j++)
    1395      536689 :       s=10*s+*t++;
    1396      339616 :     *z++=s;
    1397             :   }
    1398      343287 :   *l = lz;
    1399      343287 :   return z;
    1400             : }

Generated by: LCOV version 1.13