Code coverage tests

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

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

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

LCOV - code coverage report
Current view: top level - kernel/gmp - mp.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.8.0 lcov report (development 19352-1b11b25) Lines: 639 716 89.2 %
Date: 2016-08-25 06:11:27 Functions: 53 54 98.1 %
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      340738 : static void pari_gmp_free(void *ptr, size_t old_size){
      51      340738 :   (void)old_size; pari_free(ptr);
      52      340738 : }
      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         768 : pari_kernel_init(void)
      60             : {
      61         768 :   mp_get_memory_functions (&old_gmp_malloc, &old_gmp_realloc, &old_gmp_free);
      62         768 :   mp_set_memory_functions((void *(*)(size_t)) pari_malloc, pari_gmp_realloc, pari_gmp_free);
      63         768 : }
      64             : 
      65             : void
      66        1276 : 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        1276 :   mp_get_memory_functions (&new_gmp_malloc, &new_gmp_realloc, &new_gmp_free);
      72        1276 :   if (new_gmp_malloc==pari_malloc) new_gmp_malloc = old_gmp_malloc;
      73        1276 :   if (new_gmp_realloc==pari_gmp_realloc) new_gmp_realloc = old_gmp_realloc;
      74        1276 :   if (new_gmp_free==pari_gmp_free) new_gmp_free = old_gmp_free;
      75        1276 :   mp_set_memory_functions(new_gmp_malloc, new_gmp_realloc, new_gmp_free);
      76        1276 : }
      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     1254382 : xmpn_copy(mp_limb_t *x, mp_limb_t *y, long n)
      86             : {
      87     1254382 :   while (--n >= 0) x[n]=y[n];
      88     1254382 : }
      89             : 
      90             : INLINE void
      91   208489915 : xmpn_mirror(mp_limb_t *x, long n)
      92             : {
      93             :   long i;
      94  1025680221 :   for(i=0;i<(n>>1);i++)
      95             :   {
      96   817190306 :     ulong m=x[i];
      97   817190306 :     x[i]=x[n-1-i];
      98   817190306 :     x[n-1-i]=m;
      99             :   }
     100   208489915 : }
     101             : 
     102             : INLINE void
     103   269829604 : xmpn_mirrorcopy(mp_limb_t *z, mp_limb_t *x, long n)
     104             : {
     105             :   long i;
     106  2447333296 :   for(i=0;i<n;i++)
     107  2177503692 :     z[i]=x[n-1-i];
     108   269829604 : }
     109             : 
     110             : INLINE void
     111    55813159 : xmpn_zero(mp_ptr x, mp_size_t n)
     112             : {
     113    55813159 :   while (--n >= 0) x[n]=0;
     114    55813159 : }
     115             : 
     116             : INLINE GEN
     117    18726543 : icopy_ef(GEN x, long l)
     118             : {
     119    18726543 :   register long lx = lgefint(x);
     120    18726543 :   const GEN y = cgeti(l);
     121             : 
     122    18726521 :   while (--lx > 0) y[lx]=x[lx];
     123    18726521 :   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     2391187 : setloop(GEN a)
     141             : {
     142     2391187 :   pari_sp av = avma - 2 * sizeof(long);
     143     2391187 :   (void)cgetg(lgefint(a) + 3, t_VECSMALL);
     144     2391191 :   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      156728 : resetloop(GEN a, GEN b) {
     150      156728 :   a[0] = evaltyp(t_INT) | evallg(lgefint(b));
     151      156728 :   affii(b, a); return a;
     152             : }
     153             : 
     154             : /* assume a > 0, initialized by setloop. Do a++ */
     155             : static GEN
     156    32556969 : incpos(GEN a)
     157             : {
     158    32556969 :   long i, l = lgefint(a);
     159    32556974 :   for (i=2; i<l; i++)
     160    32556970 :     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        6184 : incneg(GEN a)
     170             : {
     171        6184 :   long i, l = lgefint(a)-1;
     172        6184 :   if (uel(a,2)--)
     173             :   {
     174        6180 :     if (!a[l]) /* implies l = 2 */
     175             :     {
     176         532 :       a[0] = evaltyp(t_INT) | _evallg(2);
     177         532 :       a[1] = evalsigne(0) | evallgefint(2);
     178             :     }
     179        6180 :     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    32576043 : incloop(GEN a)
     194             : {
     195    32576043 :   switch(signe(a))
     196             :   {
     197             :     case 0:
     198       20688 :       a[0]=evaltyp(t_INT) | _evallg(3);
     199       20688 :       a[1]=evalsigne(1) | evallgefint(3);
     200       20688 :       a[2]=1; return a;
     201        6184 :     case -1: return incneg(a);
     202    32549171 :     default: return incpos(a);
     203             :   }
     204             : }
     205             : 
     206             : INLINE GEN
     207  1307665975 : adduispec(ulong s, GEN x, long nx)
     208             : {
     209             :   GEN  zd;
     210             :   long lz;
     211             : 
     212  1307665975 :   if (nx == 1) return adduu(uel(x,0), s);
     213   300034274 :   lz = nx+3; zd = cgeti(lz);
     214   299427740 :   if (mpn_add_1(LIMBS(zd),(mp_limb_t *)x,nx,s))
     215     5047041 :     zd[lz-1]=1;
     216             :   else
     217   294473699 :     lz--;
     218   299520740 :   zd[1] = evalsigne(1) | evallgefint(lz);
     219   299520740 :   return zd;
     220             : }
     221             : 
     222             : GEN
     223   325958427 : adduispec_offset(ulong s, GEN x, long offset, long nx)
     224             : {
     225   325958427 :   GEN xd=x+2+offset;
     226   325958427 :   while (nx && *(xd+nx-1)==0) nx--;
     227   325958427 :   if (!nx) return utoi(s);
     228   306792383 :   return adduispec(s,xd,nx);
     229             : }
     230             : 
     231             : INLINE GEN
     232  1668973803 : addiispec(GEN x, GEN y, long nx, long ny)
     233             : {
     234             :   GEN zd;
     235             :   long lz;
     236             : 
     237  1668973803 :   if (nx < ny) swapspec(x,y, nx,ny);
     238  1668973803 :   if (ny == 1) return adduispec(*y,x,nx);
     239   714830554 :   lz = nx+3; zd = cgeti(lz);
     240             : 
     241   714164209 :   if (mpn_add(LIMBS(zd),(mp_limb_t *)x,nx,(mp_limb_t *)y,ny))
     242    11145905 :     zd[lz-1]=1;
     243             :   else
     244   703475722 :     lz--;
     245             : 
     246   714621627 :   zd[1] = evalsigne(1) | evallgefint(lz);
     247   714621627 :   return zd;
     248             : }
     249             : 
     250             : /* assume x >= y */
     251             : INLINE GEN
     252  1089809826 : subiuspec(GEN x, ulong s, long nx)
     253             : {
     254             :   GEN zd;
     255             :   long lz;
     256             : 
     257  1089809826 :   if (nx == 1) return utoi(x[0] - s);
     258             : 
     259    77659727 :   lz = nx + 2; zd = cgeti(lz);
     260    77646632 :   mpn_sub_1 (LIMBS(zd), (mp_limb_t *)x, nx, s);
     261    77647438 :   if (! zd[lz - 1]) { --lz; }
     262             : 
     263    77647438 :   zd[1] = evalsigne(1) | evallgefint(lz);
     264    77647438 :   return zd;
     265             : }
     266             : 
     267             : /* assume x > y */
     268             : INLINE GEN
     269  1771261911 : subiispec(GEN x, GEN y, long nx, long ny)
     270             : {
     271             :   GEN zd;
     272             :   long lz;
     273  1771261911 :   if (ny==1) return subiuspec(x,*y,nx);
     274   690821064 :   lz = nx+2; zd = cgeti(lz);
     275             : 
     276   687344059 :   mpn_sub (LIMBS(zd), (mp_limb_t *)x, nx, (mp_limb_t *)y, ny);
     277   689032331 :   while (lz >= 3 && zd[lz - 1] == 0) { lz--; }
     278             : 
     279   689032331 :   zd[1] = evalsigne(1) | evallgefint(lz);
     280   689032331 :   return zd;
     281             : }
     282             : 
     283             : static void
     284   272129474 : roundr_up_ip(GEN x, long l)
     285             : {
     286   272129474 :   long i = l;
     287             :   for(;;)
     288             :   {
     289   272383072 :     if (++((ulong*)x)[--i]) break;
     290      327013 :     if (i == 2) { x[2] = HIGHBIT; shiftr_inplace(x, 1); break; }
     291      253598 :   }
     292   272129474 : }
     293             : 
     294             : void
     295   132260270 : affir(GEN x, GEN y)
     296             : {
     297   132260270 :   const long s = signe(x), ly = lg(y);
     298             :   long lx, sh, i;
     299             : 
     300   132260270 :   if (!s)
     301             :   {
     302     1537911 :     y[1] = evalexpo(-bit_accuracy(ly));
     303     1537911 :     return;
     304             :   }
     305   130722359 :   lx = lgefint(x); sh = bfffo(*int_MSW(x));
     306   130722359 :   y[1] = evalsigne(s) | evalexpo(bit_accuracy(lx)-sh-1);
     307   130722359 :   if (sh) {
     308   129014897 :     if (lx <= ly)
     309             :     {
     310    82157220 :       for (i=lx; i<ly; i++) y[i]=0;
     311    82157220 :       mpn_lshift(LIMBS(y),LIMBS(x),lx-2,sh);
     312    82157220 :       xmpn_mirror(LIMBS(y),lx-2);
     313    82157220 :       return;
     314             :     }
     315    46857677 :     mpn_lshift(LIMBS(y),LIMBS(x)+lx-ly,ly-2,sh);
     316    46857677 :     uel(y,2) |= uel(x,lx-ly+1) >> (BITS_IN_LONG-sh);
     317    46857677 :     xmpn_mirror(LIMBS(y),ly-2);
     318             :     /* lx > ly: round properly */
     319    46857677 :     if ((uel(x,lx-ly+1)<<sh) & HIGHBIT) roundr_up_ip(y, ly);
     320             :   }
     321             :   else {
     322     1707462 :     GEN xd=int_MSW(x);
     323     1707462 :     if (lx <= ly)
     324             :     {
     325      976325 :       for (i=2; i<lx; i++,xd=int_precW(xd)) y[i]=*xd;
     326      976325 :       for (   ; i<ly; i++) y[i]=0;
     327      976325 :       return;
     328             :     }
     329      731137 :     for (i=2; i<ly; i++,xd=int_precW(xd)) y[i]=*xd;
     330             :     /* lx > ly: round properly */
     331      731137 :     if (uel(x,lx-ly+1) & HIGHBIT) roundr_up_ip(y, ly);
     332             :   }
     333             : }
     334             : 
     335             : INLINE GEN
     336   351962373 : shiftispec(GEN x, long nx, long n)
     337             : {
     338             :   long ny,m;
     339             :   GEN yd, y;
     340             : 
     341   351962373 :   if (!n) return icopyspec(x, nx);
     342             : 
     343   336343357 :   if (n > 0)
     344             :   {
     345   224698545 :     long d = dvmdsBIL(n, &m);
     346             :     long i;
     347             : 
     348   224671273 :     ny = nx + d + (m!=0);
     349   224671273 :     y = cgeti(ny + 2); yd = y + 2;
     350   224510922 :     for (i=0; i<d; i++) yd[i] = 0;
     351             : 
     352   224510922 :     if (!m) xmpn_copy((mp_limb_t *) (yd + d), (mp_limb_t *) x, nx);
     353             :     else
     354             :     {
     355   223787887 :       ulong carryd = mpn_lshift((mp_limb_t *) (yd + d), (mp_limb_t *) x, nx, m);
     356   223872062 :       if (carryd) yd[ny - 1] = carryd;
     357   211643157 :       else ny--;
     358             :     }
     359             :   }
     360             :   else
     361             :   {
     362   111644812 :     long d = dvmdsBIL(-n, &m);
     363             : 
     364   111644783 :     ny = nx - d;
     365   111644783 :     if (ny < 1) return gen_0;
     366   111305648 :     y = cgeti(ny + 2); yd = y + 2;
     367             : 
     368   111305034 :     if (!m) xmpn_copy((mp_limb_t *) yd, (mp_limb_t *) (x + d), nx - d);
     369             :     else
     370             :     {
     371   110859798 :       mpn_rshift((mp_limb_t *) yd, (mp_limb_t *) (x + d), nx - d, m);
     372   110859820 :       if (yd[ny - 1] == 0)
     373             :       {
     374     8592149 :         if (ny == 1) { avma = (pari_sp)(yd + 1); return gen_0; }
     375     7693753 :         ny--;
     376             :       }
     377             :     }
     378             :   }
     379   334978480 :   y[1] = evalsigne(1)|evallgefint(ny + 2);
     380   334978480 :   return y;
     381             : }
     382             : 
     383             : GEN
     384    23510578 : mantissa2nr(GEN x, long n)
     385             : {
     386    23510578 :   long ly, i, m, s = signe(x), lx = lg(x);
     387             :   GEN y;
     388    23510578 :   if (!s) return gen_0;
     389    23509265 :   if (!n)
     390             :   {
     391     9884850 :     y = cgeti(lx);
     392     9884850 :     y[1] = evalsigne(s) | evallgefint(lx);
     393     9884850 :     xmpn_mirrorcopy(LIMBS(y),RLIMBS(x),lx-2);
     394     9884850 :     return y;
     395             :   }
     396    13624415 :   if (n > 0)
     397             :   {
     398      568756 :     GEN z = (GEN)avma;
     399      568756 :     long d = dvmdsBIL(n, &m);
     400             : 
     401      568756 :     ly = lx+d; y = new_chunk(ly);
     402      568756 :     for ( ; d; d--) *--z = 0;
     403      568756 :     if (!m) for (i=2; i<lx; i++) y[i]=x[i];
     404             :     else
     405             :     {
     406      568198 :       register const ulong sh = BITS_IN_LONG - m;
     407      568198 :       shift_left(y,x, 2,lx-1, 0,m);
     408      568198 :       i = uel(x,2) >> sh;
     409             :       /* Extend y on the left? */
     410      568198 :       if (i) { ly++; y = new_chunk(1); y[2] = i; }
     411             :     }
     412             :   }
     413             :   else
     414             :   {
     415    13055659 :     ly = lx - dvmdsBIL(-n, &m);
     416    13055659 :     if (ly<3) return gen_0;
     417    13055659 :     y = new_chunk(ly);
     418    13055659 :     if (m) {
     419    13025127 :       shift_right(y,x, 2,ly, 0,m);
     420    13025127 :       if (y[2] == 0)
     421             :       {
     422           0 :         if (ly==3) { avma = (pari_sp)(y+3); return gen_0; }
     423           0 :         ly--; avma = (pari_sp)(++y);
     424             :       }
     425             :     } else {
     426       30532 :       for (i=2; i<ly; i++) y[i]=x[i];
     427             :     }
     428             :   }
     429    13624415 :   xmpn_mirror(LIMBS(y),ly-2);
     430    13624415 :   y[1] = evalsigne(s)|evallgefint(ly);
     431    13624415 :   y[0] = evaltyp(t_INT)|evallg(ly); return y;
     432             : }
     433             : 
     434             : GEN
     435     1452845 : truncr(GEN x)
     436             : {
     437             :   long s, e, d, m, i;
     438             :   GEN y;
     439     1452845 :   if ((s=signe(x)) == 0 || (e=expo(x)) < 0) return gen_0;
     440     1284127 :   d = nbits2prec(e+1); m = remsBIL(e);
     441     1284127 :   if (d > lg(x)) pari_err_PREC( "truncr (precision loss in truncation)");
     442             : 
     443     1284127 :   y=cgeti(d); y[1] = evalsigne(s) | evallgefint(d);
     444     1284127 :   if (++m == BITS_IN_LONG)
     445         154 :     for (i=2; i<d; i++) y[d-i+1]=x[i];
     446             :   else
     447             :   {
     448     1283973 :     GEN z=cgeti(d);
     449     1283973 :     for (i=2; i<d; i++) z[d-i+1]=x[i];
     450     1283973 :     mpn_rshift(LIMBS(y),LIMBS(z),d-2,BITS_IN_LONG-m);
     451     1283973 :     avma=(pari_sp)y;
     452             :   }
     453     1284127 :   return y;
     454             : }
     455             : 
     456             : /* integral part */
     457             : GEN
     458      401451 : floorr(GEN x)
     459             : {
     460             :   long e, d, m, i, lx;
     461             :   GEN y;
     462      401451 :   if (signe(x) >= 0) return truncr(x);
     463      160905 :   if ((e=expo(x)) < 0) return gen_m1;
     464       66827 :   d = nbits2prec(e+1); m = remsBIL(e);
     465       66827 :   lx=lg(x); if (d>lx) pari_err_PREC( "floorr (precision loss in truncation)");
     466       66827 :   y = cgeti(d+1);
     467       66827 :   if (++m == BITS_IN_LONG)
     468             :   {
     469         171 :     for (i=2; i<d; i++) y[d-i+1]=x[i];
     470         171 :     i=d; while (i<lx && !x[i]) i++;
     471         171 :     if (i==lx) goto END;
     472             :   }
     473             :   else
     474             :   {
     475       66656 :     GEN z=cgeti(d);
     476       66656 :     for (i=2; i<d; i++) z[d-i+1]=x[i];
     477       66656 :     mpn_rshift(LIMBS(y),LIMBS(z),d-2,BITS_IN_LONG-m);
     478       66656 :     if (uel(x,d-1)<<m == 0)
     479             :     {
     480       15313 :       i=d; while (i<lx && !x[i]) i++;
     481       15313 :       if (i==lx) goto END;
     482             :     }
     483             :   }
     484       60452 :   if (mpn_add_1(LIMBS(y),LIMBS(y),d-2,1))
     485           0 :     y[d++]=1;
     486             : END:
     487       66827 :   y[1] = evalsigne(-1) | evallgefint(d);
     488       66827 :   return y;
     489             : }
     490             : 
     491             : INLINE int
     492  2109478748 : cmpiispec(GEN x, GEN y, long lx, long ly)
     493             : {
     494  2109478748 :   if (lx < ly) return -1;
     495  2028781223 :   if (lx > ly) return  1;
     496  1949223698 :   return mpn_cmp((mp_limb_t*)x,(mp_limb_t*)y, lx);
     497             : }
     498             : 
     499             : INLINE int
     500    56492628 : equaliispec(GEN x, GEN y, long lx, long ly)
     501             : {
     502    56492628 :   if (lx != ly) return 0;
     503    56490159 :   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  2124540786 : muluispec(ulong x, GEN y, long ny)
     514             : {
     515  2124540786 :   if (ny == 1)
     516  1660560040 :     return muluu(x, *y);
     517             :   else
     518             :   {
     519   463980746 :     long lz = ny+3;
     520   463980746 :     GEN z = cgeti(lz);
     521   461592255 :     ulong hi = mpn_mul_1 (LIMBS(z), (mp_limb_t *)y, ny, x);
     522   462179591 :     if (hi) { z[lz - 1] = hi; } else lz--;
     523   462179591 :     z[1] = evalsigne(1) | evallgefint(lz);
     524   462179591 :     return z;
     525             :   }
     526             : }
     527             : 
     528             : /* a + b*|y| */
     529             : GEN
     530      759654 : addumului(ulong a, ulong b, GEN y)
     531             : {
     532             :   GEN z;
     533             :   long i, lz;
     534             :   ulong hi;
     535      759654 :   if (!b || !signe(y)) return utoi(a);
     536      759588 :   lz = lgefint(y)+1;
     537      759588 :   z = cgeti(lz);
     538      759588 :   z[2]=a;
     539      759588 :   for(i=3;i<lz;i++) z[i]=0;
     540      759588 :   hi=mpn_addmul_1(LIMBS(z), LIMBS(y), NLIMBS(y), b);
     541      759588 :   if (hi) z[lz-1]=hi; else lz--;
     542      759588 :   z[1] = evalsigne(1) | evallgefint(lz);
     543      759588 :   avma=(pari_sp)z; return z;
     544             : }
     545             : 
     546             : /***********************************************************************/
     547             : /**                                                                   **/
     548             : /**                          DIVISION                                 **/
     549             : /**                                                                   **/
     550             : /***********************************************************************/
     551             : 
     552             : ulong
     553   778879564 : umodiu(GEN y, ulong x)
     554             : {
     555   778879564 :   long sy=signe(y);
     556             :   ulong hi;
     557   778879564 :   if (!x) pari_err_INV("umodiu",gen_0);
     558   778900833 :   if (!sy) return 0;
     559   592593633 :   hi = mpn_mod_1(LIMBS(y),NLIMBS(y),x);
     560   592601338 :   if (!hi) return 0;
     561   584717315 :   return (sy > 0)? hi: x - hi;
     562             : }
     563             : 
     564             : /* return |y| \/ x */
     565             : GEN
     566   125759009 : diviu_rem(GEN y, ulong x, ulong *rem)
     567             : {
     568             :   long ly;
     569             :   GEN z;
     570             : 
     571   125759009 :   if (!x) pari_err_INV("diviu_rem",gen_0);
     572   125787139 :   if (!signe(y)) { *rem = 0; return gen_0; }
     573             : 
     574   125251311 :   ly = lgefint(y);
     575   125251311 :   if (ly == 3 && (ulong)x > uel(y,2)) { *rem = uel(y,2); return gen_0; }
     576             : 
     577   125200859 :   z = cgeti(ly);
     578   125166084 :   *rem = mpn_divrem_1(LIMBS(z), 0, LIMBS(y), NLIMBS(y), x);
     579   125200695 :   if (z [ly - 1] == 0) ly--;
     580   125200695 :   z[1] = evallgefint(ly) | evalsigne(1);
     581   125200695 :   return z;
     582             : }
     583             : 
     584             : GEN
     585    32204180 : divis_rem(GEN y, long x, long *rem)
     586             : {
     587    32204180 :   long sy=signe(y),ly,s;
     588             :   GEN z;
     589             : 
     590    32204180 :   if (!x) pari_err_INV("divis_rem",gen_0);
     591    32204180 :   if (!sy) { *rem = 0; return gen_0; }
     592    25495785 :   if (x<0) { s = -sy; x = -x; } else s = sy;
     593             : 
     594    25495785 :   ly = lgefint(y);
     595    25495785 :   if (ly == 3 && (ulong)x > uel(y,2)) { *rem = itos(y); return gen_0; }
     596             : 
     597    19786774 :   z = cgeti(ly);
     598    19786774 :   *rem = mpn_divrem_1(LIMBS(z), 0, LIMBS(y), NLIMBS(y), x);
     599    19786774 :   if (sy<0) *rem = -  *rem;
     600    19786774 :   if (z[ly - 1] == 0) ly--;
     601    19786774 :   z[1] = evallgefint(ly) | evalsigne(s);
     602    19786774 :   return z;
     603             : }
     604             : 
     605             : GEN
     606      461309 : divis(GEN y, long x)
     607             : {
     608      461309 :   long sy=signe(y),ly,s;
     609             :   GEN z;
     610             : 
     611      461309 :   if (!x) pari_err_INV("divis",gen_0);
     612      461309 :   if (!sy) return gen_0;
     613      461309 :   if (x<0) { s = -sy; x = -x; } else s=sy;
     614             : 
     615      461309 :   ly = lgefint(y);
     616      461309 :   if (ly == 3 && (ulong)x > uel(y,2)) return gen_0;
     617             : 
     618      461309 :   z = cgeti(ly);
     619      461309 :   (void)mpn_divrem_1(LIMBS(z), 0, LIMBS(y), NLIMBS(y), x);
     620      461309 :   if (z[ly - 1] == 0) ly--;
     621      461309 :   z[1] = evallgefint(ly) | evalsigne(s);
     622      461309 :   return z;
     623             : }
     624             : 
     625             : /* We keep llx bits of x and lly bits of y*/
     626             : static GEN
     627    31174734 : divrr_with_gmp(GEN x, GEN y)
     628             : {
     629    31174734 :   long lx=RNLIMBS(x),ly=RNLIMBS(y);
     630    31174734 :   long lw=minss(lx,ly);
     631    31174734 :   long lly=minss(lw+1,ly);
     632    31174734 :   GEN  w=cgetr(lw+2);
     633    31174734 :   long lu=lw+lly;
     634    31174734 :   long llx=minss(lu,lx);
     635    31174734 :   mp_limb_t *u=(mp_limb_t *)new_chunk(lu);
     636    31174734 :   mp_limb_t *z=(mp_limb_t *)new_chunk(lly);
     637             :   mp_limb_t *q,*r;
     638    31174734 :   long e=expo(x)-expo(y);
     639    31174734 :   long sx=signe(x),sy=signe(y);
     640    31174734 :   xmpn_mirrorcopy(z,RLIMBS(y),lly);
     641    31174734 :   xmpn_mirrorcopy(u+lu-llx,RLIMBS(x),llx);
     642    31174734 :   xmpn_zero(u,lu-llx);
     643    31174734 :   q = (mp_limb_t *)new_chunk(lw+1);
     644    31174734 :   r = (mp_limb_t *)new_chunk(lly);
     645             : 
     646    31174734 :   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    31174734 :   if (uel(r,lly-1) > (uel(z,lly-1)>>1))
     650    15423459 :     mpn_add_1(q,q,lw+1,1);
     651             : 
     652    31174734 :   xmpn_mirrorcopy(RLIMBS(w),q,lw);
     653             : 
     654    31174734 :   if (q[lw] == 0) e--;
     655    11748465 :   else if (q[lw] == 1) { shift_right(w,w, 2,lw+2, 1,1); }
     656          24 :   else { w[2] = HIGHBIT; e++; }
     657    31174734 :   if (sy < 0) sx = -sx;
     658    31174734 :   w[1] = evalsigne(sx) | evalexpo(e);
     659    31174734 :   avma=(pari_sp) w;
     660    31174734 :   return w;
     661             : }
     662             : 
     663             : /* We keep llx bits of x and lly bits of y*/
     664             : static GEN
     665     3582356 : divri_with_gmp(GEN x, GEN y)
     666             : {
     667     3582356 :   long llx=RNLIMBS(x),ly=NLIMBS(y);
     668     3582356 :   long lly=minss(llx+1,ly);
     669     3582356 :   GEN  w=cgetr(llx+2);
     670     3582356 :   long lu=llx+lly, ld=ly-lly;
     671     3582356 :   mp_limb_t *u=(mp_limb_t *)new_chunk(lu);
     672     3582356 :   mp_limb_t *z=(mp_limb_t *)new_chunk(lly);
     673             :   mp_limb_t *q,*r;
     674     3582356 :   long sh=bfffo(y[ly+1]);
     675     3582356 :   long e=expo(x)-expi(y);
     676     3582356 :   long sx=signe(x),sy=signe(y);
     677     3582356 :   if (sh) mpn_lshift(z,LIMBS(y)+ld,lly,sh);
     678       86110 :   else xmpn_copy(z,LIMBS(y)+ld,lly);
     679     3582356 :   xmpn_mirrorcopy(u+lu-llx,RLIMBS(x),llx);
     680     3582356 :   xmpn_zero(u,lu-llx);
     681     3582356 :   q = (mp_limb_t *)new_chunk(llx+1);
     682     3582356 :   r = (mp_limb_t *)new_chunk(lly);
     683             : 
     684     3582356 :   mpn_tdiv_qr(q,r,0,u,lu,z,lly);
     685             : 
     686             :   /*Round up: This is not exactly correct we should test 2*r>z*/
     687     3582356 :   if (uel(r,lly-1) > (uel(z,lly-1)>>1))
     688     1604082 :     mpn_add_1(q,q,llx+1,1);
     689             : 
     690     3582356 :   xmpn_mirrorcopy(RLIMBS(w),q,llx);
     691             : 
     692     3582356 :   if (q[llx] == 0) e--;
     693     2863956 :   else if (q[llx] == 1) { shift_right(w,w, 2,llx+2, 1,1); }
     694           0 :   else { w[2] = HIGHBIT; e++; }
     695     3582356 :   if (sy < 0) sx = -sx;
     696     3582356 :   w[1] = evalsigne(sx) | evalexpo(e);
     697     3582356 :   avma=(pari_sp) w;
     698     3582356 :   return w;
     699             : }
     700             : 
     701             : GEN
     702    17711975 : divri(GEN x, GEN y)
     703             : {
     704    17711975 :   long  s = signe(y);
     705             : 
     706    17711975 :   if (!s) pari_err_INV("divri",gen_0);
     707    17711975 :   if (!signe(x)) return real_0_bit(expo(x) - expi(y));
     708    17678506 :   if (!is_bigint(y)) {
     709    14096150 :     GEN z = divru(x, y[2]);
     710    14096150 :     if (s < 0) togglesign(z);
     711    14096150 :     return z;
     712             :   }
     713     3582356 :   return divri_with_gmp(x,y);
     714             : }
     715             : 
     716             : GEN
     717   107204821 : divrr(GEN x, GEN y)
     718             : {
     719   107204821 :   long sx=signe(x), sy=signe(y), lx,ly,lr,e,i,j;
     720             :   ulong y0,y1;
     721             :   GEN r, r1;
     722             : 
     723   107204821 :   if (!sy) pari_err_INV("divrr",y);
     724   107204821 :   e = expo(x) - expo(y);
     725   107204821 :   if (!sx) return real_0_bit(e);
     726   105296692 :   if (sy<0) sx = -sx;
     727             : 
     728   105296692 :   lx=lg(x); ly=lg(y);
     729   105296692 :   if (ly==3)
     730             :   {
     731    74121958 :     ulong k = x[2], l = (lx>3)? x[3]: 0;
     732             :     LOCAL_HIREMAINDER;
     733    74121958 :     if (k < uel(y,2)) e--;
     734             :     else
     735             :     {
     736    34224607 :       l >>= 1; if (k&1) l |= HIGHBIT;
     737    34224607 :       k >>= 1;
     738             :     }
     739    74121958 :     hiremainder = k; k = divll(l,y[2]);
     740    74121958 :     if (hiremainder & HIGHBIT)
     741             :     {
     742    20560523 :       k++;
     743    20560523 :       if (!k) { k = HIGHBIT; e++; }
     744             :     }
     745    74121958 :     r = cgetr(3);
     746    74121958 :     r[1] = evalsigne(sx) | evalexpo(e);
     747    74121958 :     r[2] = k; return r;
     748             :   }
     749             : 
     750    31174734 :   if (ly>=DIVRR_GMP_LIMIT)
     751    31174734 :     return divrr_with_gmp(x,y);
     752             : 
     753           0 :   lr = minss(lx,ly); r = new_chunk(lr);
     754           0 :   r1 = r-1;
     755           0 :   r1[1] = 0; for (i=2; i<lr; i++) r1[i]=x[i];
     756           0 :   r1[lr] = (lx>ly)? x[lr]: 0;
     757           0 :   y0 = y[2]; y1 = y[3];
     758           0 :   for (i=0; i<lr-1; i++)
     759             :   { /* r1 = r + (i-1), OK up to r1[2] (accesses at most r[lr]) */
     760             :     ulong k, qp;
     761             :     LOCAL_HIREMAINDER;
     762             :     LOCAL_OVERFLOW;
     763             : 
     764           0 :     if (uel(r1,1) == y0)
     765             :     {
     766           0 :       qp = ULONG_MAX; k = addll(y0,r1[2]);
     767             :     }
     768             :     else
     769             :     {
     770           0 :       if (uel(r1,1) > y0) /* can't happen if i=0 */
     771             :       {
     772           0 :         GEN y1 = y+1;
     773           0 :         j = lr-i; r1[j] = subll(r1[j],y1[j]);
     774           0 :         for (j--; j>0; j--) r1[j] = subllx(r1[j],y1[j]);
     775           0 :         j=i; do uel(r,--j)++; while (j && !r[j]);
     776             :       }
     777           0 :       hiremainder = r1[1]; overflow = 0;
     778           0 :       qp = divll(r1[2],y0); k = hiremainder;
     779             :     }
     780           0 :     j = lr-i+1;
     781           0 :     if (!overflow)
     782             :     {
     783             :       long k3, k4;
     784           0 :       k3 = mulll(qp,y1);
     785           0 :       if (j == 3) /* i = lr - 2 maximal, r1[3] undefined -> 0 */
     786           0 :         k4 = subll(hiremainder,k);
     787             :       else
     788             :       {
     789           0 :         k3 = subll(k3, r1[3]);
     790           0 :         k4 = subllx(hiremainder,k);
     791             :       }
     792           0 :       while (!overflow && k4) { qp--; k3=subll(k3,y1); k4=subllx(k4,y0); }
     793             :     }
     794           0 :     if (j<ly) (void)mulll(qp,y[j]); else { hiremainder = 0 ; j = ly; }
     795           0 :     for (j--; j>1; j--)
     796             :     {
     797           0 :       r1[j] = subll(r1[j], addmul(qp,y[j]));
     798           0 :       hiremainder += overflow;
     799             :     }
     800           0 :     if (uel(r1,1) != hiremainder)
     801             :     {
     802           0 :       if (uel(r1,1) < hiremainder)
     803             :       {
     804           0 :         qp--;
     805           0 :         j = lr-i-(lr-i>=ly); r1[j] = addll(r1[j], y[j]);
     806           0 :         for (j--; j>1; j--) r1[j] = addllx(r1[j], y[j]);
     807             :       }
     808             :       else
     809             :       {
     810           0 :         uel(r1,1) -= hiremainder;
     811           0 :         while (r1[1])
     812             :         {
     813           0 :           qp++; if (!qp) { j=i; do uel(r,--j)++; while (j && !r[j]); }
     814           0 :           j = lr-i-(lr-i>=ly); r1[j] = subll(r1[j],y[j]);
     815           0 :           for (j--; j>1; j--) r1[j] = subllx(r1[j],y[j]);
     816           0 :           uel(r1,1) -= overflow;
     817             :         }
     818             :       }
     819             :     }
     820           0 :     *++r1 = qp;
     821             :   }
     822             :   /* i = lr-1 */
     823             :   /* round correctly */
     824           0 :   if (uel(r1,1) > (y0>>1))
     825             :   {
     826           0 :     j=i; do uel(r,--j)++; while (j && !r[j]);
     827             :   }
     828           0 :   r1 = r-1; for (j=i; j>=2; j--) r[j]=r1[j];
     829           0 :   if (r[0] == 0) e--;
     830           0 :   else if (r[0] == 1) { shift_right(r,r, 2,lr, 1,1); }
     831             :   else { /* possible only when rounding up to 0x2 0x0 ... */
     832           0 :     r[2] = (long)HIGHBIT; e++;
     833             :   }
     834           0 :   r[0] = evaltyp(t_REAL)|evallg(lr);
     835           0 :   r[1] = evalsigne(sx) | evalexpo(e);
     836           0 :   return r;
     837             : }
     838             : 
     839             : /* Integer division x / y: such that sign(r) = sign(x)
     840             :  *   if z = ONLY_REM return remainder, otherwise return quotient
     841             :  *   if z != NULL set *z to remainder
     842             :  *   *z is the last object on stack (and thus can be disposed of with cgiv
     843             :  *   instead of gerepile)
     844             :  * If *z is zero, we put gen_0 here and no copy.
     845             :  * space needed: lx + ly
     846             :  */
     847             : GEN
     848   828070690 : dvmdii(GEN x, GEN y, GEN *z)
     849             : {
     850   828070690 :   long sx=signe(x),sy=signe(y);
     851             :   long lx, ly, lq;
     852             :   pari_sp av;
     853             :   GEN r,q;
     854             : 
     855   828070690 :   if (!sy) pari_err_INV("dvmdii",y);
     856   829741351 :   if (!sx)
     857             :   {
     858    17463655 :     if (!z || z == ONLY_REM) return gen_0;
     859     7589152 :     *z=gen_0; return gen_0;
     860             :   }
     861   812277696 :   lx=lgefint(x);
     862   812277696 :   ly=lgefint(y); lq=lx-ly;
     863   812277696 :   if (lq <= 0)
     864             :   {
     865   453888788 :     if (lq == 0)
     866             :     {
     867   436248099 :       long s=mpn_cmp(LIMBS(x),LIMBS(y),NLIMBS(x));
     868   436681620 :       if (s>0) goto DIVIDE;
     869   150435119 :       if (s==0)
     870             :       {
     871     8043274 :         if (z == ONLY_REM) return gen_0;
     872     3031269 :         if (z) *z = gen_0;
     873     3031269 :         if (sx < 0) sy = -sy;
     874     3031269 :         return stoi(sy);
     875             :       }
     876             :     }
     877   160032534 :     if (z == ONLY_REM) return icopy(x);
     878     8014695 :     if (z) *z = icopy(x);
     879     8014695 :     return gen_0;
     880             :   }
     881             : DIVIDE: /* quotient is non-zero */
     882   644635409 :   av=avma; if (sx<0) sy = -sy;
     883   644635409 :   if (ly==3)
     884             :   {
     885   333055655 :     ulong lq = lx;
     886             :     ulong si;
     887   333055655 :     q  = cgeti(lq);
     888   332875068 :     si = mpn_divrem_1(LIMBS(q), 0, LIMBS(x), NLIMBS(x), y[2]);
     889   333095564 :     if (q[lq - 1] == 0) lq--;
     890   333095564 :     if (z == ONLY_REM)
     891             :     {
     892   292466107 :       avma=av; if (!si) return gen_0;
     893   266472666 :       r=cgeti(3);
     894   266321026 :       r[1] = evalsigne(sx) | evallgefint(3);
     895   266321026 :       r[2] = si; return r;
     896             :     }
     897    40629457 :     q[1] = evalsigne(sy) | evallgefint(lq);
     898    40629457 :     if (!z) return q;
     899    36711743 :     if (!si) { *z=gen_0; return q; }
     900    16384426 :     r=cgeti(3);
     901    16384413 :     r[1] = evalsigne(sx) | evallgefint(3);
     902    16384413 :     r[2] = si; *z=r; return q;
     903             :   }
     904   311579754 :   if (z == ONLY_REM)
     905             :   {
     906   303012262 :     ulong lr = lgefint(y);
     907   303012262 :     ulong lq = lgefint(x)-lgefint(y)+3;
     908   303012262 :     GEN r = cgeti(lr);
     909   300529860 :     GEN q = cgeti(lq);
     910   299495432 :     mpn_tdiv_qr(LIMBS(q), LIMBS(r),0, LIMBS(x), NLIMBS(x), LIMBS(y), NLIMBS(y));
     911   303738490 :     if (!r[lr - 1])
     912             :     {
     913    84719154 :       while(lr>2 && !r[lr - 1]) lr--;
     914    84719154 :       if (lr == 2) {avma=av; return gen_0;} /* exact division */
     915             :     }
     916   299877066 :     r[1] = evalsigne(sx) | evallgefint(lr);
     917   299877066 :     avma = (pari_sp) r; return r;
     918             :   }
     919             :   else
     920             :   {
     921     8567492 :     ulong lq = lgefint(x)-lgefint(y)+3;
     922     8567492 :     ulong lr = lgefint(y);
     923     8567492 :     GEN q = cgeti(lq);
     924     8567492 :     GEN r = cgeti(lr);
     925     8567492 :     mpn_tdiv_qr(LIMBS(q), LIMBS(r),0, LIMBS(x), NLIMBS(x), LIMBS(y), NLIMBS(y));
     926     8567492 :     if (q[lq - 1] == 0) lq--;
     927     8567492 :     q[1] = evalsigne(sy) | evallgefint(lq);
     928     8567492 :     if (!z) { avma = (pari_sp)q; return q; }
     929     7966579 :     if (!r[lr - 1])
     930             :     {
     931     3386950 :       while(lr>2 && !r[lr - 1]) lr--;
     932     3386950 :       if (lr == 2) {avma=(pari_sp) q; *z=gen_0; return q;} /* exact division */
     933             :     }
     934     5504672 :     r[1] = evalsigne(sx) | evallgefint(lr);
     935     5504672 :     avma = (pari_sp) r; *z = r; return q;
     936             :   }
     937             : }
     938             : 
     939             : /* Montgomery reduction.
     940             :  * N has k words, assume T >= 0 has less than 2k.
     941             :  * Return res := T / B^k mod N, where B = 2^BIL
     942             :  * such that 0 <= res < T/B^k + N  and  res has less than k words */
     943             : GEN
     944     6607940 : red_montgomery(GEN T, GEN N, ulong inv)
     945             : {
     946             :   pari_sp av;
     947             :   GEN Te, Td, Ne, Nd, scratch;
     948     6607940 :   ulong i, j, m, t, d, k = NLIMBS(N);
     949             :   int carry;
     950             :   LOCAL_HIREMAINDER;
     951             :   LOCAL_OVERFLOW;
     952             : 
     953     6607940 :   if (k == 0) return gen_0;
     954     6607940 :   d = NLIMBS(T); /* <= 2*k */
     955     6607940 :   if (d == 0) return gen_0;
     956             : #ifdef DEBUG
     957             :   if (d > 2*k) pari_err_BUG("red_montgomery");
     958             : #endif
     959     6607931 :   if (k == 1)
     960             :   { /* as below, special cased for efficiency */
     961           0 :     ulong n = uel(N,2);
     962           0 :     if (d == 1) {
     963           0 :       hiremainder = uel(T,2);
     964           0 :       m = hiremainder * inv;
     965           0 :       (void)addmul(m, n); /* t + m*n = 0 */
     966           0 :       return utoi(hiremainder);
     967             :     } else { /* d = 2 */
     968           0 :       hiremainder = uel(T,2);
     969           0 :       m = hiremainder * inv;
     970           0 :       (void)addmul(m, n); /* t + m*n = 0 */
     971           0 :       t = addll(hiremainder, uel(T,3));
     972           0 :       if (overflow) t -= n; /* t > n doesn't fit in 1 word */
     973           0 :       return utoi(t);
     974             :     }
     975             :   }
     976             :   /* assume k >= 2 */
     977     6607931 :   av = avma; scratch = new_chunk(k<<1); /* >= k + 2: result fits */
     978             : 
     979             :   /* copy T to scratch space (pad with zeroes to 2k words) */
     980     6603735 :   Td = scratch;
     981     6603735 :   Te = T + 2;
     982     6603735 :   for (i=0; i < d     ; i++) *Td++ = *Te++;
     983     6603735 :   for (   ; i < (k<<1); i++) *Td++ = 0;
     984             : 
     985     6603735 :   Te = scratch - 1; /* 1 beyond end of current T mantissa (in scratch) */
     986     6603735 :   Ne = N + 1;       /* 1 beyond end of N mantissa */
     987             : 
     988     6603735 :   carry = 0;
     989    27235911 :   for (i=0; i<k; i++) /* set T := T/B nod N, k times */
     990             :   {
     991    20632176 :     Td = Te; /* one beyond end of (new) T mantissa */
     992    20632176 :     Nd = Ne;
     993    20632176 :     hiremainder = *++Td;
     994    20632176 :     m = hiremainder * inv; /* solve T + m N = O(B) */
     995             : 
     996             :     /* set T := (T + mN) / B */
     997    20632176 :     Te = Td;
     998    20632176 :     (void)addmul(m, *++Nd); /* = 0 */
     999    81210241 :     for (j=1; j<k; j++)
    1000             :     {
    1001    60578065 :       t = addll(addmul(m, *++Nd), *++Td);
    1002    60578065 :       *Td = t;
    1003    60578065 :       hiremainder += overflow;
    1004             :     }
    1005    20632176 :     t = addll(hiremainder, *++Td); *Td = t + carry;
    1006    20632176 :     carry = (overflow || (carry && *Td == 0));
    1007             :   }
    1008     6603735 :   if (carry)
    1009             :   { /* Td > N overflows (k+1 words), set Td := Td - N */
    1010       12100 :     Td = Te;
    1011       12100 :     Nd = Ne;
    1012       12100 :     t = subll(*++Td, *++Nd); *Td = t;
    1013       12100 :     while (Td < (GEN)av) { t = subllx(*++Td, *++Nd); *Td = t; }
    1014             :   }
    1015             : 
    1016             :   /* copy result */
    1017     6603735 :   Td = (GEN)av - 1; /* *Td = high word of final result */
    1018     6603735 :   while (*Td == 0 && Te < Td) Td--; /* strip leading 0s */
    1019     6603735 :   k = Td - Te; if (!k) { avma = av; return gen_0; }
    1020     6603735 :   Td = (GEN)av - k; /* will write mantissa there */
    1021     6603735 :   (void)memmove(Td, Te+1, k*sizeof(long));
    1022     6603735 :   Td -= 2;
    1023     6603735 :   Td[0] = evaltyp(t_INT) | evallg(k+2);
    1024     6600662 :   Td[1] = evalsigne(1) | evallgefint(k+2);
    1025             : #ifdef DEBUG
    1026             : {
    1027             :   long l = NLIMBS(N), s = BITS_IN_LONG*l;
    1028             :   GEN R = int2n(s);
    1029             :   GEN res = remii(mulii(T, Fp_inv(R, N)), N);
    1030             :   if (k > lgefint(N)
    1031             :     || !equalii(remii(Td,N),res)
    1032             :     || cmpii(Td, addii(shifti(T, -s), N)) >= 0) pari_err_BUG("red_montgomery");
    1033             : }
    1034             : #endif
    1035     6600662 :   avma = (pari_sp)Td; return Td;
    1036             : }
    1037             : 
    1038             : /* EXACT INTEGER DIVISION */
    1039             : 
    1040             : #if 1 && !defined(_WIN64) /* use undocumented GMP interface */
    1041             :                           /* mpz_divexact_ui is not LLP64 friendly */
    1042             : static void
    1043   242413106 : GEN2mpz(mpz_t X, GEN x)
    1044             : {
    1045   242413106 :   long l = lgefint(x)-2;
    1046   242413106 :   X->_mp_alloc = l;
    1047   242413106 :   X->_mp_size = signe(x) > 0? l: -l;
    1048   242413106 :   X->_mp_d = LIMBS(x);
    1049   242413106 : }
    1050             : static void
    1051   225093601 : mpz2GEN(GEN z, mpz_t Z)
    1052             : {
    1053   225093601 :   long l = Z->_mp_size;
    1054   225093601 :   z[1] = evalsigne(l > 0? 1: -1) | evallgefint(labs(l)+2);
    1055   225093601 : }
    1056             : 
    1057             : /* assume y != 0 and the division is exact */
    1058             : GEN
    1059   315447571 : diviuexact(GEN x, ulong y)
    1060             : {
    1061   315447571 :   if (!signe(x)) return gen_0;
    1062             :   {
    1063   207762814 :     long l = lgefint(x);
    1064             :     mpz_t X, Z;
    1065   207762814 :     GEN z = cgeti(l);
    1066   207749643 :     GEN2mpz(X, x);
    1067   207750928 :     Z->_mp_alloc = l-2;
    1068   207750928 :     Z->_mp_size  = l-2;
    1069   207750928 :     Z->_mp_d = LIMBS(z);
    1070   207750928 :     mpz_divexact_ui(Z, X, y);
    1071   207762802 :     mpz2GEN(z, Z); return z;
    1072             :   }
    1073             : }
    1074             : 
    1075             : /* Find z such that x=y*z, knowing that y | x (unchecked) */
    1076             : GEN
    1077   322073854 : diviiexact(GEN x, GEN y)
    1078             : {
    1079   322073854 :   if (!signe(y)) pari_err_INV("diviiexact",y);
    1080   322073970 :   if (lgefint(y) == 3)
    1081             :   {
    1082   289064533 :     GEN z = diviuexact(x, y[2]);
    1083   289063628 :     if (signe(y) < 0) togglesign(z);
    1084   289064750 :     return z;
    1085             :   }
    1086    33009437 :   if (!signe(x)) return gen_0;
    1087             :   {
    1088    17329698 :     long l = lgefint(x);
    1089             :     mpz_t X, Y, Z;
    1090    17329698 :     GEN z = cgeti(l);
    1091    17329698 :     GEN2mpz(X, x);
    1092    17329698 :     GEN2mpz(Y, y);
    1093    17329698 :     Z->_mp_alloc = l-2;
    1094    17329698 :     Z->_mp_size  = l-2;
    1095    17329698 :     Z->_mp_d = LIMBS(z);
    1096    17329698 :     mpz_divexact(Z, X, Y);
    1097    17329698 :     mpz2GEN(z, Z); return z;
    1098             :   }
    1099             : }
    1100             : #else
    1101             : /* assume y != 0 and the division is exact */
    1102             : GEN
    1103             : diviuexact(GEN x, ulong y)
    1104             : {
    1105             :   /*TODO: implement true exact division.*/
    1106             :   return divii(x,utoi(y));
    1107             : }
    1108             : 
    1109             : /* Find z such that x=y*z, knowing that y | x (unchecked)
    1110             :  * Method: y0 z0 = x0 mod B = 2^BITS_IN_LONG ==> z0 = 1/y0 mod B.
    1111             :  *    Set x := (x - z0 y) / B, updating only relevant words, and repeat */
    1112             : GEN
    1113             : diviiexact(GEN x, GEN y)
    1114             : {
    1115             :   /*TODO: use mpn_bdivmod instead*/
    1116             :   if (!signe(y)) pari_err_INV("diviiexact",y);
    1117             :   return divii(x,y);
    1118             : }
    1119             : #endif
    1120             : 
    1121             : /* assume yz != and yz | x */
    1122             : GEN
    1123      510975 : diviuuexact(GEN x, ulong y, ulong z)
    1124             : {
    1125             :   long tmp[4];
    1126             :   ulong t;
    1127             :   LOCAL_HIREMAINDER;
    1128      510975 :   t = mulll(y, z);
    1129      510975 :   if (!hiremainder) return diviuexact(x, t);
    1130           0 :   tmp[0] = evaltyp(t_INT)|_evallg(4);
    1131           0 :   tmp[1] = evalsigne(1)|evallgefint(4);
    1132           0 :   tmp[2] = t;
    1133           0 :   tmp[3] = hiremainder;
    1134           0 :   return diviiexact(x, tmp);
    1135             : }
    1136             : 
    1137             : 
    1138             : /********************************************************************/
    1139             : /**                                                                **/
    1140             : /**               INTEGER MULTIPLICATION                           **/
    1141             : /**                                                                **/
    1142             : /********************************************************************/
    1143             : 
    1144             : /* nx >= ny = num. of digits of x, y (not GEN, see mulii) */
    1145             : GEN
    1146  1874287133 : muliispec(GEN x, GEN y, long nx, long ny)
    1147             : {
    1148             :   GEN zd;
    1149             :   long lz;
    1150             :   ulong hi;
    1151             : 
    1152  1874287133 :   if (nx < ny) swapspec(x,y, nx,ny);
    1153  1874287133 :   if (!ny) return gen_0;
    1154  1874287133 :   if (ny == 1) return muluispec((ulong)*y, x, nx);
    1155             : 
    1156   498661952 :   lz = nx+ny+2;
    1157   498661952 :   zd = cgeti(lz);
    1158   497336028 :   hi = mpn_mul(LIMBS(zd), (mp_limb_t *)x, nx, (mp_limb_t *)y, ny);
    1159   499796872 :   if (!hi) lz--;
    1160             :   /*else zd[lz-1]=hi; GH tell me it is not necessary.*/
    1161             : 
    1162   499796872 :   zd[1] = evalsigne(1) | evallgefint(lz);
    1163   499796872 :   return zd;
    1164             : }
    1165             : GEN
    1166      533183 : muluui(ulong x, ulong y, GEN z)
    1167             : {
    1168      533183 :   long t, s = signe(z);
    1169             :   GEN r;
    1170             :   LOCAL_HIREMAINDER;
    1171             : 
    1172      533183 :   if (!x || !y || !signe(z)) return gen_0;
    1173      532787 :   t = mulll(x,y);
    1174      532787 :   if (!hiremainder)
    1175      532787 :     r = muluispec(t, z+2, lgefint(z)-2);
    1176             :   else
    1177             :   {
    1178             :     long tmp[2];
    1179           0 :     tmp[1] = hiremainder;
    1180           0 :     tmp[0] = t;
    1181           0 :     r = muliispec(z+2,tmp, lgefint(z)-2, 2);
    1182             :   }
    1183      532787 :   setsigne(r,s); return r;
    1184             : }
    1185             : 
    1186             : GEN
    1187   648443878 : sqrispec(GEN x, long nx)
    1188             : {
    1189             :   GEN zd;
    1190             :   long lz;
    1191             : 
    1192   648443878 :   if (!nx) return gen_0;
    1193   162055528 :   if (nx==1) return sqru(*x);
    1194             : 
    1195    99277719 :   lz = (nx<<1)+2;
    1196    99277719 :   zd = cgeti(lz);
    1197             : #ifdef mpn_sqr
    1198    98940138 :   mpn_sqr(LIMBS(zd), (mp_limb_t *)x, nx);
    1199             : #else
    1200             :   mpn_mul_n(LIMBS(zd), (mp_limb_t *)x, (mp_limb_t *)x, nx);
    1201             : #endif
    1202    99798958 :   if (zd[lz-1]==0) lz--;
    1203             : 
    1204    99798958 :   zd[1] = evalsigne(1) | evallgefint(lz);
    1205    99798958 :   return zd;
    1206             : }
    1207             : 
    1208             : INLINE GEN
    1209     4272539 : sqrispec_mirror(GEN x, long nx)
    1210             : {
    1211     4272539 :   GEN cx=new_chunk(nx);
    1212             :   GEN z;
    1213     4272539 :   xmpn_mirrorcopy((mp_limb_t *)cx,(mp_limb_t *)x,nx);
    1214     4272539 :   z=sqrispec(cx, nx);
    1215     4272539 :   xmpn_mirror(LIMBS(z), NLIMBS(z));
    1216     4272539 :   return z;
    1217             : }
    1218             : 
    1219             : /* leaves garbage on the stack. */
    1220             : INLINE GEN
    1221    61578064 : muliispec_mirror(GEN x, GEN y, long nx, long ny)
    1222             : {
    1223             :   GEN cx, cy, z;
    1224    61578064 :   long s = 0;
    1225    61578064 :   while (nx && x[nx-1]==0) { nx--; s++; }
    1226    61578064 :   while (ny && y[ny-1]==0) { ny--; s++; }
    1227    61578064 :   cx=new_chunk(nx); cy=new_chunk(ny);
    1228    61578064 :   xmpn_mirrorcopy((mp_limb_t *)cx,(mp_limb_t *)x,nx);
    1229    61578064 :   xmpn_mirrorcopy((mp_limb_t *)cy,(mp_limb_t *)y,ny);
    1230    61578064 :   z =  nx>=ny ? muliispec(cx, cy, nx, ny): muliispec(cy, cx, ny, nx);
    1231    61578064 :   if (s)
    1232             :   {
    1233     7703687 :     long i, lz = lgefint(z) + s;
    1234     7703687 :     (void)new_chunk(s);
    1235     7703687 :     z -= s;
    1236     7703687 :     for (i=0; i<s; i++) z[2+i]=0;
    1237     7703687 :     z[1] = evalsigne(1) | evallgefint(lz);
    1238     7703687 :     z[0] = evaltyp(t_INT) | evallg(lz);
    1239             :   }
    1240    61578064 :   xmpn_mirror(LIMBS(z), NLIMBS(z));
    1241    61578064 :   return z;
    1242             : }
    1243             : 
    1244             : /* x % (2^n), assuming n >= 0 */
    1245             : GEN
    1246    11486161 : remi2n(GEN x, long n)
    1247             : {
    1248             :   ulong hi;
    1249    11486161 :   long l, k, lx, ly, sx = signe(x);
    1250             :   GEN z, xd, zd;
    1251             : 
    1252    11486161 :   if (!sx || !n) return gen_0;
    1253             : 
    1254    11337557 :   k = dvmdsBIL(n, &l);
    1255    11356169 :   lx = lgefint(x);
    1256    11356169 :   if (lx < k+3) return icopy(x);
    1257             : 
    1258    11266702 :   xd = x + (2 + k);
    1259             :   /* x = |k|...|1|#|... : copy the last l bits of # and the first k words
    1260             :    *              ^--- initial xd  */
    1261    11266702 :   hi = ((ulong)*xd) & ((1UL<<l)-1); /* last l bits of # = top bits of result */
    1262    11266702 :   if (!hi)
    1263             :   { /* strip leading zeroes from result */
    1264     1098396 :     xd--; while (k && !*xd) { k--; xd--; }
    1265     1098396 :     if (!k) return gen_0;
    1266      483070 :     ly = k+2;
    1267             :   }
    1268             :   else
    1269    10168306 :     ly = k+3;
    1270             : 
    1271    10651376 :   zd = z = cgeti(ly);
    1272    10517826 :   *++zd = evalsigne(sx) | evallgefint(ly);
    1273    10517826 :   xd = x+1; for ( ;k; k--) *++zd = *++xd;
    1274    10517826 :   if (hi) *++zd = hi;
    1275    10517826 :   return z;
    1276             : }
    1277             : 
    1278             : /********************************************************************/
    1279             : /**                                                                **/
    1280             : /**                      INTEGER SQUARE ROOT                       **/
    1281             : /**                                                                **/
    1282             : /********************************************************************/
    1283             : 
    1284             : /* Return S (and set R) s.t S^2 + R = N, 0 <= R <= 2S.
    1285             :  * As for dvmdii, R is last on stack and guaranteed to be gen_0 in case the
    1286             :  * remainder is 0. R = NULL is allowed. */
    1287             : GEN
    1288     2605723 : sqrtremi(GEN a, GEN *r)
    1289             : {
    1290     2605723 :   long l, na = NLIMBS(a);
    1291             :   mp_size_t nr;
    1292             :   GEN S;
    1293     2605723 :   if (!na) {
    1294           0 :     if (r) *r = gen_0;
    1295           0 :     return gen_0;
    1296             :   }
    1297     2605723 :   l = (na + 5) >> 1; /* 2 + ceil(na/2) */
    1298     2605723 :   S = cgetipos(l);
    1299     2605723 :   if (r) {
    1300        3744 :     GEN R = cgeti(2 + na);
    1301        3744 :     nr = mpn_sqrtrem(LIMBS(S), LIMBS(R), LIMBS(a), na);
    1302        3744 :     if (nr) R[1] = evalsigne(1) | evallgefint(nr+2);
    1303        3363 :     else    { avma = (pari_sp)S; R = gen_0; }
    1304        3744 :     *r = R;
    1305             :   }
    1306             :   else
    1307     2601979 :     (void)mpn_sqrtrem(LIMBS(S), NULL, LIMBS(a), na);
    1308     2605723 :   return S;
    1309             : }
    1310             : 
    1311             : /* compute sqrt(|a|), assuming a != 0 */
    1312             : GEN
    1313    21056069 : sqrtr_abs(GEN a)
    1314             : {
    1315             :   GEN res;
    1316             :   mp_limb_t *b, *c;
    1317    21056069 :   long l = RNLIMBS(a), e = expo(a), er = e>>1;
    1318             :   long n;
    1319    21056069 :   res = cgetr(2 + l);
    1320    21056069 :   res[1] = evalsigne(1) | evalexpo(er);
    1321    21056069 :   if (e&1)
    1322             :   {
    1323    10771104 :     b = (mp_limb_t *) new_chunk(l<<1);
    1324    10771104 :     xmpn_zero(b,l);
    1325    10771104 :     xmpn_mirrorcopy(b+l, RLIMBS(a), l);
    1326    10771104 :     c = (mp_limb_t *) new_chunk(l);
    1327    10771104 :     n = mpn_sqrtrem(c,b,b,l<<1); /* c <- sqrt; b <- rem */
    1328    10771104 :     if (n>l || (n==l && mpn_cmp(b,c,l) > 0)) mpn_add_1(c,c,l,1);
    1329             :   }
    1330             :   else
    1331             :   {
    1332             :     ulong u;
    1333    10284965 :     b = (mp_limb_t *) mantissa2nr(a,-1);
    1334    10284965 :     b[1] = uel(a,l+1)<<(BITS_IN_LONG-1);
    1335    10284965 :     b = (mp_limb_t *) new_chunk(l);
    1336    10284965 :     xmpn_zero(b,l+1); /* overwrites the former b[0] */
    1337    10284965 :     c = (mp_limb_t *) new_chunk(l + 1);
    1338    10284965 :     n = mpn_sqrtrem(c,b,b,(l<<1)+2); /* c <- sqrt; b <- rem */
    1339    10284965 :     u = (ulong)*c++;
    1340    10284965 :     if ( u&HIGHBIT || (u == ~HIGHBIT &&
    1341           0 :              (n>l || (n==l && mpn_cmp(b,c,l)>0))))
    1342     5000022 :       mpn_add_1(c,c,l,1);
    1343             :   }
    1344    21056069 :   xmpn_mirrorcopy(RLIMBS(res),c,l);
    1345    21056069 :   avma = (pari_sp)res; return res;
    1346             : }
    1347             : 
    1348             : /* Normalize a non-negative integer */
    1349             : GEN
    1350   150449220 : int_normalize(GEN x, long known_zero_words)
    1351             : {
    1352   150449220 :   long i =  lgefint(x) - 1 - known_zero_words;
    1353  1269913444 :   for ( ; i > 1; i--)
    1354  1206895518 :     if (x[i]) { setlgefint(x, i+1); return x; }
    1355    63017926 :   x[1] = evalsigne(0) | evallgefint(2); return x;
    1356             : }
    1357             : 
    1358             : /********************************************************************
    1359             :  **                                                                **
    1360             :  **                           Base Conversion                      **
    1361             :  **                                                                **
    1362             :  ********************************************************************/
    1363             : 
    1364             : ulong *
    1365      482017 : convi(GEN x, long *l)
    1366             : {
    1367      482017 :   long n = nchar2nlong(2 + (long)(NLIMBS(x) * (BITS_IN_LONG * LOG10_2)));
    1368      482017 :   GEN str = cgetg(n+1, t_VECSMALL);
    1369      482017 :   unsigned char *res = (unsigned char*) GSTR(str);
    1370      482017 :   long llz = mpn_get_str(res, 10, LIMBS(icopy(x)), NLIMBS(x));
    1371             :   long lz;
    1372             :   ulong *z;
    1373             :   long i, j;
    1374             :   unsigned char *t;
    1375      482017 :   while (!*res) {res++; llz--;} /*Strip leading zeros*/
    1376      482017 :   lz  = (8+llz)/9;
    1377      482017 :   z = (ulong*)new_chunk(1+lz);
    1378      482017 :   t=res+llz+9;
    1379     1209851 :   for(i=0;i<llz-8;i+=9)
    1380             :   {
    1381             :     ulong s;
    1382      727834 :     t-=18;
    1383      727834 :     s=*t++;
    1384     6550506 :     for (j=1; j<9;j++)
    1385     5822672 :       s=10*s+*t++;
    1386      727834 :     *z++=s;
    1387             :   }
    1388      482017 :   if (i<llz)
    1389             :   {
    1390      466908 :     unsigned char *t = res;
    1391      466908 :     ulong s=*t++;
    1392     1976169 :     for (j=i+1; j<llz;j++)
    1393     1509261 :       s=10*s+*t++;
    1394      466908 :     *z++=s;
    1395             :   }
    1396      482017 :   *l = lz;
    1397      482017 :   return z;
    1398             : }

Generated by: LCOV version 1.11