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 - basemath - prime.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.12.1 lcov report (development 25819-e703fe1174) Lines: 645 734 87.9 %
Date: 2020-09-18 06:10:04 Functions: 69 75 92.0 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2000  The PARI group.
       2             : 
       3             : This file is part of the PARI/GP package.
       4             : 
       5             : PARI/GP is free software; you can redistribute it and/or modify it under the
       6             : terms of the GNU General Public License as published by the Free Software
       7             : Foundation. It is distributed in the hope that it will be useful, but WITHOUT
       8             : ANY WARRANTY WHATSOEVER.
       9             : 
      10             : Check the License for details. You should have received a copy of it, along
      11             : with the package; see the file 'COPYING'. If not, write to the Free Software
      12             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      13             : 
      14             : #include "pari.h"
      15             : #include "paripriv.h"
      16             : /*********************************************************************/
      17             : /**                                                                 **/
      18             : /**               PSEUDO PRIMALITY (MILLER-RABIN)                   **/
      19             : /**                                                                 **/
      20             : /*********************************************************************/
      21             : typedef struct {
      22             :   GEN n, sqrt1, sqrt2, t1, t;
      23             :   long r1;
      24             : } MR_Jaeschke_t;
      25             : 
      26             : static void
      27         226 : init_MR_Jaeschke(MR_Jaeschke_t *S, GEN n)
      28             : {
      29         226 :   S->n = n = absi_shallow(n);
      30         226 :   S->t = subiu(n,1);
      31         226 :   S->r1 = vali(S->t);
      32         226 :   S->t1 = shifti(S->t, -S->r1);
      33         226 :   S->sqrt1 = cgeti(lg(n)); S->sqrt1[1] = evalsigne(0)|evallgefint(2);
      34         226 :   S->sqrt2 = cgeti(lg(n)); S->sqrt2[1] = evalsigne(0)|evallgefint(2);
      35         226 : }
      36             : 
      37             : /* is n strong pseudo-prime for base a ? 'End matching' (check for square
      38             :  * roots of -1): if ends do mismatch, then we have factored n, and this
      39             :  * information should be made available to the factoring machinery. But so
      40             :  * exceedingly rare... besides we use BSPW now. */
      41             : static int
      42         383 : ispsp(MR_Jaeschke_t *S, ulong a)
      43             : {
      44         383 :   pari_sp av = avma;
      45         383 :   GEN c = Fp_pow(utoipos(a), S->t1, S->n);
      46             :   long r;
      47             : 
      48         383 :   if (is_pm1(c) || equalii(S->t, c)) return 1;
      49             :   /* go fishing for -1, not for 1 (saves one squaring) */
      50         431 :   for (r = S->r1 - 1; r; r--) /* r1 - 1 squarings */
      51             :   {
      52         362 :     GEN c2 = c;
      53         362 :     c = remii(sqri(c), S->n);
      54         362 :     if (equalii(S->t, c))
      55             :     {
      56         115 :       if (!signe(S->sqrt1))
      57             :       {
      58          72 :         affii(subii(S->n, c2), S->sqrt2);
      59          72 :         affii(c2, S->sqrt1); return 1;
      60             :       }
      61             :       /* saw one earlier: too many sqrt(-1)s mod n ? */
      62          43 :       return equalii(c2, S->sqrt1) || equalii(c2, S->sqrt2);
      63             :     }
      64         247 :     if (gc_needed(av,1))
      65             :     {
      66           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"ispsp, r = %ld", r);
      67           0 :       c = gerepileuptoint(av, c);
      68             :     }
      69             :   }
      70          69 :   return 0;
      71             : }
      72             : 
      73             : /* is n > 0 strong pseudo-prime for base 2 ? Only used when lgefint(n) > 3,
      74             :  * so don't test */
      75             : static int
      76      109712 : is2psp(GEN n)
      77             : {
      78      109712 :   GEN c, t = subiu(n, 1);
      79      109713 :   pari_sp av = avma;
      80      109713 :   long e = vali(t);
      81             : 
      82      109713 :   c = Fp_pow(gen_2, shifti(t, -e), n);
      83      109711 :   if (is_pm1(c) || equalii(t, c)) return 1;
      84      195487 :   while (--e)
      85             :   { /* go fishing for -1, not for 1 (e - 1 squaring) */
      86      104526 :     c = remii(sqri(c), n);
      87      104526 :     if (equalii(t, c)) return 1;
      88             :     /* can return 0 if (c == 1) but very infrequent */
      89       96043 :     if (gc_needed(av,1))
      90             :     {
      91           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"is2psp, r = %ld", e);
      92           0 :       c = gerepileuptoint(av, c);
      93             :     }
      94             :   }
      95       90961 :   return 0;
      96             : }
      97             : static int
      98     5318978 : uispsp_pre(ulong a, ulong n, ulong ni)
      99             : {
     100     5318978 :   ulong c, t = n - 1;
     101     5318978 :   long e = vals(t);
     102             : 
     103     5319656 :   c = Fl_powu_pre(a, t >> e, n, ni);
     104     5307324 :   if (c == 1 || c == t) return 1;
     105     9586475 :   while (--e)
     106             :   { /* go fishing for -1, not for 1 (saves one squaring) */
     107     5070653 :     c = Fl_sqr_pre(c, n, ni);
     108     5077420 :     if (c == t) return 1;
     109             :     /* can return 0 if (c == 1) but very infrequent */
     110             :   }
     111     4515822 :   return 0;
     112             : }
     113             : int
     114     6901971 : uispsp(ulong a, ulong n)
     115             : {
     116             :   ulong c, t;
     117             :   long e;
     118             : 
     119     6901971 :   if (n & HIGHMASK) return uispsp_pre(a, n, get_Fl_red(n));
     120     6751297 :   t = n - 1;
     121     6751297 :   e = vals(t);
     122     6767543 :   c = Fl_powu(a, t >> e, n);
     123     6809780 :   if (c == 1 || c == t) return 1;
     124     7686435 :   while (--e)
     125             :   { /* go fishing for -1, not for 1 (e - 1 squaring) */
     126     5283629 :     c = Fl_sqr(c, n);
     127     5284286 :     if (c == t) return 1;
     128             :     /* can return 0 if (c == 1) but very infrequent */
     129             :   }
     130     2402806 :   return 0;
     131             : }
     132             : int
     133           0 : uis2psp(ulong n) { return uispsp(2, n); }
     134             : 
     135             : /* Miller-Rabin test for k random bases */
     136             : long
     137          28 : millerrabin(GEN n, long k)
     138             : {
     139          28 :   pari_sp av2, av = avma;
     140             :   ulong r;
     141             :   long i;
     142             :   MR_Jaeschke_t S;
     143             : 
     144          28 :   if (typ(n) != t_INT) pari_err_TYPE("millerrabin",n);
     145          28 :   if (signe(n) <= 0) return 0;
     146             :   /* If |n| <= 3, check if n = +- 1 */
     147          28 :   if (lgefint(n) == 3 && uel(n,2) <= 3) return uel(n,2) != 1;
     148             : 
     149          14 :   if (!mod2(n)) return 0;
     150           7 :   init_MR_Jaeschke(&S, n); av2 = avma;
     151          21 :   for (i = 1; i <= k; i++)
     152             :   {
     153          20 :     do r = umodui(pari_rand(), n); while (!r);
     154          14 :     if (DEBUGLEVEL > 4) err_printf("Miller-Rabin: testing base %ld\n", r);
     155          14 :     if (!ispsp(&S, r)) return gc_long(av,0);
     156          14 :     set_avma(av2);
     157             :   }
     158           7 :   return gc_long(av,1);
     159             : }
     160             : 
     161             : GEN
     162          14 : gispseudoprime(GEN x, long flag)
     163          14 : { return flag? map_proto_lGL(millerrabin, x, flag): map_proto_lG(BPSW_psp,x); }
     164             : 
     165             : long
     166           0 : ispseudoprime(GEN x, long flag)
     167           0 : { return flag? millerrabin(x, flag): BPSW_psp(x); }
     168             : 
     169             : int
     170        2007 : MR_Jaeschke(GEN n)
     171             : {
     172             :   MR_Jaeschke_t S;
     173             :   pari_sp av;
     174             : 
     175        2007 :   if (lgefint(n) == 3) return uisprime(uel(n,2));
     176         219 :   if (!mod2(n)) return 0;
     177         219 :   av = avma; init_MR_Jaeschke(&S, n);
     178         219 :   return gc_int(av, ispsp(&S, 31) && ispsp(&S, 73));
     179             : }
     180             : 
     181             : /*********************************************************************/
     182             : /**                                                                 **/
     183             : /**                      PSEUDO PRIMALITY (LUCAS)                   **/
     184             : /**                                                                 **/
     185             : /*********************************************************************/
     186             : /* compute n-th term of Lucas sequence modulo N.
     187             :  * v_{k+2} = P v_{k+1} - v_k, v_0 = 2, v_1 = P.
     188             :  * Assume n > 0 */
     189             : static GEN
     190       18749 : LucasMod(GEN n, ulong P, GEN N)
     191             : {
     192       18749 :   pari_sp av = avma;
     193       18749 :   GEN nd = int_MSW(n);
     194       18749 :   ulong m = *nd;
     195             :   long i, j;
     196       18749 :   GEN v = utoipos(P), v1 = utoipos(P*P - 2);
     197             : 
     198       18749 :   if (m == 1)
     199        1811 :     j = 0;
     200             :   else
     201             :   {
     202       16938 :     j = 1+bfffo(m); /* < BIL */
     203       16938 :     m <<= j; j = BITS_IN_LONG - j;
     204             :   }
     205       18749 :   for (i=lgefint(n)-2;;) /* cf. leftright_pow */
     206             :   {
     207     2087847 :     for (; j; m<<=1,j--)
     208             :     { /* v = v_k, v1 = v_{k+1} */
     209     2039953 :       if (m&HIGHBIT)
     210             :       { /* set v = v_{2k+1}, v1 = v_{2k+2} */
     211      649629 :         v = subiu(mulii(v,v1), P);
     212      649404 :         v1= subiu(sqri(v1), 2);
     213             :       }
     214             :       else
     215             :       {/* set v = v_{2k}, v1 = v_{2k+1} */
     216     1390324 :         v1= subiu(mulii(v,v1), P);
     217     1389976 :         v = subiu(sqri(v), 2);
     218             :       }
     219     2039492 :       v = modii(v, N);
     220     2039456 :       v1= modii(v1,N);
     221     2039487 :       if (gc_needed(av,1))
     222             :       {
     223           0 :         if(DEBUGMEM>1) pari_warn(warnmem,"LucasMod");
     224           0 :         gerepileall(av, 2, &v,&v1);
     225             :       }
     226             :     }
     227       47894 :     if (--i == 0) return v;
     228       29666 :     j = BITS_IN_LONG;
     229       29666 :     nd=int_precW(nd);
     230       29666 :     m = *nd;
     231             :   }
     232             : }
     233             : /* compute n-th term of Lucas sequence modulo N.
     234             :  * v_{k+2} = P v_{k+1} - v_k, v_0 = 2, v_1 = P.
     235             :  * Assume n > 0 */
     236             : static ulong
     237      694166 : u_LucasMod_pre(ulong n, ulong P, ulong N, ulong NI)
     238             : {
     239             :   ulong v, v1, m;
     240             :   long j;
     241             : 
     242      694166 :   if (n == 1) return P;
     243      694154 :   j = 1 + bfffo(n); /* < BIL */
     244      694154 :   v = P; v1 = P*P - 2;
     245      694154 :   m = n<<j; j = BITS_IN_LONG - j;
     246    42385857 :   for (; j; m<<=1,j--)
     247             :   { /* v = v_k, v1 = v_{k+1} */
     248    41712524 :     if (m & HIGHBIT)
     249             :     { /* set v = v_{2k+1}, v1 = v_{2k+2} */
     250     4215397 :       v = Fl_sub(Fl_mul_pre(v,v1,N,NI), P, N);
     251     4215366 :       v1= Fl_sub(Fl_sqr_pre(v1,N,NI), 2UL, N);
     252             :     }
     253             :     else
     254             :     {/* set v = v_{2k}, v1 = v_{2k+1} */
     255    37497127 :       v1= Fl_sub(Fl_mul_pre(v,v1,N,NI),P, N);
     256    37491206 :       v = Fl_sub(Fl_sqr_pre(v,N,NI), 2UL, N);
     257             :     }
     258             :   }
     259      673333 :   return v;
     260             : }
     261             : 
     262             : /* !(N & HIGHMASK) */
     263             : static ulong
     264           0 : u_LucasMod(ulong n, ulong P, ulong N)
     265             : {
     266             :   ulong v, v1, m;
     267             :   long j;
     268             : 
     269           0 :   if (n == 1) return P;
     270           0 :   j = 1 + bfffo(n); /* < BIL */
     271           0 :   v = P; v1 = P*P - 2;
     272           0 :   m = n<<j; j = BITS_IN_LONG - j;
     273           0 :   for (; j; m<<=1,j--)
     274             :   { /* v = v_k, v1 = v_{k+1} */
     275           0 :     if (m & HIGHBIT)
     276             :     { /* set v = v_{2k+1}, v1 = v_{2k+2} */
     277           0 :       v = Fl_sub((v*v1) % N, P, N);
     278           0 :       v1= Fl_sub((v1*v1)% N, 2UL, N);
     279             :     }
     280             :     else
     281             :     {/* set v = v_{2k}, v1 = v_{2k+1} */
     282           0 :       v1= Fl_sub((v*v1) % N ,P, N);
     283           0 :       v = Fl_sub((v*v) % N, 2UL, N);
     284             :     }
     285             :   }
     286           0 :   return v;
     287             : }
     288             : 
     289             : static ulong
     290      693940 : get_disc(ulong n)
     291             : {
     292             :   ulong b;
     293             :   long i;
     294      693940 :   for (b = 3, i = 0;; b += 2, i++)
     295      793566 :   {
     296     1487506 :     ulong c = b*b - 4; /* = 1 mod 4 */
     297     1487506 :     if (krouu(n % c, c) < 0) break;
     298      793566 :     if (i == 64 && uissquareall(n, &c)) return 0; /* oo loop if N = m^2 */
     299             :   }
     300      694177 :   return b;
     301             : }
     302             : static int
     303      693937 : uislucaspsp_pre(ulong n, ulong ni)
     304             : {
     305             :   long i, v;
     306      693937 :   ulong b, z, m = n + 1;
     307             : 
     308      693937 :   if (!m) return 0; /* neither 2^32-1 nor 2^64-1 are Lucas-pp */
     309      693937 :   b = get_disc(n); if (!b) return 0;
     310      694173 :   v = vals(m); m >>= v;
     311      694170 :   z = u_LucasMod_pre(m, b, n, ni);
     312      694284 :   if (z == 2 || z == n-2) return 1;
     313      605241 :   for (i=1; i<v; i++)
     314             :   {
     315      605244 :     if (!z) return 1;
     316      318726 :     z = Fl_sub(Fl_sqr_pre(z,n,ni), 2UL, n);
     317      318727 :     if (z == 2) return 0;
     318             :   }
     319           0 :   return 0;
     320             : }
     321             : int
     322           0 : uislucaspsp(ulong n)
     323             : {
     324             :   long i, v;
     325             :   ulong b, z, m;
     326             : 
     327           0 :   if (n & HIGHMASK) return uislucaspsp_pre(n, get_Fl_red(n));
     328           0 :   m = n + 1;
     329           0 :   if (!m) return 0; /* neither 2^32-1 nor 2^64-1 are Lucas-pp */
     330           0 :   b = get_disc(n); if (!b) return 0;
     331           0 :   v = vals(m); m >>= v;
     332           0 :   z = u_LucasMod(m, b, n);
     333           0 :   if (z == 2 || z == n-2) return 1;
     334           0 :   for (i=1; i<v; i++)
     335             :   {
     336           0 :     if (!z) return 1;
     337           0 :     z = Fl_sub((z*z) % n, 2UL, n);
     338           0 :     if (z == 2) return 0;
     339             :   }
     340           0 :   return 0;
     341             : }
     342             : /* N > 3. Caller should check that N is not a square first (taken care of here,
     343             :  * but inefficient) */
     344             : static int
     345       18749 : islucaspsp(GEN N)
     346             : {
     347       18749 :   pari_sp av = avma;
     348             :   GEN m, z;
     349             :   long i, v;
     350             :   ulong b;
     351             : 
     352       18749 :   for (b=3;; b+=2)
     353       21669 :   {
     354       40418 :     ulong c = b*b - 4; /* = 1 mod 4 */
     355       40418 :     if (b == 129 && Z_issquare(N)) return 0; /* avoid oo loop if N = m^2 */
     356       40418 :     if (krouu(umodiu(N,c), c) < 0) break;
     357             :   }
     358       18749 :   m = addiu(N,1); v = vali(m); m = shifti(m,-v);
     359       18749 :   z = LucasMod(m, b, N);
     360       18749 :   if (absequaliu(z, 2)) return 1;
     361       16565 :   if (equalii(z, subiu(N,2))) return 1;
     362       18754 :   for (i=1; i<v; i++)
     363             :   {
     364       18587 :     if (!signe(z)) return 1;
     365       10952 :     z = modii(subiu(sqri(z), 2), N);
     366       10952 :     if (absequaliu(z, 2)) return 0;
     367       10952 :     if (gc_needed(av,1))
     368             :     {
     369           0 :       if(DEBUGMEM>1) pari_warn(warnmem,"islucaspsp");
     370           0 :       z = gerepileupto(av, z);
     371             :     }
     372             :   }
     373         167 :   return 0;
     374             : }
     375             : 
     376             : /* composite strong 2-pseudoprime < 1016801 whose prime divisors are > 101.
     377             :  * All have a prime divisor <= 661 */
     378             : static int
     379        3766 : is_2_prp_101(ulong n)
     380             : {
     381        3766 :   switch(n) {
     382           0 :   case 42799:
     383             :   case 49141:
     384             :   case 88357:
     385             :   case 90751:
     386             :   case 104653:
     387             :   case 130561:
     388             :   case 196093:
     389             :   case 220729:
     390             :   case 253241:
     391             :   case 256999:
     392             :   case 271951:
     393             :   case 280601:
     394             :   case 357761:
     395             :   case 390937:
     396             :   case 458989:
     397             :   case 486737:
     398             :   case 489997:
     399             :   case 514447:
     400             :   case 580337:
     401             :   case 741751:
     402             :   case 838861:
     403             :   case 873181:
     404             :   case 877099:
     405             :   case 916327:
     406             :   case 976873:
     407           0 :   case 983401: return 1;
     408        3766 :   } return 0;
     409             : }
     410             : 
     411             : static int
     412     6881816 : _uispsp(ulong a, long n) { a %= n; return !a || uispsp(a, n); }
     413             : static int
     414    11422080 : _uisprime(ulong n)
     415             : {
     416             : #ifdef LONG_IS_64BIT
     417    11318630 :   if (n < 341531)
     418     5085839 :     return _uispsp(9345883071009581737UL, n);
     419     6232791 :   if (n < 1050535501)
     420     1039016 :     return _uispsp(336781006125UL, n)
     421     1039014 :        && _uispsp(9639812373923155UL, n);
     422     5193775 :   if (n < 350269456337)
     423       28374 :     return _uispsp(4230279247111683200UL, n)
     424       11952 :         && _uispsp(14694767155120705706UL, n)
     425       40326 :         && _uispsp(16641139526367750375UL, n);
     426     5165401 :   if (n & HIGHMASK)
     427             :   {
     428     5165427 :     ulong ni = get_Fl_red(n);
     429     5168122 :     return (uispsp_pre(2, n, ni) && uislucaspsp_pre(n,ni));
     430             :   }
     431           0 :   return uispsp(2, n) && uislucaspsp(n);
     432             : #else
     433      103450 :   if (n < 360018361) return _uispsp(1143370UL, n) && _uispsp(2350307676UL, n);
     434        8185 :   return uispsp(15, n) && uispsp(176006322UL, n) && _uispsp(4221622697UL, n);
     435             : #endif
     436             : }
     437             : 
     438             : int
     439    56997037 : uisprime(ulong n)
     440             : {
     441    56997037 :   if (n < 103)
     442     2551817 :     switch(n)
     443             :     {
     444     1999528 :       case 2:
     445             :       case 3:
     446             :       case 5:
     447             :       case 7:
     448             :       case 11:
     449             :       case 13:
     450             :       case 17:
     451             :       case 19:
     452             :       case 23:
     453             :       case 29:
     454             :       case 31:
     455             :       case 37:
     456             :       case 41:
     457             :       case 43:
     458             :       case 47:
     459             :       case 53:
     460             :       case 59:
     461             :       case 61:
     462             :       case 67:
     463             :       case 71:
     464             :       case 73:
     465             :       case 79:
     466             :       case 83:
     467             :       case 89:
     468             :       case 97:
     469     1999528 :       case 101: return 1;
     470      552289 :       default: return 0;
     471             :     }
     472             :   /* gcd-extraction is much slower */
     473    85113068 :   return odd(n) && n % 3 && n % 5 && n % 7 && n % 11 && n % 13 && n % 17
     474    14675602 :                 && n % 19 && n % 23 && n % 29 && n % 31 && n % 37 && n % 41
     475    85110505 :                 && (n < 1849 || _uisprime(n));
     476             : }
     477             : 
     478             : /* assume no prime divisor <= 101 */
     479             : int
     480       16487 : uisprime_101(ulong n)
     481             : {
     482       16487 :   if (n < 1016801) return n < 10609? 1: (uispsp(2, n) && !is_2_prp_101(n));
     483       12708 :   return _uisprime(n);
     484             : }
     485             : 
     486             : /* assume no prime divisor <= 661 */
     487             : int
     488       26677 : uisprime_661(ulong n)
     489             : {
     490       26677 :   if (n < 1016801) return n < 452929? 1: uispsp(2, n);
     491       21324 :   return _uisprime(n);
     492             : }
     493             : 
     494             : static int
     495      306871 : iu_coprime(GEN N, ulong u) { return ugcd(u, umodiu(N,u)) == 1; }
     496             : long
     497    41581235 : BPSW_psp(GEN N)
     498             : {
     499             :   pari_sp av;
     500    41581235 :   if (typ(N) != t_INT) pari_err_TYPE("BPSW_psp",N);
     501    41588239 :   if (signe(N) <= 0) return 0;
     502    41588225 :   if (lgefint(N) == 3) return uisprime(uel(N,2));
     503      151543 :   if (!mod2(N)) return 0;
     504             : #ifdef LONG_IS_64BIT
     505             :   /* 16294579238595022365 = 3*5*7*11*13*17*19*23*29*31*37*41*43*47*53
     506             :    *  7145393598349078859 = 59*61*67*71*73*79*83*89*97*101 */
     507      126877 :   if (!iu_coprime(N, 16294579238595022365UL) ||
     508       77165 :       !iu_coprime(N,  7145393598349078859UL)) return 0;
     509             : #else
     510             :   /* 4127218095 = 3*5*7*11*13*17*19*23*37
     511             :    * 3948078067 = 29*31*41*43*47*53
     512             :    * 4269855901 = 59*83*89*97*101
     513             :    * 1673450759 = 61*67*71*73*79 */
     514      104594 :   if (!iu_coprime(N, 4127218095UL) ||
     515       82883 :       !iu_coprime(N, 3948078067UL) ||
     516       75322 :       !iu_coprime(N, 1673450759UL) ||
     517       62158 :       !iu_coprime(N, 4269855901UL)) return 0;
     518             : #endif
     519             :   /* no prime divisor < 103 */
     520       84684 :   av = avma;
     521       84684 :   return gc_long(av, is2psp(N) && islucaspsp(N));
     522             : }
     523             : 
     524             : /* can we write n = x^k ? Assume N has no prime divisor <= 2^14.
     525             :  * Not memory clean */
     526             : long
     527       23686 : isanypower_nosmalldiv(GEN N, GEN *px)
     528             : {
     529       23686 :   GEN x = N, y;
     530       23686 :   ulong mask = 7;
     531       23686 :   long ex, k = 1;
     532             :   forprime_t T;
     533       28279 :   while (Z_issquareall(x, &y)) { k <<= 1; x = y; }
     534       23799 :   while ( (ex = is_357_power(x, &y, &mask)) ) { k *= ex; x = y; }
     535       23686 :   (void)u_forprime_init(&T, 11, ULONG_MAX);
     536             :   /* stop when x^(1/k) < 2^14 */
     537       23686 :   while ( (ex = is_pth_power(x, &y, &T, 15)) ) { k *= ex; x = y; }
     538       23686 :   *px = x; return k;
     539             : }
     540             : 
     541             : /* no prime divisor <= 2^14 (> 661) */
     542             : long
     543       35587 : BPSW_psp_nosmalldiv(GEN N)
     544             : {
     545             :   pari_sp av;
     546       35587 :   long l = lgefint(N);
     547             : 
     548       35587 :   if (l == 3) return uisprime_661(uel(N,2));
     549       25032 :   av = avma;
     550             :   /* N large: test for pure power, rarely succeeds, but requires < 1% of
     551             :    * compositeness test times */
     552       25032 :   if (bit_accuracy(l) > 512 && isanypower_nosmalldiv(N,&N) != 1)
     553          14 :     return gc_long(av,0);
     554       25018 :   N = absi_shallow(N);
     555       25018 :   return gc_long(av, is2psp(N) && islucaspsp(N));
     556             : }
     557             : 
     558             : /***********************************************************************/
     559             : /**                                                                   **/
     560             : /**                       Pocklington-Lehmer                          **/
     561             : /**                        P-1 primality test                         **/
     562             : /**                                                                   **/
     563             : /***********************************************************************/
     564             : /* Assume x BPSW pseudoprime. Check whether it's small enough to be certified
     565             :  * prime (< 2^64). Reference for strong 2-pseudoprimes:
     566             :  *   http://www.cecm.sfu.ca/Pseudoprimes/index-2-to-64.html */
     567             : static int
     568     4096079 : BPSW_isprime_small(GEN x)
     569             : {
     570     4096079 :   long l = lgefint(x);
     571             : #ifdef LONG_IS_64BIT
     572     4030951 :   return (l == 3);
     573             : #else
     574       65128 :   return (l <= 4);
     575             : #endif
     576             : }
     577             : 
     578             : /* Assume N > 1, p^e || N-1, p prime. Find a witness a(p) such that
     579             :  *   a^(N-1) = 1 (mod N)
     580             :  *   a^(N-1)/p - 1 invertible mod N.
     581             :  * Proves that any divisor of N is 1 mod p^e. Return 0 if N is composite */
     582             : static ulong
     583       13542 : pl831(GEN N, GEN p)
     584             : {
     585       13542 :   GEN b, c, g, Nmunp = diviiexact(subiu(N,1), p);
     586       13542 :   pari_sp av = avma;
     587             :   ulong a;
     588       19916 :   for(a = 2;; a++, set_avma(av))
     589             :   {
     590       19916 :     b = Fp_pow(utoipos(a), Nmunp, N);
     591       19916 :     if (!equali1(b)) break;
     592             :   }
     593       13542 :   c = Fp_pow(b,p,N);
     594       13542 :   g = gcdii(subiu(b,1), N); /* 0 < g < N */
     595       13542 :   return (equali1(c) && equali1(g))? a: 0;
     596             : }
     597             : 
     598             : /* Brillhart, Lehmer, Selfridge test (Crandall & Pomerance, Th 4.1.5)
     599             :  * N^(1/3) <= f fully factored, f | N-1. If pl831(p) is true for
     600             :  * any prime divisor p of f, then any divisor of N is 1 mod f.
     601             :  * In that case return 1 iff N is prime */
     602             : static int
     603          49 : BLS_test(GEN N, GEN f)
     604             : {
     605             :   GEN c1, c2, r, q;
     606          49 :   q = dvmdii(N, f, &r);
     607          49 :   if (!is_pm1(r)) return 0;
     608          49 :   c2 = dvmdii(q, f, &c1);
     609             :   /* N = 1 + f c1 + f^2 c2, 0 <= c_i < f; check whether it is of the form
     610             :    * (1 + fa)(1 + fb) */
     611          49 :   return ! Z_issquare(subii(sqri(c1), shifti(c2,2)));
     612             : }
     613             : 
     614             : /* BPSW_psp(N) && !BPSW_isprime_small(N). Decide between Pocklington-Lehmer
     615             :  * and APRCL/ECPP. Return a vector of (small) primes such that PL-witnesses
     616             :  * guarantee the primality of N. Return NULL if PL is likely too expensive.
     617             :  * Return gen_0 if BLS test finds N to be composite */
     618             : static GEN
     619        4429 : BPSW_try_PL(GEN N)
     620             : {
     621        4429 :   ulong B = minuu(1UL<<19, maxprime());
     622        4429 :   GEN E, p, U, F, N_1 = subiu(N,1);
     623        4429 :   GEN fa = absZ_factor_limit_strict(N_1, B, &U), P = gel(fa,1);
     624             : 
     625        4429 :   if (!U) return P; /* N-1 fully factored */
     626        1721 :   p = gel(U,1);
     627        1721 :   E = gel(fa,2);
     628        1721 :   U = powii(p, gel(U,2)); /* unfactored part of N-1 */
     629        1721 :   F = (lg(P) == 2)? powii(gel(P,1), gel(E,1)): diviiexact(N_1,  U);
     630             : 
     631             :   /* N-1 = F U, F factored, U composite, (U,F) = 1 */
     632        1721 :   if (cmpii(F, U) > 0) return P; /* 1/2-smooth */
     633        1714 :   if (cmpii(sqri(F), U) > 0) return BLS_test(N,F)? P: gen_0; /* 1/3-smooth */
     634        1672 :   return NULL; /* not smooth enough */
     635             : }
     636             : 
     637             : static GEN isprimePL(GEN N);
     638             : /* F is a vector whose entries are primes. For each of them, find a PL
     639             :  * witness. Return 0 if caller lied and F contains a composite */
     640             : static long
     641        2757 : PL_certify(GEN N, GEN F)
     642             : {
     643        2757 :   long i, l = lg(F);
     644       16222 :   for(i = 1; i < l; i++)
     645       13465 :     if (! pl831(N, gel(F,i))) return 0;
     646        2757 :   return 1;
     647             : }
     648             : /* F is a vector whose entries are *believed* to be primes (BPSW_psp).
     649             :  * For each of them, recording a witness and recursive primality certificate */
     650             : static GEN
     651          98 : PL_certificate(GEN N, GEN F)
     652             : {
     653          98 :   long i, l = lg(F);
     654             :   GEN C;
     655          98 :   if (BPSW_isprime_small(N)) return N;
     656          98 :   C = cgetg(l, t_VEC);
     657         574 :   for (i = 1; i < l; i++)
     658             :   {
     659         476 :     GEN p = gel(F,i), C0;
     660             :     ulong w;
     661             : 
     662         476 :     if (BPSW_isprime_small(p)) { gel(C,i) = p; continue; }
     663          77 :     w = pl831(N,p); if (!w) return gen_0;
     664          77 :     C0 = isprimePL(p);
     665          77 :     if (isintzero(C0))
     666             :     { /* composite in prime factorisation ! */
     667           0 :       err_printf("Not a prime: %Ps", p);
     668           0 :       pari_err_BUG("PL_certificate [false prime number]");
     669             :     }
     670          77 :     gel(C,i) = mkvec3(p,utoipos(w), C0);
     671             :   }
     672          98 :   return mkvec2(N, C);
     673             : }
     674             : /* M a t_MAT */
     675             : static int
     676          84 : PL_isvalid(GEN v)
     677             : {
     678             :   GEN C, F, F2, N, N1, U;
     679             :   long i, l;
     680          84 :   switch(typ(v))
     681             :   {
     682           0 :     case t_INT: return BPSW_isprime_small(v) && BPSW_psp(v);
     683          84 :     case t_VEC: if (lg(v) == 3) break;/*fall through */
     684           0 :     default: return 0;
     685             :   }
     686          84 :   N = gel(v,1);
     687          84 :   C = gel(v,2);
     688          84 :   if (typ(N) != t_INT || signe(N) <= 0 || typ(C) != t_VEC) return 0;
     689          84 :   U = N1 = subiu(N,1);
     690          84 :   l = lg(C);
     691         427 :   for (i = 1; i < l; i++)
     692             :   {
     693         350 :     GEN p = gel(C,i), a = NULL, C0 = NULL, ap;
     694             :     long vp;
     695         350 :     if (typ(p) != t_INT)
     696             :     {
     697          70 :       if (lg(p) != 4) return 0;
     698          70 :       a = gel(p,2); C0 = gel(p,3); p = gel(p,1);
     699          70 :       if (typ(p) != t_INT || typ(a) != t_INT || !PL_isvalid(C0)) return 0;
     700             :     }
     701         350 :     vp = Z_pvalrem(U, p, &U); if (!vp) return 0;
     702         343 :     if (!a)
     703             :     {
     704         280 :       if (!BPSW_isprime_small(p)) return 0;
     705         280 :       continue;
     706             :     }
     707          63 :     if (!equalii(gel(C0,1), p)) return 0;
     708          63 :     ap = Fp_pow(a, diviiexact(N1,p), N);
     709          63 :     if (!equali1(gcdii(subiu(ap,1), N)) || !equali1(Fp_pow(ap, p, N))) return 0;
     710             :   }
     711          77 :   F = diviiexact(N1, U); /* factored part of N-1 */
     712          77 :   F2= sqri(F);
     713          77 :   if (cmpii(F2, N) > 0) return 1;
     714           0 :   if (cmpii(mulii(F2,F), N) <= 0) return 0;
     715           0 :   return BLS_test(N,F);
     716             : }
     717             : 
     718             : /* Assume N is a strong BPSW pseudoprime, Pocklington-Lehmer primality proof.
     719             :  * Return gen_0 (non-prime), N (small prime), matrix (large prime)
     720             :  *
     721             :  * The matrix has 3 columns, [a,b,c] with
     722             :  * a[i] prime factor of N-1,
     723             :  * b[i] witness for a[i] as in pl831
     724             :  * c[i] check_prime(a[i]) */
     725             : static GEN
     726         119 : isprimePL(GEN N)
     727             : {
     728             :   GEN cbrtN, N_1, F, f;
     729         119 :   if (BPSW_isprime_small(N)) return N;
     730          98 :   cbrtN = sqrtnint(N,3);
     731          98 :   N_1 = subiu(N,1);
     732          98 :   F = Z_factor_until(N_1, sqri(cbrtN));
     733          98 :   f = factorback(F); /* factored part of N-1, f^3 > N */
     734          98 :   if (DEBUGLEVEL>3)
     735             :   {
     736           0 :     GEN r = divri(itor(f,LOWDEFAULTPREC), N);
     737           0 :     err_printf("Pocklington-Lehmer: proving primality of N = %Ps\n", N);
     738           0 :     err_printf("Pocklington-Lehmer: N-1 factored up to %Ps! (%.3Ps%%)\n", f, r);
     739             :   }
     740             :   /* if N-1 is only N^(1/3)-smooth, BLS test */
     741          98 :   if (!equalii(f,N_1) && cmpii(sqri(f),N) <= 0 && !BLS_test(N,f))
     742           0 :     return gen_0; /* Failed, N is composite */
     743          98 :   F = gel(F,1); settyp(F, t_VEC);
     744          98 :   return PL_certificate(N, F);
     745             : }
     746             : 
     747             : /* assume N a BPSW pseudoprime, in particular, it is odd > 2. Prove N prime */
     748             : long
     749     4096214 : BPSW_isprime(GEN N)
     750             : {
     751             :   pari_sp av;
     752             :   long t;
     753             :   GEN P;
     754     4096214 :   if (BPSW_isprime_small(N)) return 1;
     755        4349 :   av = avma; P = BPSW_try_PL(N);
     756        4429 :   if (!P) /* not smooth enough */
     757        1672 :     t = expi(N) < 768? isprimeAPRCL(N): isprimeECPP(N);
     758             :   else
     759        2757 :     t = (typ(P) == t_INT)? 0: PL_certify(N,P);
     760        4429 :   return gc_long(av,t);
     761             : }
     762             : 
     763             : static long
     764          42 : _isprimePL(GEN x)
     765             : {
     766          42 :   pari_sp av = avma;
     767          42 :   if (!BPSW_psp(x)) return 0;
     768          35 :   return gc_long(av, !isintzero(isprimePL(x)));
     769             : }
     770             : GEN
     771    38978330 : gisprime(GEN x, long flag)
     772             : {
     773    38978330 :   switch (flag)
     774             :   {
     775    38979385 :     case 0: return map_proto_lG(isprime,x);
     776          28 :     case 1: return map_proto_lG(_isprimePL,x);
     777          14 :     case 2: return map_proto_lG(isprimeAPRCL,x);
     778          21 :     case 3: return map_proto_lG(isprimeECPP,x);
     779             :   }
     780           0 :   pari_err_FLAG("gisprime");
     781             :   return NULL;/*LCOV_EXCL_LINE*/
     782             : }
     783             : 
     784             : long
     785    39666209 : isprime(GEN x) { return BPSW_psp(x) && BPSW_isprime(x); }
     786             : 
     787             : GEN
     788          70 : primecert(GEN x, long flag)
     789             : {
     790          70 :   if (!BPSW_psp(x)) return gen_0;
     791          63 :   switch(flag)
     792             :   {
     793          56 :     case 0: return ecpp(x);
     794           7 :     case 1: { pari_sp av = avma; return gerepilecopy(av, isprimePL(x)); }
     795             :   }
     796           0 :   pari_err_FLAG("primecert");
     797             :   return NULL;/*LCOV_EXCL_LINE*/
     798             : }
     799             : 
     800             : enum { c_VOID = 0, c_ECPP, c_N1 };
     801             : static long
     802          91 : cert_type(GEN c)
     803             : {
     804          91 :   switch(typ(c))
     805             :   {
     806           0 :     case t_INT: return c_ECPP;
     807          91 :     case t_VEC:
     808          91 :       if (lg(c) == 3 && typ(gel(c,1)) == t_INT) return c_N1;
     809          77 :       return c_ECPP;
     810             :   }
     811           0 :   return c_VOID;
     812             : }
     813             : 
     814             : long
     815          49 : primecertisvalid(GEN c)
     816             : {
     817          49 :   switch(typ(c))
     818             :   {
     819           7 :     case t_INT: return BPSW_isprime_small(c) && BPSW_psp(c);
     820          42 :     case t_VEC:
     821          42 :       return cert_type(c) == c_ECPP? ecppisvalid(c): PL_isvalid(c);
     822             :   }
     823           0 :   return 0;
     824             : }
     825             : 
     826             : static long
     827         329 : check_eccpcertentry(GEN c)
     828             : {
     829             :   GEN v;
     830         329 :   long i,l = lg(c);
     831         329 :   if (typ(c)!=t_VEC || l!=6) return 0;
     832        1610 :   for(i=1; i<=4; i++)
     833        1288 :     if (typ(gel(c,i))!=t_INT) return 0;
     834         322 :   v = gel(c,5);
     835         322 :   if(typ(v)!=t_VEC) return 0;
     836         966 :   for(i=1; i<=2; i++)
     837         644 :     if (typ(gel(v,i))!=t_INT) return 0;
     838         322 :   return 1;
     839             : }
     840             : 
     841             : static long
     842          49 : check_eccpcert(GEN c)
     843             : {
     844             :   long i, l;
     845          49 :   switch(typ(c))
     846             :   {
     847           0 :     case t_INT: return signe(c) >= 0;
     848          49 :     case t_VEC: break;
     849           0 :     default: return 0;
     850             :   }
     851          49 :   l = lg(c); if (l == 1) return 0;
     852         364 :   for (i = 1; i < l; i++)
     853         329 :     if (check_eccpcertentry(gel(c,i))==0) return 0;
     854          35 :   return 1;
     855             : }
     856             : 
     857             : GEN
     858          49 : primecertexport(GEN c, long flag)
     859             : {
     860          49 :   if (cert_type(c) != c_ECPP) pari_err_IMPL("N-1 certificate");
     861          49 :   if (!check_eccpcert(c))
     862          14 :     pari_err_TYPE("primecertexport - invalid certificate", c);
     863          35 :   return ecppexport(c, flag);
     864             : }
     865             : 
     866             : /***********************************************************************/
     867             : /**                                                                   **/
     868             : /**                          PRIME NUMBERS                            **/
     869             : /**                                                                   **/
     870             : /***********************************************************************/
     871             : 
     872             : static struct {
     873             :   ulong p;
     874             :   long n;
     875             : } prime_table[] = {
     876             :   {           0,          0},
     877             :   {        7919,       1000},
     878             :   {       17389,       2000},
     879             :   {       27449,       3000},
     880             :   {       37813,       4000},
     881             :   {       48611,       5000},
     882             :   {       59359,       6000},
     883             :   {       70657,       7000},
     884             :   {       81799,       8000},
     885             :   {       93179,       9000},
     886             :   {      104729,      10000},
     887             :   {      224737,      20000},
     888             :   {      350377,      30000},
     889             :   {      479909,      40000},
     890             :   {      611953,      50000},
     891             :   {      746773,      60000},
     892             :   {      882377,      70000},
     893             :   {     1020379,      80000},
     894             :   {     1159523,      90000},
     895             :   {     1299709,     100000},
     896             :   {     2750159,     200000},
     897             :   {     7368787,     500000},
     898             :   {    15485863,    1000000},
     899             :   {    32452843,    2000000},
     900             :   {    86028121,    5000000},
     901             :   {   179424673,   10000000},
     902             :   {   373587883,   20000000},
     903             :   {   982451653,   50000000},
     904             :   {  2038074743,  100000000},
     905             :   {  4000000483UL,189961831},
     906             :   {  4222234741UL,200000000},
     907             : #if BITS_IN_LONG == 64
     908             :   { 5336500537UL,   250000000L},
     909             :   { 6461335109UL,   300000000L},
     910             :   { 7594955549UL,   350000000L},
     911             :   { 8736028057UL,   400000000L},
     912             :   { 9883692017UL,   450000000L},
     913             :   { 11037271757UL,  500000000L},
     914             :   { 13359555403UL,  600000000L},
     915             :   { 15699342107UL,  700000000L},
     916             :   { 18054236957UL,  800000000L},
     917             :   { 20422213579UL,  900000000L},
     918             :   { 22801763489UL, 1000000000L},
     919             :   { 47055833459UL, 2000000000L},
     920             :   { 71856445751UL, 3000000000L},
     921             :   { 97011687217UL, 4000000000L},
     922             :   {122430513841UL, 5000000000L},
     923             :   {148059109201UL, 6000000000L},
     924             :   {173862636221UL, 7000000000L},
     925             :   {200000000507UL, 8007105083L},
     926             :   {225898512559UL, 9000000000L},
     927             :   {252097800623UL,10000000000L},
     928             :   {384489816343UL,15000000000L},
     929             :   {518649879439UL,20000000000L},
     930             :   {654124187867UL,25000000000L},
     931             :   {790645490053UL,30000000000L},
     932             :   {928037044463UL,35000000000L},
     933             :  {1066173339601UL,40000000000L},
     934             :  {1344326694119UL,50000000000L},
     935             :  {1624571841097UL,60000000000L},
     936             :  {1906555030411UL,70000000000L},
     937             :  {2190026988349UL,80000000000L},
     938             :  {2474799787573UL,90000000000L},
     939             :  {2760727302517UL,100000000000L}
     940             : #endif
     941             : };
     942             : static const int prime_table_len = numberof(prime_table);
     943             : 
     944             : /* find prime closest to n in prime_table. */
     945             : static long
     946    18097751 : prime_table_closest_p(ulong n)
     947             : {
     948             :   long i;
     949    18211580 :   for (i = 1; i < prime_table_len; i++)
     950             :   {
     951    18211614 :     ulong p = prime_table[i].p;
     952    18211614 :     if (p > n)
     953             :     {
     954    18097785 :       ulong u = n - prime_table[i-1].p;
     955    18097785 :       if (p - n > u) i--;
     956    18097785 :       break;
     957             :     }
     958             :   }
     959    18097751 :   if (i == prime_table_len) i = prime_table_len - 1;
     960    18097751 :   return i;
     961             : }
     962             : 
     963             : /* return the n-th successor of prime p > 2; n > 0 */
     964             : static GEN
     965          77 : prime_successor(ulong p, ulong n)
     966             : {
     967          77 :   GEN Q = utoipos(p), P = NULL;
     968             :   ulong i;
     969          77 : RESET:
     970             :   for(;;)
     971             :   {
     972             :     forprime_t S;
     973          77 :     GEN L = dbltor(4*n * log(gtodouble(Q) + 1)), R = addii(Q, ceil_safe(L));
     974             : 
     975          77 :     forprime_init(&S, addiu(Q, 1), R);
     976          77 :     Q = R;
     977     2402708 :     for (i = 1; i <= n; i++)
     978             :     {
     979     2402631 :       P = forprime_next(&S);
     980     2402631 :       if (!P) { n -= i-1; goto RESET; }
     981             :     }
     982          77 :     return P;
     983             :   }
     984             : }
     985             : /* find the N-th prime */
     986             : static GEN
     987         273 : prime_table_find_n(ulong N)
     988             : {
     989             :   byteptr d;
     990         273 :   ulong n, p, maxp = maxprime();
     991             :   long i;
     992        2240 :   for (i = 1; i < prime_table_len; i++)
     993             :   {
     994        2240 :     n = prime_table[i].n;
     995        2240 :     if (n > N)
     996             :     {
     997         273 :       ulong u = N - prime_table[i-1].n;
     998         273 :       if (n - N > u) i--;
     999         273 :       break;
    1000             :     }
    1001             :   }
    1002         273 :   if (i == prime_table_len) i = prime_table_len - 1;
    1003         273 :   p = prime_table[i].p;
    1004         273 :   n = prime_table[i].n;
    1005         273 :   if (n > N && p > maxp)
    1006             :   {
    1007          14 :     i--;
    1008          14 :     p = prime_table[i].p;
    1009          14 :     n = prime_table[i].n;
    1010             :   }
    1011             :   /* if beyond prime table, then n <= N */
    1012         273 :   d = diffptr + n;
    1013         273 :   if (n > N)
    1014             :   {
    1015          14 :     n -= N;
    1016       50624 :     do { n--; PREC_PRIME_VIADIFF(p,d); } while (n) ;
    1017             :   }
    1018         259 :   else if (n < N)
    1019             :   {
    1020         259 :     ulong maxpN = maxprimeN();
    1021         259 :     if (N >= maxpN)
    1022             :     {
    1023          77 :       if (N == maxpN) return utoipos(maxp);
    1024          77 :       if (p < maxp) { p = maxp; n = maxpN; }
    1025          77 :       return prime_successor(p, N - n);
    1026             :     }
    1027             :     /* can find prime(N) in table */
    1028         182 :     n = N - n;
    1029       45234 :     do { n--; NEXT_PRIME_VIADIFF(p,d); } while (n) ;
    1030             :   }
    1031         196 :   return utoipos(p);
    1032             : }
    1033             : 
    1034             : ulong
    1035           0 : uprime(long N)
    1036             : {
    1037           0 :   pari_sp av = avma;
    1038             :   GEN p;
    1039           0 :   if (N <= 0) pari_err_DOMAIN("prime", "n", "<=",gen_0, stoi(N));
    1040           0 :   p = prime_table_find_n(N);
    1041           0 :   if (lgefint(p) != 3) pari_err_OVERFLOW("uprime");
    1042           0 :   return gc_ulong(av, p[2]);
    1043             : }
    1044             : GEN
    1045         280 : prime(long N)
    1046             : {
    1047         280 :   pari_sp av = avma;
    1048             :   GEN p;
    1049         280 :   if (N <= 0) pari_err_DOMAIN("prime", "n", "<=",gen_0, stoi(N));
    1050         273 :   new_chunk(4); /*HACK*/
    1051         273 :   p = prime_table_find_n(N);
    1052         273 :   set_avma(av); return icopy(p);
    1053             : }
    1054             : 
    1055             : static void
    1056          98 : prime_interval(GEN N, GEN *pa, GEN *pb, GEN *pd)
    1057             : {
    1058             :   GEN a, b, d;
    1059          98 :   switch(typ(N))
    1060             :   {
    1061          56 :     case t_INT:
    1062          56 :       a = gen_2;
    1063          56 :       b = subiu(N,1); /* between 2 and N-1 */
    1064          56 :       d = subiu(N,2);
    1065          56 :       if (signe(d) <= 0)
    1066           7 :         pari_err_DOMAIN("randomprime","N", "<=", gen_2, N);
    1067          49 :       break;
    1068          42 :     case t_VEC:
    1069          42 :       if (lg(N) != 3) pari_err_TYPE("randomprime",N);
    1070          42 :       a = gel(N,1);
    1071          42 :       b = gel(N,2);
    1072          42 :       if (gcmp(b, a) < 0)
    1073           7 :         pari_err_DOMAIN("randomprime","b-a", "<", gen_0, mkvec2(a,b));
    1074          35 :       if (typ(a) != t_INT)
    1075             :       {
    1076           7 :         a = gceil(a);
    1077           7 :         if (typ(a) != t_INT) pari_err_TYPE("randomprime",a);
    1078             :       }
    1079          35 :       if (typ(b) != t_INT)
    1080             :       {
    1081           7 :         b = gfloor(b);
    1082           7 :         if (typ(b) != t_INT) pari_err_TYPE("randomprime",b);
    1083             :       }
    1084          35 :       if (cmpiu(a, 2) < 0)
    1085             :       {
    1086           7 :         a = gen_2;
    1087           7 :         d = subiu(b,1);
    1088             :       }
    1089             :       else
    1090          28 :         d = addiu(subii(b,a), 1);
    1091          35 :       if (signe(d) <= 0)
    1092          14 :         pari_err_DOMAIN("randomprime","floor(b) - max(ceil(a),2)", "<",
    1093             :                         gen_0, mkvec2(a,b));
    1094          21 :       break;
    1095           0 :     default:
    1096           0 :       pari_err_TYPE("randomprime", N);
    1097             :       a = b = d = NULL; /*LCOV_EXCL_LINE*/
    1098             :   }
    1099          70 :   *pa = a; *pb = b; *pd = d;
    1100          70 : }
    1101             : 
    1102             : /* random b-bit prime */
    1103             : GEN
    1104          56 : randomprime(GEN N)
    1105             : {
    1106          56 :   pari_sp av = avma, av2;
    1107             :   GEN a, b, d;
    1108          56 :   if (!N)
    1109             :     for(;;)
    1110          56 :     {
    1111          63 :       ulong p = random_bits(31);
    1112          63 :       if (uisprime(p)) return utoipos(p);
    1113             :     }
    1114          49 :   prime_interval(N, &a, &b, &d); av2 = avma;
    1115             :   for (;;)
    1116        5889 :   {
    1117        5910 :     GEN p = addii(a, randomi(d));
    1118        5910 :     if (BPSW_psp(p)) return gerepileuptoint(av, p);
    1119        5889 :     set_avma(av2);
    1120             :   }
    1121             : }
    1122             : GEN
    1123         105 : randomprime0(GEN N, GEN q)
    1124             : {
    1125         105 :   pari_sp av = avma, av2;
    1126             :   GEN C, D, a, b, d, r;
    1127         105 :   if (!q) return randomprime(N);
    1128          49 :   switch(typ(q))
    1129             :   {
    1130          28 :     case t_INT: C = gen_1; D = q; break;
    1131          21 :     case t_INTMOD: C = gel(q,2); D = gel(q,1); break;
    1132           0 :     default:
    1133           0 :       pari_err_TYPE("randomprime", q);
    1134             :       return NULL;/*LCOV_EXCL_LINE*/
    1135             :   }
    1136          49 :   if (!N) N = int2n(31);
    1137          49 :   prime_interval(N, &a, &b, &d);
    1138          49 :   r = modii(subii(C, a), D);
    1139          49 :   if (signe(r)) { a = addii(a, r); d = subii(d, r); }
    1140          49 :   if (!equali1(gcdii(C,D)))
    1141             :   {
    1142          14 :     if (isprime(a)) return gerepilecopy(av, a);
    1143           7 :     pari_err_COPRIME("randomprime", C, D);
    1144             :   }
    1145          35 :   d = divii(d, D); if (!signe(d)) d = gen_1;
    1146          35 :   av2 = avma;
    1147             :   for (;;)
    1148        3584 :   {
    1149        3619 :     GEN p = addii(a, mulii(D, randomi(d)));
    1150        3619 :     if (BPSW_psp(p)) return gerepileuptoint(av, p);
    1151        3584 :     set_avma(av2);
    1152             :   }
    1153             :   return NULL;
    1154             : }
    1155             : 
    1156             : /* set *pp = nextprime(a) = p
    1157             :  *     *pd so that NEXT_PRIME_VIADIFF(d, p) = nextprime(p+1)
    1158             :  *     *pn so that p = the n-th prime
    1159             :  * error if nextprime(a) is out of primetable bounds */
    1160             : void
    1161    18097652 : prime_table_next_p(ulong a, byteptr *pd, ulong *pp, ulong *pn)
    1162             : {
    1163             :   byteptr d;
    1164    18097652 :   ulong p, n, maxp = maxprime();
    1165    18097657 :   long i = prime_table_closest_p(a);
    1166    18097663 :   p = prime_table[i].p;
    1167    18097663 :   if (p > a && p > maxp)
    1168             :   {
    1169           0 :     i--;
    1170           0 :     p = prime_table[i].p;
    1171             :   }
    1172             :   /* if beyond prime table, then p <= a */
    1173    18097663 :   n = prime_table[i].n;
    1174    18097663 :   d = diffptr + n;
    1175    18097663 :   if (p < a)
    1176             :   {
    1177    18064582 :     if (a > maxp) pari_err_MAXPRIME(a);
    1178    44851892 :     do { n++; NEXT_PRIME_VIADIFF(p,d); } while (p < a);
    1179             :   }
    1180       33081 :   else if (p != a)
    1181             :   {
    1182    14273170 :     do { n--; PREC_PRIME_VIADIFF(p,d); } while (p > a) ;
    1183       33081 :     if (p < a) { NEXT_PRIME_VIADIFF(p,d); n++; }
    1184             :   }
    1185    18097663 :   *pn = n;
    1186    18097663 :   *pp = p;
    1187    18097663 :   *pd = d;
    1188    18097663 : }
    1189             : 
    1190             : ulong
    1191       17331 : uprimepi(ulong a)
    1192             : {
    1193       17331 :   ulong p, n, maxp = maxprime();
    1194       17331 :   if (a <= maxp)
    1195             :   {
    1196             :     byteptr d;
    1197       17206 :     prime_table_next_p(a, &d, &p, &n);
    1198       17206 :     return p == a? n: n-1;
    1199             :   }
    1200             :   else
    1201             :   {
    1202         125 :     long i = prime_table_closest_p(a);
    1203             :     forprime_t S;
    1204         125 :     p = prime_table[i].p;
    1205         125 :     if (p > a)
    1206             :     {
    1207          28 :       i--;
    1208          28 :       p = prime_table[i].p;
    1209             :     }
    1210             :     /* p = largest prime in table <= a */
    1211         125 :     n = prime_table[i].n;
    1212         125 :     (void)u_forprime_init(&S, p+1, a);
    1213    52168646 :     for (; p; n++) p = u_forprime_next(&S);
    1214         125 :     return n-1;
    1215             :   }
    1216             : }
    1217             : 
    1218             : GEN
    1219         252 : primepi(GEN x)
    1220             : {
    1221         252 :   pari_sp av = avma;
    1222         252 :   GEN pp, nn, N = typ(x) == t_INT? x: gfloor(x);
    1223             :   forprime_t S;
    1224             :   ulong n, p;
    1225             :   long i;
    1226         252 :   if (typ(N) != t_INT) pari_err_TYPE("primepi",N);
    1227         252 :   if (signe(N) <= 0) return gen_0;
    1228         252 :   if (lgefint(N) == 3) { n = N[2]; set_avma(av); return utoi(uprimepi(n)); }
    1229           1 :   i = prime_table_len-1;
    1230           1 :   p = prime_table[i].p;
    1231           1 :   n = prime_table[i].n;
    1232           1 :   (void)forprime_init(&S, utoipos(p+1), N);
    1233           1 :   nn = setloop(utoipos(n));
    1234           1 :   pp = gen_0;
    1235     3280223 :   for (; pp; incloop(nn)) pp = forprime_next(&S);
    1236           1 :   return gerepileuptoint(av, subiu(nn,1));
    1237             : }
    1238             : 
    1239             : /* pi(x) < x/log x * (1 + 1/log x + 2.51/log^2 x)), x>=355991 [ Dusart ]
    1240             :  * pi(x) < x/(log x - 1.1), x >= 60184 [ Dusart ]
    1241             :  * ? \p9
    1242             :  * ? M = 0; for(x = 4, 60184, M = max(M, log(x) - x/primepi(x))); M
    1243             :  * %1 = 1.11196252 */
    1244             : double
    1245       78728 : primepi_upper_bound(double x)
    1246             : {
    1247       78728 :   if (x >= 355991)
    1248             :   {
    1249          77 :     double L = 1/log(x);
    1250          77 :     return x * L * (1 + L + 2.51*L*L);
    1251             :   }
    1252       78651 :   if (x >= 60184) return x / (log(x) - 1.1);
    1253       78651 :   if (x < 5) return 2; /* don't bother */
    1254       54613 :   return x / (log(x) - 1.111963);
    1255             : }
    1256             : /* pi(x) > x/log x (1 + 1/log x), x >= 599 [ Dusart ]
    1257             :  * pi(x) > x / (log x + 2), x >= 55 [ Rosser ] */
    1258             : double
    1259          14 : primepi_lower_bound(double x)
    1260             : {
    1261          14 :   if (x >= 599)
    1262             :   {
    1263          14 :     double L = 1/log(x);
    1264          14 :     return x * L * (1 + L);
    1265             :   }
    1266           0 :   if (x < 55) return 0; /* don't bother */
    1267           0 :   return x / (log(x) + 2.);
    1268             : }
    1269             : GEN
    1270           1 : gprimepi_upper_bound(GEN x)
    1271             : {
    1272           1 :   pari_sp av = avma;
    1273             :   double L;
    1274           1 :   if (typ(x) != t_INT) x = gfloor(x);
    1275           1 :   if (expi(x) <= 1022)
    1276             :   {
    1277           1 :     set_avma(av);
    1278           1 :     return dbltor(primepi_upper_bound(gtodouble(x)));
    1279             :   }
    1280           0 :   x = itor(x, LOWDEFAULTPREC);
    1281           0 :   L = 1 / rtodbl(logr_abs(x));
    1282           0 :   x = mulrr(x, dbltor(L * (1 + L + 2.51*L*L)));
    1283           0 :   return gerepileuptoleaf(av, x);
    1284             : }
    1285             : GEN
    1286           1 : gprimepi_lower_bound(GEN x)
    1287             : {
    1288           1 :   pari_sp av = avma;
    1289             :   double L;
    1290           1 :   if (typ(x) != t_INT) x = gfloor(x);
    1291           1 :   if (abscmpiu(x, 55) <= 0) return gen_0;
    1292           1 :   if (expi(x) <= 1022)
    1293             :   {
    1294           1 :     set_avma(av);
    1295           1 :     return dbltor(primepi_lower_bound(gtodouble(x)));
    1296             :   }
    1297           0 :   x = itor(x, LOWDEFAULTPREC);
    1298           0 :   L = 1 / rtodbl(logr_abs(x));
    1299           0 :   x = mulrr(x, dbltor(L * (1 + L)));
    1300           0 :   return gerepileuptoleaf(av, x);
    1301             : }
    1302             : 
    1303             : GEN
    1304          77 : primes(long n)
    1305             : {
    1306             :   forprime_t S;
    1307             :   long i;
    1308             :   GEN y;
    1309          77 :   if (n <= 0) return cgetg(1, t_VEC);
    1310          77 :   y = cgetg(n+1, t_VEC);
    1311          77 :   (void)new_chunk(3*n); /*HACK*/
    1312          77 :   u_forprime_init(&S, 2, ULONG_MAX);
    1313          77 :   set_avma((pari_sp)y);
    1314        3871 :   for (i = 1; i <= n; i++) gel(y, i) = utoipos( u_forprime_next(&S) );
    1315          77 :   return y;
    1316             : }
    1317             : GEN
    1318           0 : primes_zv(long n)
    1319             : {
    1320             :   forprime_t S;
    1321             :   long i;
    1322             :   GEN y;
    1323           0 :   if (n <= 0) return cgetg(1, t_VECSMALL);
    1324           0 :   y = cgetg(n+1,t_VECSMALL);
    1325           0 :   u_forprime_init(&S, 2, ULONG_MAX);
    1326           0 :   for (i = 1; i <= n; i++) y[i] =  u_forprime_next(&S);
    1327           0 :   set_avma((pari_sp)y); return y;
    1328             : }
    1329             : GEN
    1330         175 : primes0(GEN N)
    1331             : {
    1332         175 :   switch(typ(N))
    1333             :   {
    1334          77 :     case t_INT: return primes(itos(N));
    1335          98 :     case t_VEC:
    1336          98 :       if (lg(N) == 3) return primes_interval(gel(N,1),gel(N,2));
    1337             :   }
    1338           0 :   pari_err_TYPE("primes", N);
    1339             :   return NULL;/*LCOV_EXCL_LINE*/
    1340             : }
    1341             : 
    1342             : GEN
    1343         154 : primes_interval(GEN a, GEN b)
    1344             : {
    1345         154 :   pari_sp av = avma;
    1346             :   forprime_t S;
    1347             :   long i, n;
    1348             :   GEN y, d, p;
    1349         154 :   if (typ(a) != t_INT)
    1350             :   {
    1351           0 :     a = gceil(a);
    1352           0 :     if (typ(a) != t_INT) pari_err_TYPE("primes_interval",a);
    1353             :   }
    1354         154 :   if (typ(b) != t_INT)
    1355             :   {
    1356           7 :     b = gfloor(b);
    1357           7 :     if (typ(b) != t_INT) pari_err_TYPE("primes_interval",b);
    1358             :   }
    1359         147 :   if (signe(a) < 0) a = gen_2;
    1360         147 :   d = subii(b, a);
    1361         147 :   if (signe(d) < 0 || signe(b) <= 0) { set_avma(av); return cgetg(1, t_VEC); }
    1362         147 :   if (lgefint(b) == 3)
    1363             :   {
    1364         129 :     set_avma(av);
    1365         129 :     y = primes_interval_zv(itou(a), itou(b));
    1366         129 :     n = lg(y); settyp(y, t_VEC);
    1367      469848 :     for (i = 1; i < n; i++) gel(y,i) = utoipos(y[i]);
    1368         129 :     return y;
    1369             :   }
    1370             :   /* at most d+1 primes in [a,b]. If d large, try better bound to lower
    1371             :    * memory use */
    1372          18 :   if (abscmpiu(d,100000) > 0)
    1373             :   {
    1374           1 :     GEN D = gsub(gprimepi_upper_bound(b), gprimepi_lower_bound(a));
    1375           1 :     D = ceil_safe(D);
    1376           1 :     if (cmpii(D, d) < 0) d = D;
    1377             :   }
    1378          18 :   n = itos(d)+1;
    1379          18 :   forprime_init(&S, a, b);
    1380          18 :   y = cgetg(n+1, t_VEC); i = 1;
    1381        5900 :   while ((p = forprime_next(&S))) gel(y, i++) = icopy(p);
    1382          18 :   setlg(y, i); return gerepileupto(av, y);
    1383             : }
    1384             : 
    1385             : /* a <= b, at most d primes in [a,b]. Return them */
    1386             : static GEN
    1387        7465 : primes_interval_i(ulong a, ulong b, ulong d)
    1388             : {
    1389        7465 :   ulong p, i = 1, n = d + 1;
    1390             :   forprime_t S;
    1391        7465 :   GEN y = cgetg(n+1, t_VECSMALL);
    1392        7465 :   pari_sp av = avma;
    1393        7465 :   u_forprime_init(&S, a, b);
    1394    13277377 :   while ((p = u_forprime_next(&S))) y[i++] = p;
    1395        7465 :   set_avma(av); setlg(y, i); stackdummy((pari_sp)(y + i), (pari_sp)(y + n+1));
    1396        7465 :   return y;
    1397             : }
    1398             : GEN
    1399        7269 : primes_interval_zv(ulong a, ulong b)
    1400             : {
    1401             :   ulong d;
    1402        7269 :   if (!a) return primes_upto_zv(b);
    1403        7269 :   if (b < a) return cgetg(1, t_VECSMALL);
    1404        7269 :   d = b - a;
    1405        7269 :   if (d > 100000UL)
    1406             :   {
    1407          13 :     ulong D = (ulong)ceil(primepi_upper_bound(b)-primepi_lower_bound(a));
    1408          13 :     if (D < d) d = D;
    1409             :   }
    1410        7269 :   return primes_interval_i(a, b, d);
    1411             : }
    1412             : GEN
    1413         196 : primes_upto_zv(ulong b)
    1414             : {
    1415             :   ulong d;
    1416         196 :   if (b < 2) return cgetg(1, t_VECSMALL);
    1417         196 :   d = (b > 100000UL)? (ulong)primepi_upper_bound(b): b;
    1418         196 :   return primes_interval_i(2, b, d);
    1419             : }
    1420             : 
    1421             : /***********************************************************************/
    1422             : /**                                                                   **/
    1423             : /**                       PRIVATE PRIME TABLE                         **/
    1424             : /**                                                                   **/
    1425             : /***********************************************************************/
    1426             : 
    1427             : void
    1428      176668 : pari_set_primetab(GEN global_primetab)
    1429             : {
    1430      176668 :   if (global_primetab)
    1431             :   {
    1432      174987 :     long i, l = lg(global_primetab);
    1433      174987 :     primetab = cgetg_block(l, t_VEC);
    1434      203374 :     for (i = 1; i < l; i++)
    1435       28458 :       gel(primetab,i) = gclone(gel(global_primetab,i));
    1436        1681 :   } else primetab = cgetg_block(1, t_VEC);
    1437      176606 : }
    1438             : 
    1439             : /* delete dummy NULL entries */
    1440             : static void
    1441          21 : cleanprimetab(GEN T)
    1442             : {
    1443          21 :   long i,j, l = lg(T);
    1444          70 :   for (i = j = 1; i < l; i++)
    1445          49 :     if (T[i]) T[j++] = T[i];
    1446          21 :   setlg(T,j);
    1447          21 : }
    1448             : /* remove p from T */
    1449             : static void
    1450          28 : rmprime(GEN T, GEN p)
    1451             : {
    1452             :   long i;
    1453          28 :   if (typ(p) != t_INT) pari_err_TYPE("removeprimes",p);
    1454          28 :   i = ZV_search(T, p);
    1455          28 :   if (!i)
    1456           7 :     pari_err_DOMAIN("removeprimes","prime","not in",
    1457             :                     strtoGENstr("primetable"), p);
    1458          21 :   gunclone(gel(T,i)); gel(T,i) = NULL;
    1459          21 :   cleanprimetab(T);
    1460          21 : }
    1461             : 
    1462             : /*stolen from ZV_union_shallow() : clone entries from y */
    1463             : static GEN
    1464          56 : addp_union(GEN x, GEN y)
    1465             : {
    1466          56 :   long i, j, k, lx = lg(x), ly = lg(y);
    1467          56 :   GEN z = cgetg(lx + ly - 1, t_VEC);
    1468          56 :   i = j = k = 1;
    1469          84 :   while (i<lx && j<ly)
    1470             :   {
    1471          28 :     int s = cmpii(gel(x,i), gel(y,j));
    1472          28 :     if (s < 0)
    1473          21 :       gel(z,k++) = gel(x,i++);
    1474           7 :     else if (s > 0)
    1475           0 :       gel(z,k++) = gclone(gel(y,j++));
    1476             :     else {
    1477           7 :       gel(z,k++) = gel(x,i++);
    1478           7 :       j++;
    1479             :     }
    1480             :   }
    1481          56 :   while (i<lx) gel(z,k++) = gel(x,i++);
    1482         147 :   while (j<ly) gel(z,k++) = gclone(gel(y,j++));
    1483          56 :   setlg(z, k); return z;
    1484             : }
    1485             : 
    1486             : /* p is NULL, or a single element or a row vector with "primes" to add to
    1487             :  * prime table. */
    1488             : static GEN
    1489         189 : addp(GEN *T, GEN p)
    1490             : {
    1491         189 :   pari_sp av = avma;
    1492             :   long i, l;
    1493             :   GEN v;
    1494             : 
    1495         189 :   if (!p || lg(p) == 1) return *T;
    1496          70 :   if (!is_vec_t(typ(p))) p = mkvec(p);
    1497             : 
    1498          70 :   RgV_check_ZV(p, "addprimes");
    1499          63 :   v = gen_indexsort_uniq(p, (void*)&cmpii, &cmp_nodata);
    1500          63 :   p = vecpermute(p, v);
    1501          63 :   if (abscmpiu(gel(p,1), 2) < 0) pari_err_DOMAIN("addprimes", "p", "<", gen_2,p);
    1502          56 :   p = addp_union(*T, p);
    1503          56 :   l = lg(p);
    1504          56 :   if (l != lg(*T))
    1505             :   {
    1506          56 :     GEN old = *T, t = cgetg_block(l, t_VEC);
    1507         175 :     for (i = 1; i < l; i++) gel(t,i) = gel(p,i);
    1508          56 :     *T = t; gunclone(old);
    1509             :   }
    1510          56 :   set_avma(av); return *T;
    1511             : }
    1512             : GEN
    1513         189 : addprimes(GEN p) { return addp(&primetab, p); }
    1514             : 
    1515             : static GEN
    1516          28 : rmprimes(GEN T, GEN prime)
    1517             : {
    1518             :   long i,tx;
    1519             : 
    1520          28 :   if (!prime) return T;
    1521          28 :   tx = typ(prime);
    1522          28 :   if (is_vec_t(tx))
    1523             :   {
    1524          14 :     if (prime == T)
    1525             :     {
    1526          14 :       for (i=1; i < lg(prime); i++) gunclone(gel(prime,i));
    1527           7 :       setlg(prime, 1);
    1528             :     }
    1529             :     else
    1530             :     {
    1531          21 :       for (i=1; i < lg(prime); i++) rmprime(T, gel(prime,i));
    1532             :     }
    1533          14 :     return T;
    1534             :   }
    1535          14 :   rmprime(T, prime); return T;
    1536             : }
    1537             : GEN
    1538          28 : removeprimes(GEN prime) { return rmprimes(primetab, prime); }

Generated by: LCOV version 1.13