Code coverage tests

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

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

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

LCOV - code coverage report
Current view: top level - kernel/none - level1.h (source / functions) Hit Total Coverage
Test: PARI/GP v2.14.0 lcov report (development 26712-590d837a1c) Lines: 588 750 78.4 %
Date: 2021-06-22 07:13:04 Functions: 211 278 75.9 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : #line 2 "../src/kernel/none/level1.h"
       2             : /* Copyright (C) 2000  The PARI group.
       3             : 
       4             : This file is part of the PARI/GP package.
       5             : 
       6             : PARI/GP is free software; you can redistribute it and/or modify it under the
       7             : terms of the GNU General Public License as published by the Free Software
       8             : Foundation; either version 2 of the License, or (at your option) any later
       9             : version. It is distributed in the hope that it will be useful, but WITHOUT
      10             : ANY WARRANTY WHATSOEVER.
      11             : 
      12             : Check the License for details. You should have received a copy of it, along
      13             : with the package; see the file 'COPYING'. If not, write to the Free Software
      14             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      15             : 
      16             : /* This file defines "level 1" kernel functions.
      17             :  * These functions can be inline; they are also defined externally in
      18             :  * mpinl.c, which includes this file and never needs to be changed */
      19             : 
      20             : INLINE long
      21 74222229856 : evallg(long x)
      22             : {
      23 74222229856 :   if (x & ~LGBITS) pari_err_OVERFLOW("lg()");
      24 74225527008 :   return _evallg(x);
      25             : }
      26             : INLINE long
      27    54655891 : evalvalp(long x)
      28             : {
      29    54655891 :   long v = _evalvalp(x);
      30    54655891 :   if (v & ~VALPBITS) pari_err_OVERFLOW("valp()");
      31    54655891 :   return v;
      32             : }
      33             : INLINE long
      34 10948100941 : evalexpo(long x)
      35             : {
      36 10948100941 :   long v = _evalexpo(x);
      37 10948100941 :   if (v & ~EXPOBITS) pari_err_OVERFLOW("expo()");
      38 10942600359 :   return v;
      39             : }
      40             : INLINE long
      41    21870171 : evalprecp(long x)
      42             : {
      43    21870171 :   long v = _evalprecp(x);
      44    21870171 :   if (x & ~((1UL<<(BITS_IN_LONG-VALPnumBITS))-1)) pari_err_OVERFLOW("precp()");
      45    21870174 :   return v;
      46             : }
      47             : 
      48             : INLINE int
      49   161875164 : varncmp(long x, long y)
      50             : {
      51   161875164 :   if (varpriority[x] < varpriority[y]) return  1;
      52   124736368 :   if (varpriority[x] > varpriority[y]) return -1;
      53    68435670 :   return 0;
      54             : }
      55             : INLINE long
      56           0 : varnmin(long x, long y)
      57           0 : { return (varpriority[x] <= varpriority[y])? x: y; }
      58             : INLINE long
      59         203 : varnmax(long x, long y)
      60         203 : { return (varpriority[x] >= varpriority[y])? x: y; }
      61             : 
      62             : /* Inhibit some area gerepile-wise: declare it to be a non recursive
      63             :  * type, of length l. Thus gerepile won't inspect the zone, just copy it.
      64             :  * For the following situation:
      65             :  *   z = cgetg(t,a); av = avma; garbage(); ltop = avma;
      66             :  *   for (i=1; i<HUGE; i++) gel(z,i) = blah();
      67             :  *   stackdummy(av,ltop);
      68             :  * loses (av-ltop) words but save a costly gerepile. */
      69             : INLINE void
      70  2817650241 : stackdummy(pari_sp av, pari_sp ltop) {
      71  2817650241 :   long l = ((GEN)av) - ((GEN)ltop);
      72  2817650241 :   if (l > 0) {
      73   901207443 :     GEN z = (GEN)ltop;
      74   901207443 :     z[0] = evaltyp(t_VECSMALL) | evallg(l);
      75             : #ifdef DEBUG
      76             :     { long i; for (i = 1; i < l; i++) z[i] = 0; }
      77             : #endif
      78             :   }
      79  2817312369 : }
      80             : INLINE void
      81    71922006 : fixlg(GEN x, long ly) {
      82    71922006 :   long lx = lg(x), l = lx - ly;
      83    71922006 :   if (l > 0)
      84             :   { /* stackdummy(x+lx, x+ly) */
      85    45405587 :     GEN z = x + ly;
      86    45405587 :     z[0] = evaltyp(t_VECSMALL) | evallg(l);
      87    45405648 :     setlg(x, ly);
      88             : #ifdef DEBUG
      89             :     { long i; for (i = 1; i < l; i++) z[i] = 0; }
      90             : #endif
      91             :   }
      92    71922140 : }
      93             : /* update lg(z) before affrr(y, z)  [ to cater for precision loss ]*/
      94             : INLINE void
      95    26524262 : affrr_fixlg(GEN y, GEN z) { fixlg(z, lg(y)); affrr(y, z); }
      96             : 
      97             : /*******************************************************************/
      98             : /*                                                                 */
      99             : /*                       ALLOCATE ON STACK                         */
     100             : /*                                                                 */
     101             : /*******************************************************************/
     102             : INLINE void
     103 98417541698 : set_avma(ulong av) { avma = av; }
     104             : 
     105             : INLINE double
     106   296028476 : gc_double(pari_sp av, double d) { set_avma(av); return d; }
     107             : INLINE long
     108   137859507 : gc_long(pari_sp av, long s) { set_avma(av); return s; }
     109             : INLINE ulong
     110    20157842 : gc_ulong(pari_sp av, ulong s) { set_avma(av); return s; }
     111             : INLINE int
     112    33167251 : gc_bool(pari_sp av, int s) { set_avma(av); return s; }
     113             : INLINE int
     114     1177550 : gc_int(pari_sp av, int s) { set_avma(av); return s; }
     115             : INLINE GEN
     116     6050012 : gc_NULL(pari_sp av) { set_avma(av); return NULL; }
     117             : INLINE GEN
     118  9849174320 : gc_const(pari_sp av, GEN x) { set_avma(av); return x; }
     119             : 
     120             : INLINE GEN
     121 72192595909 : new_chunk(size_t x) /* x is a number of longs */
     122             : {
     123 72192595909 :   GEN z = ((GEN) avma) - x;
     124             :   CHECK_CTRLC
     125 72192595909 :   if (x > (avma-pari_mainstack->bot) / sizeof(long))
     126           7 :     new_chunk_resize(x);
     127 72192595902 :   set_avma((pari_sp)z);
     128             : #ifdef MEMSTEP
     129             :   if (DEBUGMEM>1 && pari_mainstack->memused != DISABLE_MEMUSED) {
     130             :     long d = (long)pari_mainstack->memused - (long)z;
     131             :     if (labs(d) > 4*MEMSTEP)
     132             :     {
     133             :       pari_mainstack->memused = (pari_sp)z;
     134             :       err_printf("...%4.0lf Mbytes used\n",
     135             :                 (pari_mainstack->top-pari_mainstack->memused)/1048576.);
     136             :     }
     137             :   }
     138             : #endif
     139 72167302014 :   return z;
     140             : }
     141             : 
     142             : INLINE char *
     143     8032737 : stack_malloc(size_t N)
     144             : {
     145     8032737 :   long n = nchar2nlong(N);
     146     8032732 :   return (char*)new_chunk(n);
     147             : }
     148             : 
     149             : INLINE char *
     150    81678116 : stack_malloc_align(size_t N, long k)
     151             : {
     152    81678116 :   ulong d = ((ulong)avma) % k, e = ((ulong)N) % k;
     153    81678116 :   if (d) (void)new_chunk(d/sizeof(long));
     154    81678114 :   if (e) N += k-e;
     155    81678114 :   return (char*) new_chunk(nchar2nlong(N));
     156             : }
     157             : 
     158             : INLINE char *
     159       83564 : stack_calloc(size_t N)
     160             : {
     161       83564 :   char *p = stack_malloc(N);
     162       83562 :   memset(p, 0, N); return p;
     163             : }
     164             : 
     165             : INLINE char *
     166         489 : stack_calloc_align(size_t N, long k)
     167             : {
     168         489 :   ulong d = ((ulong)avma) % k, e = ((ulong)N) % k;
     169         489 :   if (d) (void)new_chunk(d/sizeof(long));
     170         489 :   if (e) N += k-e;
     171         489 :   return stack_calloc(N);
     172             : }
     173             : 
     174             : /* cgetg(lg(x), typ(x)), set *lx. Implicit unsetisclone() */
     175             : INLINE GEN
     176   936089844 : cgetg_copy(GEN x, long *plx) {
     177             :   GEN y;
     178   936089844 :   *plx = lg(x); y = new_chunk((size_t)*plx);
     179   936107855 :   y[0] = x[0] & (TYPBITS|LGBITS); return y;
     180             : }
     181             : INLINE GEN
     182      242119 : cgetg_block(long x, long y)
     183             : {
     184      242119 :   GEN z = newblock((size_t)x);
     185      242401 :   z[0] = CLONEBIT | evaltyp(y) | evallg(x);
     186      242280 :   return z;
     187             : }
     188             : INLINE GEN
     189 12360897253 : cgetg(long x, long y)
     190             : {
     191 12360897253 :   GEN z = new_chunk((size_t)x);
     192 12358253891 :   z[0] = evaltyp(y) | evallg(x);
     193 12355226882 :   return z;
     194             : }
     195             : INLINE GEN
     196 21835805348 : cgeti(long x)
     197             : {
     198 21835805348 :   GEN z = new_chunk((size_t)x);
     199 21791999342 :   z[0] = evaltyp(t_INT) | evallg(x);
     200 21775102536 :   return z;
     201             : }
     202             : INLINE GEN
     203 13536964262 : cgetipos(long x)
     204             : {
     205 13536964262 :   GEN z = cgeti(x);
     206 13514834645 :   z[1] = evalsigne(1) | evallgefint(x);
     207 13514834645 :   return z;
     208             : }
     209             : INLINE GEN
     210   304544735 : cgetineg(long x)
     211             : {
     212   304544735 :   GEN z = cgeti(x);
     213   304545200 :   z[1] = evalsigne(-1) | evallgefint(x);
     214   304545200 :   return z;
     215             : }
     216             : INLINE GEN
     217       33604 : cgetr_block(long x)
     218             : {
     219       33604 :   GEN z = newblock((size_t)x);
     220       33612 :   z[0] = CLONEBIT | evaltyp(t_REAL) | evallg(x);
     221       33610 :   return z;
     222             : }
     223             : INLINE GEN
     224 10884478026 : cgetr(long x)
     225             : {
     226 10884478026 :   GEN z = new_chunk((size_t)x);
     227 10868363400 :   z[0] = evaltyp(t_REAL) | evallg(x);
     228 10862545109 :   return z;
     229             : }
     230             : 
     231             : /*******************************************************************/
     232             : /*                                                                 */
     233             : /*                     COPY, NEGATION, ABSOLUTE VALUE              */
     234             : /*                                                                 */
     235             : /*******************************************************************/
     236             : /* cannot do memcpy because sometimes x and y overlap */
     237             : INLINE GEN
     238  3483105787 : leafcopy(GEN x)
     239             : {
     240  3483105787 :   register long lx = lg(x);
     241  3483105787 :   GEN y = new_chunk(lx); /* can't use cgetg_copy, in case x,y overlap */
     242 18373675322 :   while (--lx > 0) y[lx] = x[lx];
     243  3482924797 :   y[0] = x[0] & (TYPBITS|LGBITS); return y;
     244             : }
     245             : INLINE GEN
     246  6581233799 : icopy(GEN x)
     247             : {
     248  6581233799 :   long i = lgefint(x), lx = i;
     249  6581233799 :   GEN y = new_chunk(lx); /* can't use cgeti, in case x,y overlap */
     250 27215866106 :   while (--i > 0) y[i] = x[i];
     251  6575802896 :   y[0] = evaltyp(t_INT) | evallg(lx);
     252  6578944326 :   return y;
     253             : }
     254             : INLINE GEN
     255    60913461 : icopyspec(GEN x, long nx)
     256             : {
     257    60913461 :   long i = nx+2, lx = i;
     258    60913461 :   GEN y = new_chunk(lx); /* can't use cgeti, in case x,y overlap */
     259   896489215 :   x -= 2; while (--i >= 2) y[i] = x[i];
     260    60912867 :   y[1] = evalsigne(1) | evallgefint(lx);
     261    60912867 :   y[0] = evaltyp(t_INT) | evallg(lx);
     262    60912650 :   return y;
     263             : }
     264   594887111 : INLINE GEN rcopy(GEN x) { return leafcopy(x); }
     265         294 : INLINE GEN mpcopy(GEN x) { return leafcopy(x); }
     266             : 
     267             : INLINE GEN
     268   629399698 : mpabs(GEN x) { GEN y = leafcopy(x); setabssign(y); return y; }
     269             : INLINE GEN
     270     2613808 : mpabs_shallow(GEN x) { return signe(x) < 0? mpabs(x): x; }
     271   589125246 : INLINE GEN absi(GEN x) { return mpabs(x); }
     272    54033859 : INLINE GEN absi_shallow(GEN x) { return signe(x) < 0? negi(x): x; }
     273      232946 : INLINE GEN absr(GEN x) { return mpabs(x); }
     274             : 
     275             : INLINE GEN
     276  1444978975 : mpneg(GEN x) { GEN y = leafcopy(x); togglesign(y); return y; }
     277   590388018 : INLINE GEN negi(GEN x) { return mpneg(x); }
     278     1873973 : INLINE GEN negr(GEN x) { return mpneg(x); }
     279             : 
     280             : /* negate in place */
     281             : INLINE void
     282  2208472007 : togglesign(GEN x) { if (x[1] & SIGNBITS) { x[1] ^= HIGHBIT; } }
     283             : INLINE void
     284   672884003 : setabssign(GEN x) { x[1] &= ~HIGHBIT; }
     285             : /* negate in place, except universal constants */
     286             : INLINE void
     287   113146378 : togglesign_safe(GEN *px)
     288             : {
     289   113146378 :   switch(*px - gen_1) /* gen_1, gen_2, gen_m1, gen_m2 */
     290             :   {
     291     2381499 :     case 0: *px = gen_m1; break;
     292           4 :     case 3: *px = gen_m2;  break;
     293      556922 :     case 6: *px = gen_1; break;
     294           0 :     case 9: *px = gen_2;  break;
     295   110207953 :     default: togglesign(*px);
     296             :   }
     297   113149172 : }
     298             : /* setsigne(y, signe(x)) */
     299             : INLINE void
     300           0 : affectsign(GEN x, GEN y)
     301             : {
     302           0 :   y[1] = (x[1] & SIGNBITS) | (y[1] & ~SIGNBITS);
     303           0 : }
     304             : /* copies sign in place, except for universal constants */
     305             : INLINE void
     306     9914275 : affectsign_safe(GEN x, GEN *py)
     307             : {
     308     9914275 :   if (((*py)[1] ^ x[1]) & HIGHBIT) togglesign_safe(py);
     309     9914275 : }
     310             : /*******************************************************************/
     311             : /*                                                                 */
     312             : /*                     GEN -> LONG, LONG -> GEN                    */
     313             : /*                                                                 */
     314             : /*******************************************************************/
     315             : /* assume x != 0, return -x as a t_INT */
     316             : INLINE GEN
     317   301722947 : utoineg(ulong x) { GEN y = cgetineg(3); y[2] = x; return y; }
     318             : /* assume x != 0, return utoi(x) */
     319             : INLINE GEN
     320 12016571630 : utoipos(ulong x) { GEN y = cgetipos(3); y[2] = x; return y; }
     321             : INLINE GEN
     322 10231785421 : utoi(ulong x) { return x? utoipos(x): gen_0; }
     323             : INLINE GEN
     324   453514486 : stoi(long x)
     325             : {
     326   453514486 :   if (!x) return gen_0;
     327   285530061 :   return x > 0? utoipos((ulong)x): utoineg((ulong)-x);
     328             : }
     329             : 
     330             : /* x 2^BIL + y */
     331             : INLINE GEN
     332  7290546747 : uutoi(ulong x, ulong y)
     333             : {
     334             :   GEN z;
     335  7290546747 :   if (!x) return utoi(y);
     336   623874336 :   z = cgetipos(4);
     337   626778714 :   *int_W_lg(z, 1, 4) = x;
     338   626778714 :   *int_W_lg(z, 0, 4) = y; return z;
     339             : }
     340             : /* - (x 2^BIL + y) */
     341             : INLINE GEN
     342      243986 : uutoineg(ulong x, ulong y)
     343             : {
     344             :   GEN z;
     345      243986 :   if (!x) return y? utoineg(y): gen_0;
     346       10365 :   z = cgetineg(4);
     347       10365 :   *int_W_lg(z, 1, 4) = x;
     348       10365 :   *int_W_lg(z, 0, 4) = y; return z;
     349             : }
     350             : 
     351             : INLINE long
     352   408552648 : itos(GEN x)
     353             : {
     354   408552648 :   long s = signe(x);
     355             :   long u;
     356             : 
     357   408552648 :   if (!s) return 0;
     358   394374979 :   u = x[2];
     359   394374979 :   if (lgefint(x) > 3 || u < 0)
     360          21 :     pari_err_OVERFLOW("t_INT-->long assignment");
     361   394373282 :   return (s>0) ? u : -u;
     362             : }
     363             : /* as itos, but return 0 if too large. Cf is_bigint */
     364             : INLINE long
     365    46568581 : itos_or_0(GEN x) {
     366             :   long n;
     367    46568581 :   if (lgefint(x) != 3 || (n = x[2]) & HIGHBIT) return 0;
     368    45445354 :   return signe(x) > 0? n: -n;
     369             : }
     370             : INLINE ulong
     371   125254782 : itou(GEN x)
     372             : {
     373   125254782 :   switch(lgefint(x)) {
     374     9615983 :     case 2: return 0;
     375   115638886 :     case 3: return x[2];
     376           0 :     default:
     377           0 :       pari_err_OVERFLOW("t_INT-->ulong assignment");
     378             :       return 0; /* LCOV_EXCL_LINE */
     379             :   }
     380             : }
     381             : 
     382             : /* as itou, but return 0 if too large. Cf is_bigint */
     383             : INLINE ulong
     384     4375180 : itou_or_0(GEN x) {
     385     4375180 :   if (lgefint(x) != 3) return 0;
     386     4364301 :   return (ulong)x[2];
     387             : }
     388             : 
     389             : INLINE ulong
     390     6619925 : umuluu_or_0(ulong x, ulong y)
     391             : {
     392             :   ulong z;
     393             :   LOCAL_HIREMAINDER;
     394     6619925 :   z = mulll(x, y);
     395     6619925 :   return hiremainder? 0: z;
     396             : }
     397             : /* return x*y if <= n, else 0. Beware overflow */
     398             : INLINE ulong
     399     3170739 : umuluu_le(ulong x, ulong y, ulong n)
     400             : {
     401             :   ulong z;
     402             :   LOCAL_HIREMAINDER;
     403     3170739 :   z = mulll(x, y);
     404     3170739 :   return (hiremainder || z > n)? 0: z;
     405             : }
     406             : 
     407             : INLINE GEN
     408   165988305 : real_0_bit(long bitprec) { GEN x=cgetr(2); x[1]=evalexpo(bitprec); return x; }
     409             : INLINE GEN
     410      562269 : real_0(long prec) { return real_0_bit(-prec2nbits(prec)); }
     411             : INLINE GEN
     412     3183484 : real_1_bit(long bit) { return real_1(nbits2prec(bit)); }
     413             : INLINE GEN
     414    76134349 : real_1(long prec) {
     415    76134349 :   GEN x = cgetr(prec);
     416             :   long i;
     417    76125819 :   x[1] = evalsigne(1) | _evalexpo(0);
     418   391280580 :   x[2] = (long)HIGHBIT; for (i=3; i<prec; i++) x[i] = 0;
     419    76125819 :   return x;
     420             : }
     421             : INLINE GEN
     422         182 : real_m1(long prec) {
     423         182 :   GEN x = cgetr(prec);
     424             :   long i;
     425         182 :   x[1] = evalsigne(-1) | _evalexpo(0);
     426        1096 :   x[2] = (long)HIGHBIT; for (i=3; i<prec; i++) x[i] = 0;
     427         182 :   return x;
     428             : }
     429             : 
     430             : /* 2.^n */
     431             : INLINE GEN
     432      696327 : real2n(long n, long prec) { GEN z = real_1(prec); setexpo(z, n); return z; }
     433             : INLINE GEN
     434           0 : real_m2n(long n, long prec) { GEN z = real_m1(prec); setexpo(z, n); return z; }
     435             : INLINE GEN
     436   193757773 : stor(long s, long prec) { GEN z = cgetr(prec); affsr(s,z); return z; }
     437             : INLINE GEN
     438    10126259 : utor(ulong s, long prec){ GEN z = cgetr(prec); affur(s,z); return z; }
     439             : INLINE GEN
     440   654453556 : itor(GEN x, long prec) { GEN z = cgetr(prec); affir(x,z); return z; }
     441             : INLINE GEN
     442   189750325 : rtor(GEN x, long prec) { GEN z = cgetr(prec); affrr(x,z); return z; }
     443             : 
     444    12848482 : INLINE ulong int_bit(GEN x, long n)
     445             : {
     446    12848482 :   long r, q = dvmdsBIL(n, &r);
     447    12844778 :   return q < lgefint(x)-2?((ulong)*int_W(x,q) >> r) & 1UL:0;
     448             : }
     449             : 
     450             : /*******************************************************************/
     451             : /*                                                                 */
     452             : /*                           COMPARISON                            */
     453             : /*                                                                 */
     454             : /*******************************************************************/
     455             : INLINE int
     456     1160304 : cmpss(long a, long b)
     457     1160304 : { return a>b? 1: (a<b? -1: 0); }
     458             : 
     459             : INLINE int
     460  1386643492 : cmpuu(ulong a, ulong b)
     461  1386643492 : { return a>b? 1: (a<b? -1: 0); }
     462             : 
     463             : INLINE int
     464      891109 : cmpir(GEN x, GEN y)
     465             : {
     466             :   pari_sp av;
     467             :   GEN z;
     468             : 
     469      891109 :   if (!signe(x)) return -signe(y);
     470      466846 :   if (!signe(y))
     471             :   {
     472        2213 :     if (expo(y) >= expi(x)) return 0;
     473        2185 :     return signe(x);
     474             :   }
     475      464633 :   av=avma; z = itor(x, realprec(y)); set_avma(av);
     476      464616 :   return cmprr(z,y); /* cmprr does no memory adjustment */
     477             : }
     478             : INLINE int
     479      437210 : cmpri(GEN x, GEN y) { return -cmpir(y,x); }
     480             : INLINE int
     481      100002 : cmpsr(long x, GEN y)
     482             : {
     483             :   pari_sp av;
     484             :   GEN z;
     485             : 
     486      100002 :   if (!x) return -signe(y);
     487      100002 :   av=avma; z = stor(x, LOWDEFAULTPREC); set_avma(av);
     488      100002 :   return cmprr(z,y);
     489             : }
     490             : INLINE int
     491       40131 : cmprs(GEN x, long y) { return -cmpsr(y,x); }
     492             : /* compare x and y */
     493             : INLINE int
     494     4555424 : cmpui(ulong x, GEN y)
     495             : {
     496             :   ulong p;
     497     4555424 :   if (!x) return -signe(y);
     498     4555424 :   if (signe(y) <= 0) return 1;
     499     4555403 :   if (lgefint(y) > 3) return -1;
     500     2072994 :   p = y[2]; if (p == x) return 0;
     501     2034743 :   return p < x ? 1 : -1;
     502             : }
     503             : INLINE int
     504     4555423 : cmpiu(GEN x, ulong y) { return -cmpui(y,x); }
     505             : /* compare x and |y| */
     506             : INLINE int
     507    39931196 : abscmpui(ulong x, GEN y)
     508             : {
     509    39931196 :   long l = lgefint(y);
     510             :   ulong p;
     511             : 
     512    39931196 :   if (!x) return (l > 2)? -1: 0;
     513    39931182 :   if (l == 2) return 1;
     514    39677160 :   if (l > 3) return -1;
     515    39656368 :   p = y[2]; if (p == x) return 0;
     516    38950674 :   return p < x ? 1 : -1;
     517             : }
     518             : INLINE int
     519    39085541 : abscmpiu(GEN x, ulong y) { return -abscmpui(y,x); }
     520             : INLINE int
     521     3538168 : cmpsi(long x, GEN y)
     522             : {
     523             :   ulong p;
     524             : 
     525     3538168 :   if (!x) return -signe(y);
     526             : 
     527     3537034 :   if (x > 0)
     528             :   {
     529     3536103 :     if (signe(y)<=0) return 1;
     530     3535774 :     if (lgefint(y)>3) return -1;
     531     3522022 :     p = y[2]; if (p == (ulong)x) return 0;
     532     3456657 :     return p < (ulong)x ? 1 : -1;
     533             :   }
     534             : 
     535         931 :   if (signe(y)>=0) return -1;
     536         147 :   if (lgefint(y)>3) return 1;
     537         147 :   p = y[2]; if (p == (ulong)-x) return 0;
     538          42 :   return p < (ulong)(-x) ? -1 : 1;
     539             : }
     540             : INLINE int
     541     3527277 : cmpis(GEN x, long y) { return -cmpsi(y,x); }
     542             : INLINE int
     543     1549164 : mpcmp(GEN x, GEN y)
     544             : {
     545     1549164 :   if (typ(x)==t_INT)
     546       60208 :     return (typ(y)==t_INT) ? cmpii(x,y) : cmpir(x,y);
     547     1488956 :   return (typ(y)==t_INT) ? -cmpir(y,x) : cmprr(x,y);
     548             : }
     549             : 
     550             : /* x == y ? */
     551             : INLINE int
     552     1987560 : equalui(ulong x, GEN y)
     553             : {
     554     1987560 :   if (!x) return !signe(y);
     555     1987560 :   if (signe(y) <= 0 || lgefint(y) != 3) return 0;
     556     1977350 :   return ((ulong)y[2] == (ulong)x);
     557             : }
     558             : /* x == y ? */
     559             : INLINE int
     560      580768 : equalsi(long x, GEN y)
     561             : {
     562      580768 :   if (!x) return !signe(y);
     563      580768 :   if (x > 0)
     564             :   {
     565      578668 :     if (signe(y) <= 0 || lgefint(y) != 3) return 0;
     566      528473 :     return ((ulong)y[2] == (ulong)x);
     567             :   }
     568        2100 :   if (signe(y) >= 0 || lgefint(y) != 3) return 0;
     569        1994 :   return ((ulong)y[2] == (ulong)-x);
     570             : }
     571             : /* x == |y| ? */
     572             : INLINE int
     573    37228480 : absequalui(ulong x, GEN y)
     574             : {
     575    37228480 :   if (!x) return !signe(y);
     576    37228480 :   return (lgefint(y) == 3 && (ulong)y[2] == x);
     577             : }
     578             : INLINE int
     579    35707535 : absequaliu(GEN x, ulong y) { return absequalui(y,x); }
     580             : INLINE int
     581      580538 : equalis(GEN x, long y) { return equalsi(y,x); }
     582             : INLINE int
     583     1987560 : equaliu(GEN x, ulong y) { return equalui(y,x); }
     584             : 
     585             : /* assume x != 0, is |x| == 2^n ? */
     586             : INLINE int
     587      717352 : absrnz_equal2n(GEN x) {
     588      717352 :   if ((ulong)x[2]==HIGHBIT)
     589             :   {
     590       33016 :     long i, lx = lg(x);
     591      103208 :     for (i = 3; i < lx; i++)
     592       77376 :       if (x[i]) return 0;
     593       25832 :     return 1;
     594             :   }
     595      684336 :   return 0;
     596             : }
     597             : /* assume x != 0, is |x| == 1 ? */
     598             : INLINE int
     599     2459854 : absrnz_equal1(GEN x) { return !expo(x) && absrnz_equal2n(x); }
     600             : 
     601             : INLINE long
     602  8273593498 : maxss(long x, long y) { return x>y?x:y; }
     603             : INLINE long
     604  1229051827 : minss(long x, long y) { return x<y?x:y; }
     605             : INLINE long
     606     8264667 : minuu(ulong x, ulong y) { return x<y?x:y; }
     607             : INLINE long
     608     9448608 : maxuu(ulong x, ulong y) { return x>y?x:y; }
     609             : INLINE double
     610     2693184 : maxdd(double x, double y) { return x>y?x:y; }
     611             : INLINE double
     612     2451017 : mindd(double x, double y) { return x<y?x:y; }
     613             : 
     614             : /*******************************************************************/
     615             : /*                                                                 */
     616             : /*                             ADD / SUB                           */
     617             : /*                                                                 */
     618             : /*******************************************************************/
     619             : INLINE GEN
     620       25046 : subuu(ulong x, ulong y)
     621             : {
     622             :   ulong z;
     623             :   LOCAL_OVERFLOW;
     624       25046 :   z = subll(x, y);
     625       25046 :   return overflow? utoineg(-z): utoi(z);
     626             : }
     627             : INLINE GEN
     628  3260486389 : adduu(ulong x, ulong y) { ulong t = x+y; return uutoi((t < x), t); }
     629             : 
     630             : INLINE GEN
     631       25046 : addss(long x, long y)
     632             : {
     633       25046 :   if (!x) return stoi(y);
     634       25046 :   if (!y) return stoi(x);
     635       25046 :   if (x > 0) return y > 0? adduu(x,y): subuu(x, -y);
     636             : 
     637       25046 :   if (y > 0) return subuu(y, -x);
     638             :   else { /* - adduu(-x, -y) */
     639           0 :     ulong t = (-x)+(-y); return uutoineg((t < (ulong)(-x)), t);
     640             :   }
     641             : }
     642       25046 : INLINE GEN subss(long x, long y) { return addss(-y,x); }
     643             : 
     644             : INLINE GEN
     645  6756845315 : subii(GEN x, GEN y)
     646             : {
     647  6756845315 :   if (x==y) return gen_0; /* frequent with x = y = gen_0 */
     648  4999786944 :   return addii_sign(x, signe(x), y, -signe(y));
     649             : }
     650             : INLINE GEN
     651  8605749535 : addii(GEN x, GEN y) { return addii_sign(x, signe(x), y, signe(y)); }
     652             : INLINE GEN
     653  2406187989 : addrr(GEN x, GEN y) { return addrr_sign(x, signe(x), y, signe(y)); }
     654             : INLINE GEN
     655   928197329 : subrr(GEN x, GEN y) { return addrr_sign(x, signe(x), y, -signe(y)); }
     656             : INLINE GEN
     657   338418191 : addir(GEN x, GEN y) { return addir_sign(x, signe(x), y, signe(y)); }
     658             : INLINE GEN
     659     2917436 : subir(GEN x, GEN y) { return addir_sign(x, signe(x), y, -signe(y)); }
     660             : INLINE GEN
     661     5843707 : subri(GEN x, GEN y) { return addir_sign(y, -signe(y), x, signe(x)); }
     662             : INLINE GEN
     663   116522370 : addsi(long x, GEN y) { return addsi_sign(x, y, signe(y)); }
     664             : INLINE GEN
     665    96856522 : addui(ulong x, GEN y) { return addui_sign(x, y, signe(y)); }
     666             : INLINE GEN
     667     6963279 : subsi(long x, GEN y) { return addsi_sign(x, y, -signe(y)); }
     668             : INLINE GEN
     669    62429991 : subui(ulong x, GEN y) { return addui_sign(x, y, -signe(y)); }
     670             : 
     671             : /*******************************************************************/
     672             : /*                                                                 */
     673             : /*                           MOD, REM, DIV                         */
     674             : /*                                                                 */
     675             : /*******************************************************************/
     676   125841420 : INLINE ulong mod2BIL(GEN x) { return *int_LSW(x); }
     677           0 : INLINE long mod64(GEN x) { return mod2BIL(x) & 63; }
     678         259 : INLINE long mod32(GEN x) { return mod2BIL(x) & 31; }
     679      237152 : INLINE long mod16(GEN x) { return mod2BIL(x) & 15; }
     680    52433185 : INLINE long mod8(GEN x)  { return mod2BIL(x) & 7; }
     681     7326438 : INLINE long mod4(GEN x)  { return mod2BIL(x) & 3; }
     682    50494264 : INLINE long mod2(GEN x)  { return mod2BIL(x) & 1; }
     683             : INLINE int
     684    79522683 : mpodd(GEN x) { return signe(x) && mod2(x); }
     685             : /* x mod 2^n, n < BITS_IN_LONG */
     686             : INLINE ulong
     687    22058715 : umodi2n(GEN x, long n)
     688             : {
     689    22058715 :   long s = signe(x);
     690    22058715 :   const ulong _2n = 1UL << n;
     691             :   ulong m;
     692    22058715 :   if (!s) return 0;
     693    22004881 :   m = *int_LSW(x) & (_2n - 1);
     694    22004881 :   if (s < 0 && m) m = _2n - m;
     695    22004881 :   return m;
     696             : }
     697           0 : INLINE ulong Mod64(GEN x){ return umodi2n(x,6); }
     698      167559 : INLINE ulong Mod32(GEN x){ return umodi2n(x,5); }
     699      239722 : INLINE ulong Mod16(GEN x){ return umodi2n(x,4); }
     700     1770832 : INLINE ulong Mod8(GEN x) { return umodi2n(x,3); }
     701    18330282 : INLINE ulong Mod4(GEN x) { return umodi2n(x,2); }
     702     1550073 : INLINE ulong Mod2(GEN x) { return umodi2n(x,1); }
     703             : 
     704             : INLINE GEN
     705    57766583 : truedivii(GEN a,GEN b) { return truedvmdii(a,b,NULL); }
     706             : INLINE GEN
     707      243909 : truedivis(GEN a, long b) { return truedvmdis(a,b,NULL); }
     708             : INLINE GEN
     709     6162526 : truedivsi(long a, GEN b) { return truedvmdsi(a,b,NULL); }
     710             : 
     711             : INLINE GEN
     712    12524791 : divii(GEN a, GEN b) { return dvmdii(a,b,NULL); }
     713             : INLINE GEN
     714  2027649952 : remii(GEN a, GEN b) { return dvmdii(a,b,ONLY_REM); }
     715             : 
     716             : INLINE GEN
     717           0 : divss(long x, long y) { return stoi(x / y); }
     718             : INLINE GEN
     719           0 : modss(long x, long y) { return utoi(smodss(x, y)); }
     720             : INLINE GEN
     721           0 : remss(long x, long y) { return stoi(x % y); }
     722             : INLINE long
     723        6422 : smodss(long x, long y)
     724             : {
     725        6422 :   long r = x%y;
     726        6422 :   return (r >= 0)? r: labs(y) + r;
     727             : }
     728             : INLINE ulong
     729   576549591 : umodsu(long x, ulong y)
     730             : {
     731   576549591 :   return x>=0 ? x%y: Fl_neg((-x)%y, y);
     732             : }
     733             : 
     734             : INLINE long
     735           0 : sdivss_rem(long x, long y, long *r)
     736             : {
     737             :   long q;
     738             :   LOCAL_HIREMAINDER;
     739           0 :   if (!y) pari_err_INV("sdivss_rem",gen_0);
     740           0 :   hiremainder = 0; q = divll((ulong)labs(x),(ulong)labs(y));
     741           0 :   if (x < 0) { hiremainder = -((long)hiremainder); q = -q; }
     742           0 :   if (y < 0) q = -q;
     743           0 :   *r = hiremainder; return q;
     744             : }
     745             : INLINE GEN
     746           0 : divss_rem(long x, long y, long *r) { return stoi(sdivss_rem(x,y,r)); }
     747             : INLINE ulong
     748   137544744 : udivuu_rem(ulong x, ulong y, ulong *r)
     749             : {
     750   137544744 :   if (!y) pari_err_INV("udivuu_rem",gen_0);
     751   137544744 :   *r = x % y; return x / y;
     752             : }
     753             : INLINE ulong
     754      670434 : ceildivuu(ulong a, ulong b)
     755             : {
     756      670434 :   ulong c = a/b;
     757      670434 :   return (a%b)? c+1: c;
     758             : }
     759             : 
     760             : INLINE ulong
     761       11400 : uabsdivui_rem(ulong x, GEN y, ulong *r)
     762             : {
     763       11400 :   long q, s = signe(y);
     764             :   LOCAL_HIREMAINDER;
     765             : 
     766       11400 :   if (!s) pari_err_INV("uabsdivui_rem",gen_0);
     767       11400 :   if (!x || lgefint(y)>3) { *r = x; return 0; }
     768       11092 :   hiremainder=0; q = (long)divll(x, (ulong)y[2]);
     769       11092 :   if (s < 0) q = -q;
     770       11092 :   *r = hiremainder; return q;
     771             : }
     772             : 
     773             : /* assume d != 0 and |n| / d can be represented as an ulong.
     774             :  * Return |n|/d, set *r = |n| % d */
     775             : INLINE ulong
     776    10187366 : uabsdiviu_rem(GEN n, ulong d, ulong *r)
     777             : {
     778    10187366 :   switch(lgefint(n))
     779             :   {
     780           0 :     case 2: *r = 0; return 0;
     781    10187366 :     case 3:
     782             :     {
     783    10187366 :       ulong nn = n[2];
     784    10187366 :       *r = nn % d; return nn / d;
     785             :     }
     786           0 :     default: /* 4 */
     787             :     {
     788             :       ulong n1, n0, q;
     789             :       LOCAL_HIREMAINDER;
     790           0 :       n0 = *int_W(n,0);
     791           0 :       n1 = *int_W(n,1);
     792           0 :       hiremainder = n1;
     793           0 :       q = divll(n0, d);
     794           0 :       *r = hiremainder; return q;
     795             :     }
     796             :   }
     797             : }
     798             : 
     799             : INLINE long
     800    14961591 : sdivsi_rem(long x, GEN y, long *r)
     801             : {
     802    14961591 :   long q, s = signe(y);
     803             :   LOCAL_HIREMAINDER;
     804             : 
     805    14961591 :   if (!s) pari_err_INV("sdivsi_rem",gen_0);
     806    14961591 :   if (!x || lgefint(y)>3 || ((long)y[2]) < 0) { *r = x; return 0; }
     807    13021166 :   hiremainder=0; q = (long)divll(labs(x), (ulong)y[2]);
     808    13021166 :   if (x < 0) { hiremainder = -((long)hiremainder); q = -q; }
     809    13021166 :   if (s < 0) q = -q;
     810    13021166 :   *r = hiremainder; return q;
     811             : }
     812             : INLINE GEN
     813           0 : divsi_rem(long s, GEN y, long *r) { return stoi(sdivsi_rem(s,y,r)); }
     814             : 
     815             : INLINE long
     816         140 : sdivsi(long x, GEN y)
     817             : {
     818         140 :   long q, s = signe(y);
     819             : 
     820         140 :   if (!s) pari_err_INV("sdivsi",gen_0);
     821         140 :   if (!x || lgefint(y)>3 || ((long)y[2]) < 0) return 0;
     822          98 :   q = labs(x) / y[2];
     823          98 :   if (x < 0) q = -q;
     824          98 :   if (s < 0) q = -q;
     825          98 :   return q;
     826             : }
     827             : 
     828             : INLINE GEN
     829           0 : dvmdss(long x, long y, GEN *z)
     830             : {
     831             :   long r;
     832           0 :   GEN q = divss_rem(x,y, &r);
     833           0 :   *z = stoi(r); return q;
     834             : }
     835             : INLINE long
     836  4949802841 : dvmdsBIL(long n, long *r) { *r = remsBIL(n); return divsBIL(n); }
     837             : INLINE ulong
     838   161098197 : dvmduBIL(ulong n, ulong *r) { *r = remsBIL(n); return divsBIL(n); }
     839             : INLINE GEN
     840           0 : dvmdsi(long x, GEN y, GEN *z)
     841             : {
     842             :   long r;
     843           0 :   GEN q = divsi_rem(x,y, &r);
     844           0 :   *z = stoi(r); return q;
     845             : }
     846             : INLINE GEN
     847           0 : dvmdis(GEN x, long y, GEN *z)
     848             : {
     849             :   long r;
     850           0 :   GEN q = divis_rem(x,y, &r);
     851           0 :   *z = stoi(r); return q;
     852             : }
     853             : 
     854             : INLINE long
     855    19633892 : smodis(GEN x, long y)
     856             : {
     857    19633892 :   pari_sp av = avma;
     858    19633892 :   long r; (void)divis_rem(x,y, &r);
     859    19633892 :   return gc_long(av, (r >= 0)? r: labs(y) + r);
     860             : }
     861             : INLINE GEN
     862    19596902 : modis(GEN x, long y) { return stoi(smodis(x,y)); }
     863             : INLINE GEN
     864     8798820 : modsi(long x, GEN y) {
     865     8798820 :   long r; (void)sdivsi_rem(x, y, &r);
     866     8798818 :   return (r >= 0)? stoi(r): addsi_sign(r, y, 1);
     867             : }
     868             : 
     869             : INLINE ulong
     870      814747 : umodui(ulong x, GEN y)
     871             : {
     872      814747 :   if (!signe(y)) pari_err_INV("umodui",gen_0);
     873      814747 :   if (!x || lgefint(y) > 3) return x;
     874      206589 :   return x % (ulong)y[2];
     875             : }
     876             : 
     877             : INLINE ulong
     878      210670 : ugcdiu(GEN x, ulong y) { return ugcd(umodiu(x,y), y); }
     879             : INLINE ulong
     880        4792 : ugcdui(ulong y, GEN x) { return ugcd(umodiu(x,y), y); }
     881             : 
     882             : INLINE GEN
     883           0 : remsi(long x, GEN y)
     884           0 : { long r; (void)sdivsi_rem(x,y, &r); return stoi(r); }
     885             : INLINE GEN
     886           0 : remis(GEN x, long y)
     887             : {
     888           0 :   pari_sp av = avma;
     889             :   long r;
     890           0 :   (void)divis_rem(x,y, &r); set_avma(av); return stoi(r);
     891             : }
     892             : 
     893             : INLINE GEN
     894           0 : rdivis(GEN x, long y, long prec)
     895             : {
     896           0 :   GEN z = cgetr(prec);
     897           0 :   pari_sp av = avma;
     898           0 :   affrr(divrs(itor(x,prec), y),z);
     899           0 :   set_avma(av); return z;
     900             : }
     901             : INLINE GEN
     902           0 : rdivsi(long x, GEN y, long prec)
     903             : {
     904           0 :   GEN z = cgetr(prec);
     905           0 :   pari_sp av = avma;
     906           0 :   affrr(divsr(x, itor(y,prec)), z);
     907           0 :   set_avma(av); return z;
     908             : }
     909             : INLINE GEN
     910      839650 : rdivss(long x, long y, long prec)
     911             : {
     912      839650 :   GEN z = cgetr(prec);
     913      839650 :   pari_sp av = avma;
     914      839650 :   affrr(divrs(stor(x, prec), y), z);
     915      839650 :   set_avma(av); return z;
     916             : }
     917             : 
     918             : INLINE void
     919     7860120 : rdiviiz(GEN x, GEN y, GEN z)
     920             : {
     921     7860120 :   long prec = realprec(z), lx = lgefint(x), ly = lgefint(y);
     922     7860120 :   if (lx == 2) { affur(0, z); return; }
     923     7860120 :   if (ly == 3)
     924             :   {
     925     2084275 :     affir(x, z); if (signe(y) < 0) togglesign(z);
     926     2084207 :     affrr(divru(z, y[2]), z);
     927             :   }
     928     5775845 :   else if (lx > prec + 1 || ly > prec + 1)
     929             :   {
     930     2021187 :     affir(x,z); affrr(divri(z, y), z);
     931             :   }
     932             :   else
     933             :   {
     934     3754659 :     long b = bit_accuracy(prec) + expi(y) - expi(x) + 1;
     935     3754705 :     GEN q = divii(b > 0? shifti(x, b): x, y);
     936     3754805 :     affir(q, z); if (b > 0) shiftr_inplace(z, -b);
     937             :   }
     938     7860271 :   set_avma((ulong)z);
     939             : }
     940             : INLINE GEN
     941     7823168 : rdivii(GEN x, GEN y, long prec)
     942     7823168 : { GEN z = cgetr(prec); rdiviiz(x, y, z); return z; }
     943             : INLINE GEN
     944     7478489 : fractor(GEN x, long prec)
     945     7478489 : { return rdivii(gel(x,1), gel(x,2), prec); }
     946             : 
     947             : INLINE int
     948    13201324 : dvdii(GEN x, GEN y)
     949             : {
     950    13201324 :   pari_sp av = avma;
     951             :   GEN r;
     952    13201324 :   if (!signe(x)) return 1;
     953    11939520 :   if (!signe(y)) return 0;
     954    11939520 :   r = remii(x,y);
     955    11946513 :   return gc_bool(av, r == gen_0);
     956             : }
     957             : INLINE int
     958         154 : dvdsi(long x, GEN y)
     959             : {
     960         154 :   if (x == 0) return 1;
     961          49 :   if (!signe(y)) return 0;
     962          49 :   if (lgefint(y) != 3) return 0;
     963          49 :   return x % y[2] == 0;
     964             : }
     965             : INLINE int
     966         161 : dvdui(ulong x, GEN y)
     967             : {
     968         161 :   if (x == 0) return 1;
     969         161 :   if (!signe(y)) return 0;
     970         161 :   if (lgefint(y) != 3) return 0;
     971         161 :   return x % y[2] == 0;
     972             : }
     973             : INLINE int
     974        5663 : dvdis(GEN x, long y)
     975        5663 : { return y? smodis(x, y) == 0: signe(x) == 0; }
     976             : INLINE int
     977      413814 : dvdiu(GEN x, ulong y)
     978      413814 : { return y? umodiu(x, y) == 0: signe(x) == 0; }
     979             : 
     980             : INLINE int
     981           0 : dvdisz(GEN x, long y, GEN z)
     982             : {
     983           0 :   const pari_sp av = avma;
     984             :   long r;
     985           0 :   GEN p1 = divis_rem(x,y, &r);
     986           0 :   set_avma(av); if (r) return 0;
     987           0 :   affii(p1,z); return 1;
     988             : }
     989             : INLINE int
     990           0 : dvdiuz(GEN x, ulong y, GEN z)
     991             : {
     992           0 :   const pari_sp av = avma;
     993             :   ulong r;
     994           0 :   GEN p1 = absdiviu_rem(x,y, &r);
     995           0 :   set_avma(av); if (r) return 0;
     996           0 :   affii(p1,z); return 1;
     997             : }
     998             : INLINE int
     999        4742 : dvdiiz(GEN x, GEN y, GEN z)
    1000             : {
    1001        4742 :   const pari_sp av=avma;
    1002        4742 :   GEN p2, p1 = dvmdii(x,y,&p2);
    1003        4742 :   if (signe(p2)) return gc_bool(av,0);
    1004        4003 :   affii(p1,z); return gc_bool(av,1);
    1005             : }
    1006             : 
    1007             : INLINE ulong
    1008    77952765 : remlll_pre(ulong u2, ulong u1, ulong u0, ulong n, ulong ninv)
    1009             : {
    1010    77952765 :   u1 = remll_pre(u2, u1, n, ninv);
    1011    78975844 :   return remll_pre(u1, u0, n, ninv);
    1012             : }
    1013             : 
    1014             : INLINE ulong
    1015  1830506015 : Fl_sqr_pre(ulong a, ulong p, ulong pi)
    1016             : {
    1017             :   register ulong x;
    1018             :   LOCAL_HIREMAINDER;
    1019  1830506015 :   x = mulll(a,a);
    1020  1830506015 :   return remll_pre(hiremainder, x, p, pi);
    1021             : }
    1022             : 
    1023             : INLINE ulong
    1024  3055741138 : Fl_mul_pre(ulong a, ulong b, ulong p, ulong pi)
    1025             : {
    1026             :   register ulong x;
    1027             :   LOCAL_HIREMAINDER;
    1028  3055741138 :   x = mulll(a,b);
    1029  3055741138 :   return remll_pre(hiremainder, x, p, pi);
    1030             : }
    1031             : 
    1032             : INLINE ulong
    1033  6056853496 : Fl_addmul_pre(ulong y0, ulong x0, ulong x1, ulong p, ulong pi)
    1034             : {
    1035             :   ulong l0, h0;
    1036             :   LOCAL_HIREMAINDER;
    1037  6056853496 :   hiremainder = y0;
    1038  6056853496 :   l0 = addmul(x0, x1); h0 = hiremainder;
    1039  6056853496 :   return remll_pre(h0, l0, p, pi);
    1040             : }
    1041             : 
    1042             : INLINE ulong
    1043    52619395 : Fl_addmulmul_pre(ulong x0, ulong y0, ulong x1, ulong y1, ulong p, ulong pi)
    1044             : {
    1045             :   ulong l0, l1, h0, h1;
    1046             :   LOCAL_OVERFLOW;
    1047             :   LOCAL_HIREMAINDER;
    1048    52619395 :   l0 = mulll(x0, y0); h0 = hiremainder;
    1049    52619395 :   l1 = mulll(x1, y1); h1 = hiremainder;
    1050    52619395 :   l0 = addll(l0, l1); h0 = addllx(h0, h1);
    1051    52619395 :   return overflow ? remlll_pre(1, h0, l0, p, pi): remll_pre(h0, l0, p, pi);
    1052             : }
    1053             : 
    1054             : INLINE ulong
    1055      196161 : Fl_ellj_pre(ulong a4, ulong a6, ulong p, ulong pi)
    1056             : {
    1057             :   /* a43 = 4 a4^3 */
    1058      196161 :   ulong a43 = Fl_double(Fl_double(
    1059             :               Fl_mul_pre(a4, Fl_sqr_pre(a4, p, pi), p, pi), p), p);
    1060             :   /* a62 = 27 a6^2 */
    1061      196168 :   ulong a62 = Fl_mul_pre(Fl_sqr_pre(a6, p, pi), 27 % p, p, pi);
    1062      196168 :   ulong z1 = Fl_mul_pre(a43, 1728 % p, p, pi);
    1063      196166 :   ulong z2 = Fl_add(a43, a62, p);
    1064      196164 :   return Fl_div(z1, z2, p);
    1065             : }
    1066             : 
    1067             : /*******************************************************************/
    1068             : /*                                                                 */
    1069             : /*                        MP (INT OR REAL)                         */
    1070             : /*                                                                 */
    1071             : /*******************************************************************/
    1072             : INLINE GEN
    1073          21 : mptrunc(GEN x) { return typ(x)==t_INT? icopy(x): truncr(x); }
    1074             : INLINE GEN
    1075           0 : mpfloor(GEN x) { return typ(x)==t_INT? icopy(x): floorr(x); }
    1076             : INLINE GEN
    1077           0 : mpceil(GEN x) { return typ(x)==t_INT? icopy(x): ceilr(x); }
    1078             : INLINE GEN
    1079     1343328 : mpround(GEN x) { return typ(x) == t_INT? icopy(x): roundr(x); }
    1080             : 
    1081             : INLINE long
    1082      638363 : mpexpo(GEN x) { return typ(x) == t_INT? expi(x): expo(x); }
    1083             : 
    1084             : INLINE GEN
    1085   114832192 : mpadd(GEN x, GEN y)
    1086             : {
    1087   114832192 :   if (typ(x)==t_INT)
    1088    20337074 :     return (typ(y)==t_INT) ? addii(x,y) : addir(x,y);
    1089    94495118 :   return (typ(y)==t_INT) ? addir(y,x) : addrr(x,y);
    1090             : }
    1091             : INLINE GEN
    1092    68457502 : mpsub(GEN x, GEN y)
    1093             : {
    1094    68457502 :   if (typ(x)==t_INT)
    1095      616722 :     return (typ(y)==t_INT) ? subii(x,y) : subir(x,y);
    1096    67840780 :   return (typ(y)==t_INT) ? subri(x,y) : subrr(x,y);
    1097             : }
    1098             : INLINE GEN
    1099   191286032 : mpmul(GEN x, GEN y)
    1100             : {
    1101   191286032 :   if (typ(x)==t_INT)
    1102    36379137 :     return (typ(y)==t_INT) ? mulii(x,y) : mulir(x,y);
    1103   154906895 :   return (typ(y)==t_INT) ? mulir(y,x) : mulrr(x,y);
    1104             : }
    1105             : INLINE GEN
    1106    23910168 : mpsqr(GEN x) { return (typ(x)==t_INT) ? sqri(x) : sqrr(x); }
    1107             : INLINE GEN
    1108      660527 : mpdiv(GEN x, GEN y)
    1109             : {
    1110      660527 :   if (typ(x)==t_INT)
    1111      254302 :     return (typ(y)==t_INT) ? divii(x,y) : divir(x,y);
    1112      406225 :   return (typ(y)==t_INT) ? divri(x,y) : divrr(x,y);
    1113             : }
    1114             : 
    1115             : /*******************************************************************/
    1116             : /*                                                                 */
    1117             : /*                          Z/nZ, n ULONG                          */
    1118             : /*                                                                 */
    1119             : /*******************************************************************/
    1120             : INLINE ulong
    1121   392877094 : Fl_double(ulong a, ulong p)
    1122             : {
    1123   392877094 :   ulong res = a << 1;
    1124   392877094 :   return (res >= p || res < a) ? res - p : res;
    1125             : }
    1126             : INLINE ulong
    1127    73857385 : Fl_triple(ulong a, ulong p)
    1128             : {
    1129    73857385 :   ulong res = a << 1;
    1130    73857385 :   if (res >= p || res < a) res -= p;
    1131    73857385 :   res += a;
    1132    73857385 :   return (res >= p || res < a)? res - p: res;
    1133             : }
    1134             : INLINE ulong
    1135    15596015 : Fl_halve(ulong a, ulong p)
    1136             : {
    1137             :   ulong ap, ap2;
    1138    15596015 :   if ((a&1UL)==0) return a>>1;
    1139     7854632 :   ap = a + p; ap2 = ap>>1;
    1140     7854632 :   return ap>=a ? ap2: (ap2|HIGHBIT);
    1141             : }
    1142             : 
    1143             : INLINE ulong
    1144  3448636120 : Fl_add(ulong a, ulong b, ulong p)
    1145             : {
    1146  3448636120 :   ulong res = a + b;
    1147  3448636120 :   return (res >= p || res < a) ? res - p : res;
    1148             : }
    1149             : INLINE ulong
    1150   521011557 : Fl_neg(ulong x, ulong p) { return x ? p - x: 0; }
    1151             : 
    1152             : INLINE ulong
    1153  6374413459 : Fl_sub(ulong a, ulong b, ulong p)
    1154             : {
    1155  6374413459 :   ulong res = a - b;
    1156  6374413459 :   return (res > a) ? res + p: res;
    1157             : }
    1158             : 
    1159             : /* centerlift(u mod p) */
    1160             : INLINE long
    1161     3486617 : Fl_center(ulong u, ulong p, ulong ps2) { return (long) (u > ps2)? u - p: u; }
    1162             : 
    1163             : INLINE ulong
    1164  1195781099 : Fl_mul(ulong a, ulong b, ulong p)
    1165             : {
    1166             :   register ulong x;
    1167             :   LOCAL_HIREMAINDER;
    1168  1195781099 :   x = mulll(a,b);
    1169  1195781099 :   if (!hiremainder) return x % p;
    1170   231722196 :   (void)divll(x,p); return hiremainder;
    1171             : }
    1172             : INLINE ulong
    1173    96525606 : Fl_sqr(ulong a, ulong p)
    1174             : {
    1175             :   register ulong x;
    1176             :   LOCAL_HIREMAINDER;
    1177    96525606 :   x = mulll(a,a);
    1178    96525606 :   if (!hiremainder) return x % p;
    1179    55115939 :   (void)divll(x,p); return hiremainder;
    1180             : }
    1181             : /* don't assume that p is prime: can't special case a = 0 */
    1182             : INLINE ulong
    1183    23695269 : Fl_div(ulong a, ulong b, ulong p)
    1184    23695269 : { return Fl_mul(a, Fl_inv(b, p), p); }
    1185             : 
    1186             : /*******************************************************************/
    1187             : /*                                                                 */
    1188             : /*        DEFINED FROM EXISTING ONE EXPLOITING COMMUTATIVITY       */
    1189             : /*                                                                 */
    1190             : /*******************************************************************/
    1191             : INLINE GEN
    1192     1429639 : addri(GEN x, GEN y) { return addir(y,x); }
    1193             : INLINE GEN
    1194   103039866 : addis(GEN x, long s) { return addsi(s,x); }
    1195             : INLINE GEN
    1196    95862130 : addiu(GEN x, ulong s) { return addui(s,x); }
    1197             : INLINE GEN
    1198     6120355 : addrs(GEN x, long s) { return addsr(s,x); }
    1199             : 
    1200             : INLINE GEN
    1201    58741001 : subiu(GEN x, long y) { GEN z = subui(y, x); togglesign(z); return z; }
    1202             : INLINE GEN
    1203        2619 : subis(GEN x, long y) { return addsi(-y,x); }
    1204             : INLINE GEN
    1205     8829953 : subrs(GEN x, long y) { return addsr(-y,x); }
    1206             : 
    1207             : INLINE GEN
    1208   277412586 : mulis(GEN x, long s) { return mulsi(s,x); }
    1209             : INLINE GEN
    1210   166920378 : muliu(GEN x, ulong s) { return mului(s,x); }
    1211             : INLINE GEN
    1212     3187270 : mulru(GEN x, ulong s) { return mulur(s,x); }
    1213             : INLINE GEN
    1214    20189972 : mulri(GEN x, GEN s) { return mulir(s,x); }
    1215             : INLINE GEN
    1216    10271218 : mulrs(GEN x, long s) { return mulsr(s,x); }
    1217             : 
    1218             : /*******************************************************************/
    1219             : /*                                                                 */
    1220             : /*                  VALUATION, EXPONENT, SHIFTS                    */
    1221             : /*                                                                 */
    1222             : /*******************************************************************/
    1223             : INLINE long
    1224   124333418 : vali(GEN x)
    1225             : {
    1226             :   long i;
    1227             :   GEN xp;
    1228             : 
    1229   124333418 :   if (!signe(x)) return -1;
    1230   124253521 :   xp=int_LSW(x);
    1231   129502002 :   for (i=0; !*xp; i++) xp=int_nextW(xp);
    1232   124253521 :   return vals(*xp) + i * BITS_IN_LONG;
    1233             : }
    1234             : 
    1235             : /* assume x > 0 */
    1236             : INLINE long
    1237   474795317 : expu(ulong x) { return (BITS_IN_LONG-1) - (long)bfffo(x); }
    1238             : 
    1239             : INLINE long
    1240  1220343653 : expi(GEN x)
    1241             : {
    1242  1220343653 :   const long lx=lgefint(x);
    1243  1220343653 :   return lx==2? -(long)HIGHEXPOBIT: bit_accuracy(lx)-(long)bfffo(*int_MSW(x))-1;
    1244             : }
    1245             : 
    1246             : INLINE GEN
    1247    99931276 : shiftr(GEN x, long n)
    1248             : {
    1249    99931276 :   const long e = evalexpo(expo(x)+n);
    1250    99927445 :   const GEN y = rcopy(x);
    1251             : 
    1252    99916443 :   if (e & ~EXPOBITS) pari_err_OVERFLOW("expo()");
    1253    99916677 :   y[1] = (y[1]&~EXPOBITS) | e; return y;
    1254             : }
    1255             : INLINE GEN
    1256   132233636 : mpshift(GEN x,long s) { return (typ(x)==t_INT)?shifti(x,s):shiftr(x,s); }
    1257             : 
    1258             : /* FIXME: adapt/use mpn_[lr]shift instead */
    1259             : /* z2[imin..imax] := z1[imin..imax].f shifted left sh bits
    1260             :  * (feeding f from the right). Assume sh > 0 */
    1261             : INLINE void
    1262  5797943022 : shift_left(GEN z2, GEN z1, long imin, long imax, ulong f,  ulong sh)
    1263             : {
    1264  5797943022 :   GEN sb = z1 + imin, se = z1 + imax, te = z2 + imax;
    1265  5797943022 :   ulong l, m = BITS_IN_LONG - sh, k = f >> m;
    1266 26363276639 :   while (se > sb) {
    1267 20565333617 :     l     = *se--;
    1268 20565333617 :     *te-- = (l << sh) | k;
    1269 20565333617 :     k     = l >> m;
    1270             :   }
    1271  5797943022 :   *te = (((ulong)*se) << sh) | k;
    1272  5797943022 : }
    1273             : /* z2[imin..imax] := f.z1[imin..imax-1] shifted right sh bits
    1274             :  * (feeding f from the left). Assume sh > 0 */
    1275             : INLINE void
    1276  4716514213 : shift_right(GEN z2, GEN z1, long imin, long imax, ulong f, ulong sh)
    1277             : {
    1278  4716514213 :   GEN sb = z1 + imin, se = z1 + imax, tb = z2 + imin;
    1279  4716514213 :   ulong k, l = *sb++, m = BITS_IN_LONG - sh;
    1280  4716514213 :   *tb++ = (l >> sh) | (f << m);
    1281 25681778835 :   while (sb < se) {
    1282 20965264622 :     k     = l << m;
    1283 20965264622 :     l     = *sb++;
    1284 20965264622 :     *tb++ = (l >> sh) | k;
    1285             :   }
    1286  4716514213 : }
    1287             : 
    1288             : /* Backward compatibility. Inefficient && unused */
    1289             : extern ulong hiremainder;
    1290             : INLINE ulong
    1291           0 : shiftl(ulong x, ulong y)
    1292           0 : { hiremainder = x>>(BITS_IN_LONG-y); return (x<<y); }
    1293             : 
    1294             : INLINE ulong
    1295           0 : shiftlr(ulong x, ulong y)
    1296           0 : { hiremainder = x<<(BITS_IN_LONG-y); return (x>>y); }
    1297             : 
    1298             : INLINE void
    1299   187136688 : shiftr_inplace(GEN z, long d)
    1300             : {
    1301   187136688 :   setexpo(z, expo(z)+d);
    1302   187136798 : }
    1303             : 
    1304             : /*******************************************************************/
    1305             : /*                                                                 */
    1306             : /*                           ASSIGNMENT                            */
    1307             : /*                                                                 */
    1308             : /*******************************************************************/
    1309             : INLINE void
    1310   178664447 : affii(GEN x, GEN y)
    1311             : {
    1312   178664447 :   long lx = lgefint(x);
    1313   178664447 :   if (lg(y)<lx) pari_err_OVERFLOW("t_INT-->t_INT assignment");
    1314  5865812490 :   while (--lx) y[lx] = x[lx];
    1315   178668323 : }
    1316             : INLINE void
    1317     3406869 : affsi(long s, GEN x)
    1318             : {
    1319     3406869 :   if (!s) x[1] = evalsigne(0) | evallgefint(2);
    1320             :   else
    1321             :   {
    1322     3252321 :     if (s > 0) { x[1] = evalsigne( 1) | evallgefint(3); x[2] =  s; }
    1323     1304307 :     else       { x[1] = evalsigne(-1) | evallgefint(3); x[2] = -s; }
    1324             :   }
    1325     3406869 : }
    1326             : INLINE void
    1327     7920609 : affui(ulong u, GEN x)
    1328             : {
    1329     7920609 :   if (!u) x[1] = evalsigne(0) | evallgefint(2);
    1330     7881466 :   else  { x[1] = evalsigne(1) | evallgefint(3); x[2] = u; }
    1331     7920609 : }
    1332             : 
    1333             : INLINE void
    1334   193539153 : affsr(long x, GEN y)
    1335             : {
    1336   193539153 :   long sh, i, ly = lg(y);
    1337             : 
    1338   193539153 :   if (!x)
    1339             :   {
    1340       70966 :     y[1] = evalexpo(-prec2nbits(ly));
    1341       70966 :     return;
    1342             :   }
    1343   193468187 :   if (x < 0) {
    1344        3213 :     x = -x; sh = bfffo(x);
    1345        3213 :     y[1] = evalsigne(-1) | _evalexpo((BITS_IN_LONG-1)-sh);
    1346             :   }
    1347             :   else
    1348             :   {
    1349   193464974 :     sh = bfffo(x);
    1350   193464974 :     y[1] = evalsigne(1) | _evalexpo((BITS_IN_LONG-1)-sh);
    1351             :   }
    1352  1157100569 :   y[2] = ((ulong)x)<<sh; for (i=3; i<ly; i++) y[i]=0;
    1353             : }
    1354             : 
    1355             : INLINE void
    1356    10128307 : affur(ulong x, GEN y)
    1357             : {
    1358    10128307 :   long sh, i, ly = lg(y);
    1359             : 
    1360    10128307 :   if (!x)
    1361             :   {
    1362     1299417 :     y[1] = evalexpo(-prec2nbits(ly));
    1363     1299417 :     return;
    1364             :   }
    1365     8828890 :   sh = bfffo(x);
    1366     8828890 :   y[1] = evalsigne(1) | _evalexpo((BITS_IN_LONG-1)-sh);
    1367    29619387 :   y[2] = x<<sh; for (i=3; i<ly; i++) y[i] = 0;
    1368             : }
    1369             : 
    1370             : INLINE void
    1371      266662 : affiz(GEN x, GEN y) { if (typ(y)==t_INT) affii(x,y); else affir(x,y); }
    1372             : INLINE void
    1373           0 : affsz(long x, GEN y) { if (typ(y)==t_INT) affsi(x,y); else affsr(x,y); }
    1374             : INLINE void
    1375      679720 : mpaff(GEN x, GEN y) { if (typ(x)==t_INT) affiz(x, y); else affrr(x,y); }
    1376             : 
    1377             : /*******************************************************************/
    1378             : /*                                                                 */
    1379             : /*                    OPERATION + ASSIGNMENT                       */
    1380             : /*                                                                 */
    1381             : /*******************************************************************/
    1382             : 
    1383           0 : INLINE void addiiz(GEN x, GEN y, GEN z)
    1384           0 : { pari_sp av = avma; affii(addii(x,y),z); set_avma(av); }
    1385           0 : INLINE void addirz(GEN x, GEN y, GEN z)
    1386           0 : { pari_sp av = avma; affrr(addir(x,y),z); set_avma(av); }
    1387           0 : INLINE void addriz(GEN x, GEN y, GEN z)
    1388           0 : { pari_sp av = avma; affrr(addri(x,y),z); set_avma(av); }
    1389     1307080 : INLINE void addrrz(GEN x, GEN y, GEN z)
    1390     1307080 : { pari_sp av = avma; affrr(addrr(x,y),z); set_avma(av); }
    1391           0 : INLINE void addsiz(long s, GEN y, GEN z)
    1392           0 : { pari_sp av = avma; affii(addsi(s,y),z); set_avma(av); }
    1393           0 : INLINE void addsrz(long s, GEN y, GEN z)
    1394           0 : { pari_sp av = avma; affrr(addsr(s,y),z); set_avma(av); }
    1395           0 : INLINE void addssz(long s, long y, GEN z)
    1396           0 : { pari_sp av = avma; affii(addss(s,y),z); set_avma(av); }
    1397             : 
    1398           0 : INLINE void diviiz(GEN x, GEN y, GEN z)
    1399           0 : { pari_sp av = avma; affii(divii(x,y),z); set_avma(av); }
    1400           0 : INLINE void divirz(GEN x, GEN y, GEN z)
    1401           0 : { pari_sp av = avma; mpaff(divir(x,y),z); set_avma(av); }
    1402           0 : INLINE void divisz(GEN x, long y, GEN z)
    1403           0 : { pari_sp av = avma; affii(divis(x,y),z); set_avma(av); }
    1404           0 : INLINE void divriz(GEN x, GEN y, GEN z)
    1405           0 : { pari_sp av = avma; affrr(divri(x,y),z); set_avma(av); }
    1406         417 : INLINE void divrrz(GEN x, GEN y, GEN z)
    1407         417 : { pari_sp av = avma; affrr(divrr(x,y),z); set_avma(av); }
    1408           0 : INLINE void divrsz(GEN y, long s, GEN z)
    1409           0 : { pari_sp av = avma; affrr(divrs(y,s),z); set_avma(av); }
    1410           0 : INLINE void divsiz(long x, GEN y, GEN z)
    1411           0 : { long junk; affsi(sdivsi_rem(x,y,&junk), z); }
    1412           0 : INLINE void divsrz(long s, GEN y, GEN z)
    1413           0 : { pari_sp av = avma; mpaff(divsr(s,y),z); set_avma(av); }
    1414           0 : INLINE void divssz(long x, long y, GEN z)
    1415           0 : { affsi(x/y, z); }
    1416             : 
    1417           0 : INLINE void modisz(GEN y, long s, GEN z)
    1418           0 : { affsi(smodis(y,s),z); }
    1419           0 : INLINE void modsiz(long s, GEN y, GEN z)
    1420           0 : { pari_sp av = avma; affii(modsi(s,y),z); set_avma(av); }
    1421           0 : INLINE void modssz(long s, long y, GEN z)
    1422           0 : { affsi(smodss(s,y),z); }
    1423             : 
    1424           0 : INLINE void mpaddz(GEN x, GEN y, GEN z)
    1425           0 : { pari_sp av = avma; mpaff(mpadd(x,y),z); set_avma(av); }
    1426           0 : INLINE void mpsubz(GEN x, GEN y, GEN z)
    1427           0 : { pari_sp av = avma; mpaff(mpsub(x,y),z); set_avma(av); }
    1428           0 : INLINE void mpmulz(GEN x, GEN y, GEN z)
    1429           0 : { pari_sp av = avma; mpaff(mpmul(x,y),z); set_avma(av); }
    1430             : 
    1431           0 : INLINE void muliiz(GEN x, GEN y, GEN z)
    1432           0 : { pari_sp av = avma; affii(mulii(x,y),z); set_avma(av); }
    1433           0 : INLINE void mulirz(GEN x, GEN y, GEN z)
    1434           0 : { pari_sp av = avma; mpaff(mulir(x,y),z); set_avma(av); }
    1435           0 : INLINE void mulriz(GEN x, GEN y, GEN z)
    1436           0 : { pari_sp av = avma; mpaff(mulri(x,y),z); set_avma(av); }
    1437      260344 : INLINE void mulrrz(GEN x, GEN y, GEN z)
    1438      260344 : { pari_sp av = avma; affrr(mulrr(x,y),z); set_avma(av); }
    1439           0 : INLINE void mulsiz(long s, GEN y, GEN z)
    1440           0 : { pari_sp av = avma; affii(mulsi(s,y),z); set_avma(av); }
    1441           0 : INLINE void mulsrz(long s, GEN y, GEN z)
    1442           0 : { pari_sp av = avma; mpaff(mulsr(s,y),z); set_avma(av); }
    1443           0 : INLINE void mulssz(long s, long y, GEN z)
    1444           0 : { pari_sp av = avma; affii(mulss(s,y),z); set_avma(av); }
    1445             : 
    1446           0 : INLINE void remiiz(GEN x, GEN y, GEN z)
    1447           0 : { pari_sp av = avma; affii(remii(x,y),z); set_avma(av); }
    1448           0 : INLINE void remisz(GEN y, long s, GEN z)
    1449           0 : { pari_sp av = avma; affii(remis(y,s),z); set_avma(av); }
    1450           0 : INLINE void remsiz(long s, GEN y, GEN z)
    1451           0 : { pari_sp av = avma; affii(remsi(s,y),z); set_avma(av); }
    1452           0 : INLINE void remssz(long s, long y, GEN z)
    1453           0 : { pari_sp av = avma; affii(remss(s,y),z); set_avma(av); }
    1454             : 
    1455           0 : INLINE void subiiz(GEN x, GEN y, GEN z)
    1456           0 : { pari_sp av = avma; affii(subii(x,y),z); set_avma(av); }
    1457           0 : INLINE void subirz(GEN x, GEN y, GEN z)
    1458           0 : { pari_sp av = avma; affrr(subir(x,y),z); set_avma(av); }
    1459           0 : INLINE void subisz(GEN y, long s, GEN z)
    1460           0 : { pari_sp av = avma; affii(addsi(-s,y),z); set_avma(av); }
    1461           0 : INLINE void subriz(GEN x, GEN y, GEN z)
    1462           0 : { pari_sp av = avma; affrr(subri(x,y),z); set_avma(av); }
    1463     1296708 : INLINE void subrrz(GEN x, GEN y, GEN z)
    1464     1296708 : { pari_sp av = avma; affrr(subrr(x,y),z); set_avma(av); }
    1465           0 : INLINE void subrsz(GEN y, long s, GEN z)
    1466           0 : { pari_sp av = avma; affrr(addsr(-s,y),z); set_avma(av); }
    1467           0 : INLINE void subsiz(long s, GEN y, GEN z)
    1468           0 : { pari_sp av = avma; affii(subsi(s,y),z); set_avma(av); }
    1469           0 : INLINE void subsrz(long s, GEN y, GEN z)
    1470           0 : { pari_sp av = avma; affrr(subsr(s,y),z); set_avma(av); }
    1471           0 : INLINE void subssz(long x, long y, GEN z) { addssz(x,-y,z); }
    1472             : 
    1473             : INLINE void
    1474           0 : dvmdssz(long x, long y, GEN z, GEN t) {
    1475           0 :   pari_sp av = avma;
    1476             :   long r;
    1477           0 :   affii(divss_rem(x,y, &r), z); set_avma(av); affsi(r,t);
    1478           0 : }
    1479             : INLINE void
    1480           0 : dvmdsiz(long x, GEN y, GEN z, GEN t) {
    1481           0 :   pari_sp av = avma;
    1482             :   long r;
    1483           0 :   affii(divsi_rem(x,y, &r), z); set_avma(av); affsi(r,t);
    1484           0 : }
    1485             : INLINE void
    1486           0 : dvmdisz(GEN x, long y, GEN z, GEN t) {
    1487           0 :   pari_sp av = avma;
    1488             :   long r;
    1489           0 :   affii(divis_rem(x,y, &r),z); set_avma(av); affsi(r,t);
    1490           0 : }
    1491             : INLINE void
    1492           0 : dvmdiiz(GEN x, GEN y, GEN z, GEN t) {
    1493           0 :   pari_sp av = avma;
    1494             :   GEN r;
    1495           0 :   affii(dvmdii(x,y,&r),z); affii(r,t); set_avma(av);
    1496           0 : }

Generated by: LCOV version 1.13