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 - language - forprime.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.16.2 lcov report (development 29419-8afb0ed749) Lines: 439 501 87.6 %
Date: 2024-07-02 09:03:41 Functions: 36 38 94.7 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2016  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; either version 2 of the License, or (at your option) any later
       8             : version. It is distributed in the hope that it will be useful, but WITHOUT
       9             : ANY WARRANTY WHATSOEVER.
      10             : 
      11             : Check the License for details. You should have received a copy of it, along
      12             : with the package; see the file 'COPYING'. If not, write to the Free Software
      13             : Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. */
      14             : 
      15             : #include "pari.h"
      16             : #include "paripriv.h"
      17             : 
      18             : /**********************************************************************/
      19             : /***                                                                ***/
      20             : /***                     Public prime table                         ***/
      21             : /***                                                                ***/
      22             : /**********************************************************************/
      23             : 
      24             : static ulong _maxprimelim = 0;
      25             : static GEN _prodprimes,_prodprimes_addr;
      26             : typedef unsigned char *byteptr;
      27             : 
      28             : /* Build/Rebuild table of prime differences. The actual work is done by the
      29             :  * following two subroutines;  the user entry point is the function
      30             :  * initprimes() below; initprimes1() is the basecase, called when
      31             :  * maxnum (size) is moderate. Must be called after pari_init_stack() )*/
      32             : static void
      33        1841 : initprimes1(ulong size, long *lenp, pari_prime *p1)
      34             : {
      35        1841 :   pari_sp av = avma;
      36             :   long k;
      37        1841 :   byteptr q, r, s, p = (byteptr)stack_calloc(size+2), fin = p + size;
      38             :   pari_prime *re;
      39             : 
      40       22092 :   for (r=q=p,k=1; r<=fin; )
      41             :   {
      42       33138 :     do { r+=k; k+=2; r+=k; } while (*++q);
      43      889203 :     for (s=r; s<=fin; s+=k) *s = 1;
      44             :   }
      45        1841 :   re = p1; *re++ = 2; *re++ = 3; /* 2 and 3 */
      46        1841 :   for (s=q=p+1; ; s=q)
      47             :   {
      48      940751 :     do q++; while (*q);
      49      314811 :     if (q > fin) break;
      50      312970 :     *re++ = (pari_prime) 2*(q-p)+1;
      51             :   }
      52        1841 :   *re++ = 0;
      53        1841 :   *lenp = re - p1;
      54        1841 :   set_avma(av);
      55        1841 : }
      56             : 
      57             : /*  Timing in ms (Athlon/850; reports 512K of secondary cache; looks
      58             :     like there is 64K of quickier cache too).
      59             : 
      60             :       arena|    30m     100m    300m    1000m    2000m  <-- primelimit
      61             :       =================================================
      62             :       16K       1.1053  1.1407  1.2589  1.4368   1.6086
      63             :       24K       1.0000  1.0625  1.1320  1.2443   1.3095
      64             :       32K       1.0000  1.0469  1.0761  1.1336   1.1776
      65             :       48K       1.0000  1.0000  1.0254  1.0445   1.0546
      66             :       50K       1.0000  1.0000  1.0152  1.0345   1.0464
      67             :       52K       1.0000  1.0000  1.0203  1.0273   1.0362
      68             :       54K       1.0000  1.0000  1.0812  1.0216   1.0281
      69             :       56K       1.0526  1.0000  1.0051  1.0144   1.0205
      70             :       58K       1.0000  1.0000  1.0000  1.0086   1.0123
      71             :       60K       0.9473  0.9844  1.0051  1.0014   1.0055
      72             :       62K       1.0000  0.9844  0.9949  0.9971   0.9993
      73             :       64K       1.0000  1.0000  1.0000  1.0000   1.0000
      74             :       66K       1.2632  1.2187  1.2183  1.2055   1.1953
      75             :       68K       1.4211  1.4844  1.4721  1.4425   1.4188
      76             :       70K       1.7368  1.7188  1.7107  1.6767   1.6421
      77             :       72K       1.9474  1.9531  1.9594  1.9023   1.8573
      78             :       74K       2.2105  2.1875  2.1827  2.1207   2.0650
      79             :       76K       2.4211  2.4219  2.4010  2.3305   2.2644
      80             :       78K       2.5789  2.6250  2.6091  2.5330   2.4571
      81             :       80K       2.8421  2.8125  2.8223  2.7213   2.6380
      82             :       84K       3.1053  3.1875  3.1776  3.0819   2.9802
      83             :       88K       3.5263  3.5312  3.5228  3.4124   3.2992
      84             :       92K       3.7895  3.8438  3.8375  3.7213   3.5971
      85             :       96K       4.0000  4.1093  4.1218  3.9986   3.9659
      86             :       112K      4.3684  4.5781  4.5787  4.4583   4.6115
      87             :       128K      4.7368  4.8750  4.9188  4.8075   4.8997
      88             :       192K      5.5263  5.7188  5.8020  5.6911   5.7064
      89             :       256K      6.0000  6.2187  6.3045  6.1954   6.1033
      90             :       384K      6.7368  6.9531  7.0405  6.9181   6.7912
      91             :       512K      7.3158  7.5156  7.6294  7.5000   7.4654
      92             :       768K      9.1579  9.4531  9.6395  9.5014   9.1075
      93             :       1024K    10.368  10.7497 10.9999 10.878   10.8201
      94             :       1536K    12.579  13.3124 13.7660 13.747   13.4739
      95             :       2048K    13.737  14.4839 15.0509 15.151   15.1282
      96             :       3076K    14.789  15.5780 16.2993 16.513   16.3365
      97             : 
      98             :     Now the same number relative to the model
      99             : 
     100             :     (1 + 0.36*sqrt(primelimit)/arena) * (arena <= 64 ? 1.05 : (arena-64)**0.38)
     101             : 
     102             :      [SLOW2_IN_ROOTS = 0.36, ALPHA = 0.38]
     103             : 
     104             :       arena|    30m     100m    300m    1000m    2000m  <-- primelimit
     105             :       =================================================
     106             :         16K    1.014    0.9835  0.9942  0.9889  1.004
     107             :         24K    0.9526   0.9758  0.9861  0.9942  0.981
     108             :         32K    0.971    0.9939  0.9884  0.9849  0.9806
     109             :         48K    0.9902   0.9825  0.996   0.9945  0.9885
     110             :         50K    0.9917   0.9853  0.9906  0.9926  0.9907
     111             :         52K    0.9932   0.9878  0.9999  0.9928  0.9903
     112             :         54K    0.9945   0.9902  1.064   0.9939  0.9913
     113             :         56K    1.048    0.9924  0.9925  0.993   0.9921
     114             :         58K    0.9969   0.9945  0.9909  0.9932  0.9918
     115             :         60K    0.9455   0.9809  0.9992  0.9915  0.9923
     116             :         62K    0.9991   0.9827  0.9921  0.9924  0.9929
     117             :         64K    1        1       1       1       1
     118             :         66K    1.02     0.9849  0.9857  0.9772  0.9704
     119             :         68K    0.8827   0.9232  0.9176  0.9025  0.8903
     120             :         70K    0.9255   0.9177  0.9162  0.9029  0.8881
     121             :         72K    0.9309   0.936   0.9429  0.9219  0.9052
     122             :         74K    0.9715   0.9644  0.967   0.9477  0.9292
     123             :         76K    0.9935   0.9975  0.9946  0.9751  0.9552
     124             :         78K    0.9987   1.021   1.021   1.003   0.9819
     125             :         80K    1.047    1.041   1.052   1.027   1.006
     126             :         84K    1.052    1.086   1.092   1.075   1.053
     127             :         88K    1.116    1.125   1.133   1.117   1.096
     128             :         92K    1.132    1.156   1.167   1.155   1.134
     129             :         96K    1.137    1.177   1.195   1.185   1.196
     130             :        112K    1.067    1.13    1.148   1.15    1.217
     131             :        128K    1.04     1.083   1.113   1.124   1.178
     132             :        192K    0.9368   0.985   1.025   1.051   1.095
     133             :        256K    0.8741   0.9224  0.9619  0.995   1.024
     134             :        384K    0.8103   0.8533  0.8917  0.9282  0.9568
     135             :        512K    0.7753   0.8135  0.8537  0.892   0.935
     136             :        768K    0.8184   0.8638  0.9121  0.9586  0.9705
     137             :       1024K    0.8241   0.8741  0.927   0.979   1.03
     138             :       1536K    0.8505   0.9212  0.9882  1.056   1.096
     139             :       2048K    0.8294   0.8954  0.9655  1.041   1.102
     140             : 
     141             : */
     142             : 
     143             : #ifndef SLOW2_IN_ROOTS
     144             :   /* SLOW2_IN_ROOTS below 3: some slowdown starts to be noticable
     145             :    * when things fit into the cache on Sparc.
     146             :    * The choice of 2.6 gives a slowdown of 1-2% on UltraSparcII,
     147             :    * but makes calculations for "maximum" of 436273009
     148             :    * fit into 256K cache (still common for some architectures).
     149             :    *
     150             :    * One may change it when small caches become uncommon, but the gain
     151             :    * is not going to be very noticable... */
     152             : #  ifdef i386           /* gcc defines this? */
     153             : #    define SLOW2_IN_ROOTS      0.36
     154             : #  else
     155             : #    define SLOW2_IN_ROOTS      2.6
     156             : #  endif
     157             : #endif
     158             : #ifndef CACHE_ARENA
     159             : #  ifdef i386           /* gcc defines this? */
     160             :    /* Due to smaller SLOW2_IN_ROOTS, smaller arena is OK; fit L1 cache */
     161             : #    define CACHE_ARENA (63 * 1024UL) /* No slowdown even with 64K L1 cache */
     162             : #  else
     163             : #    define CACHE_ARENA (200 * 1024UL) /* No slowdown even with 256K L2 cache */
     164             : #  endif
     165             : #endif
     166             : 
     167             : #define CACHE_ALPHA     (0.38)          /* Cache performance model parameter */
     168             : #define CACHE_CUTOFF    (0.018)         /* Cache performance not smooth here */
     169             : 
     170             : static double slow2_in_roots = SLOW2_IN_ROOTS;
     171             : 
     172             : typedef struct {
     173             :     ulong arena;
     174             :     double power;
     175             :     double cutoff;
     176             : } cache_model_t;
     177             : 
     178             : static cache_model_t cache_model = { CACHE_ARENA, CACHE_ALPHA, CACHE_CUTOFF };
     179             : 
     180             : /* Assume that some calculation requires a chunk of memory to be
     181             :    accessed often in more or less random fashion (as in sieving).
     182             :    Assume that the calculation can be done in steps by subdividing the
     183             :    chunk into smaller subchunks (arenas) and treating them
     184             :    separately.  Assume that the overhead of subdivision is equivalent
     185             :    to the number of arenas.
     186             : 
     187             :    Find an optimal size of the arena taking into account the overhead
     188             :    of subdivision, and the overhead of arena not fitting into the
     189             :    cache.  Assume that arenas of size slow2_in_roots slows down the
     190             :    calculation 2x (comparing to very big arenas; when cache hits do
     191             :    not matter).  Since cache performance varies wildly with
     192             :    architecture, load, and wheather (especially with cache coloring
     193             :    enabled), use an idealized cache model based on benchmarks above.
     194             : 
     195             :    Assume that an independent region of FIXED_TO_CACHE bytes is accessed
     196             :    very often concurrently with the arena access.
     197             :  */
     198             : static ulong
     199        1841 : good_arena_size(ulong slow2_size, ulong total, ulong fixed_to_cache,
     200             :                 cache_model_t *cache_model)
     201             : {
     202        1841 :   ulong asize, cache_arena = cache_model->arena;
     203             :   double Xmin, Xmax, A, B, C1, C2, D, V;
     204        1841 :   double alpha = cache_model->power, cut_off = cache_model->cutoff;
     205             : 
     206             :   /* Estimated relative slowdown,
     207             :      with overhead = max((fixed_to_cache+arena)/cache_arena - 1, 0):
     208             : 
     209             :      1 + slow2_size/arena due to initialization overhead;
     210             : 
     211             :      max(1, 4.63 * overhead^0.38 ) due to footprint > cache size.
     212             : 
     213             :      [The latter is hard to substantiate theoretically, but this
     214             :      function describes benchmarks pretty close; it does not hurt that
     215             :      one can minimize it explicitly too ;-).  The switch between
     216             :      different choices of max() happens when overhead=0.018.]
     217             : 
     218             :      Thus the problem is minimizing (1 + slow2_size/arena)*overhead**0.29.
     219             :      This boils down to F=((X+A)/(X+B))X^alpha, X=overhead,
     220             :      B = (1 - fixed_to_cache/cache_arena), A = B + slow2_size/cache_arena,
     221             :      alpha = 0.38, and X>=0.018, X>-B.
     222             : 
     223             :      We need to find the rightmost root of (X+A)*(X+B) - alpha(A-B)X to the
     224             :      right of 0.018 (if such exists and is below Xmax).  Then we manually
     225             :      check the remaining region [0, 0.018].
     226             : 
     227             :      Since we cannot trust the purely-experimental cache-hit slowdown
     228             :      function, as a sanity check always prefer fitting into the
     229             :      cache (or "almost fitting") if F-law predicts that the larger
     230             :      value of the arena provides less than 10% speedup.
     231             :    */
     232             : 
     233             :   /* The simplest case: we fit into cache */
     234        1841 :   asize = cache_arena - fixed_to_cache;
     235        1841 :   if (total <= asize) return total;
     236             :   /* The simple case: fitting into cache doesn't slow us down more than 10% */
     237        1841 :   if (asize > 10 * slow2_size) return asize;
     238             :   /* Slowdown of not fitting into cache is significant.  Try to optimize.
     239             :      Do not be afraid to spend some time on optimization - in trivial
     240             :      cases we do not reach this point; any gain we get should
     241             :      compensate the time spent on optimization.  */
     242             : 
     243           0 :   B = (1 - ((double)fixed_to_cache)/cache_arena);
     244           0 :   A = B + ((double)slow2_size)/cache_arena;
     245           0 :   C2 = A*B;
     246           0 :   C1 = (A + B - 1/alpha*(A - B))/2;
     247           0 :   D = C1*C1 - C2;
     248           0 :   if (D > 0)
     249           0 :     V = cut_off*cut_off + 2*C1*cut_off + C2; /* Value at CUT_OFF */
     250             :   else
     251           0 :     V = 0; /* Peacify the warning */
     252           0 :   Xmin = cut_off;
     253           0 :   Xmax = ((double)total - fixed_to_cache)/cache_arena; /* Two candidates */
     254             : 
     255           0 :   if ( D <= 0 || (V >= 0 && C1 + cut_off >= 0) ) /* slowdown increasing */
     256           0 :     Xmax = cut_off; /* Only one candidate */
     257           0 :   else if (V >= 0 && /* slowdown concave down */
     258           0 :            ((Xmax + C1) <= 0 || (Xmax*Xmax + 2*C1*Xmax + C2) <= 0))
     259             :       /* DO NOTHING */;  /* Keep both candidates */
     260           0 :   else if (V <= 0 && (Xmax*Xmax + 2*C1*Xmax + C2) <= 0) /*slowdown decreasing*/
     261           0 :       Xmin = cut_off; /* Only one candidate */
     262             :   else /* Now we know: 2 roots, the largest is in CUT_OFF..Xmax */
     263           0 :       Xmax = sqrt(D) - C1;
     264           0 :   if (Xmax != Xmin) { /* Xmin == CUT_OFF; Check which one is better */
     265           0 :     double v1 = (cut_off + A)/(cut_off + B);
     266           0 :     double v2 = 2.33 * (Xmax + A)/(Xmax + B) * pow(Xmax, alpha);
     267             : 
     268           0 :     if (1.1 * v2 >= v1) /* Prefer fitting into the cache if slowdown < 10% */
     269           0 :       V = v1;
     270             :     else
     271           0 :     { Xmin = Xmax; V = v2; }
     272           0 :   } else if (B > 0) /* We need V */
     273           0 :     V = 2.33 * (Xmin + A)/(Xmin + B) * pow(Xmin, alpha);
     274           0 :   if (B > 0 && 1.1 * V > A/B)  /* Now Xmin is the minumum.  Compare with 0 */
     275           0 :     Xmin = 0;
     276             : 
     277           0 :   asize = (ulong)((1 + Xmin)*cache_arena - fixed_to_cache);
     278           0 :   if (asize > total) asize = total; /* May happen due to approximations */
     279           0 :   return asize;
     280             : }
     281             : 
     282             : /* Use as in
     283             :     install(set_optimize,lLDG)          \\ Through some M too?
     284             :     set_optimize(2,1) \\ disable dependence on limit
     285             :     \\ 1: how much cache usable, 2: slowdown of setup, 3: alpha, 4: cutoff
     286             :     \\ 2,3,4 are in units of 0.001
     287             : 
     288             :     { time_primes_arena(ar,limit) =     \\ ar = arena size in K
     289             :         set_optimize(1,floor(ar*1024));
     290             :         default(primelimit, 200 000);   \\ 100000 results in *larger* malloc()!
     291             :         gettime;
     292             :         default(primelimit, floor(limit));
     293             :         if(ar >= 1, ar=floor(ar));
     294             :         print("arena "ar"K => "gettime"ms");
     295             :     }
     296             : */
     297             : long
     298           0 : set_optimize(long what, GEN g)
     299             : {
     300           0 :   long ret = 0;
     301             : 
     302           0 :   switch (what) {
     303           0 :   case 1:
     304           0 :     ret = (long)cache_model.arena;
     305           0 :     break;
     306           0 :   case 2:
     307           0 :     ret = (long)(slow2_in_roots * 1000);
     308           0 :     break;
     309           0 :   case 3:
     310           0 :     ret = (long)(cache_model.power * 1000);
     311           0 :     break;
     312           0 :   case 4:
     313           0 :     ret = (long)(cache_model.cutoff * 1000);
     314           0 :     break;
     315           0 :   default:
     316           0 :     pari_err_BUG("set_optimize");
     317           0 :     break;
     318             :   }
     319           0 :   if (g != NULL) {
     320           0 :     ulong val = itou(g);
     321             : 
     322           0 :     switch (what) {
     323           0 :     case 1: cache_model.arena = val; break;
     324           0 :     case 2: slow2_in_roots     = (double)val / 1000.; break;
     325           0 :     case 3: cache_model.power  = (double)val / 1000.; break;
     326           0 :     case 4: cache_model.cutoff = (double)val / 1000.; break;
     327             :     }
     328             :   }
     329           0 :   return ret;
     330             : }
     331             : 
     332             : /* s is odd; prime (starting from 3 = known_primes[2]), terminated by a 0 byte.
     333             :  * Checks n odd numbers starting at 'start', setting bytes to 0 (composite)
     334             :  * or 1 (prime), starting at data */
     335             : static void
     336        7053 : sieve_chunk(pari_prime *known_primes, ulong s, byteptr data, ulong n)
     337             : {
     338        7053 :   ulong p, cnt = n-1, start = s;
     339             :   pari_prime *q;
     340             : 
     341        7053 :   memset(data, 0, n);
     342        7053 :   start >>= 1;  /* (start - 1)/2 */
     343        7053 :   start += n; /* Corresponds to the end */
     344             :   /* data corresponds to start, q runs over primediffs */
     345     1013244 :   for (q = known_primes + 1, p = 3; p; p = *++q)
     346             :   { /* first odd number >= start > p and divisible by p
     347             :        = last odd number <= start + 2p - 2 and 0 (mod p)
     348             :        = p + last number <= start + p - 2 and 0 (mod 2p)
     349             :        = p + start+p-2 - (start+p-2) % 2p
     350             :        = start + 2(p - 1 - ((start-1)/2 + (p-1)/2) % p). */
     351     1006191 :     long off = cnt - ((start+(p>>1)) % p);
     352  1607084161 :     while (off >= 0) { data[off] = 1; off -= p; }
     353             :   }
     354        7053 : }
     355             : 
     356             : static void
     357        1841 : set_prodprimes(void)
     358             : {
     359        1841 :   pari_sp ltop = avma, av;
     360        1841 :   ulong b = 1UL << 8, M = minuu(maxprime(), GP_DATA->factorlimit);
     361        1841 :   GEN W, w, v = primes_interval_zv(3, M);
     362        1841 :   long s, u, j, jold, l = lg(v);
     363             : 
     364        1841 :   W = cgetg(64+1, t_VEC);
     365   151008025 :   for (jold = j = u = 1; j < l; j++)
     366   151006184 :     if (uel(v,j) >= b)
     367             :     {
     368       23933 :       long lw = (j == l-1? l: j) - jold + 1;
     369       23933 :       w = v+jold-1; w[0] = evaltyp(t_VECSMALL) | _evallg(lw);
     370       23933 :       gel(W,u++) = zv_prod_Z(w); /* p_jold ... p_{j-1} */
     371       23933 :       jold = j; b *= 2;
     372       23933 :       if (b > M) b = M; /* truncate last run */
     373             :     }
     374        1841 :   setlg(W, u);
     375       23933 :   for (j = 2; j < u; j++) gel(W,j) = mulii(gel(W,j-1), gel(W,j));
     376        1841 :   s = gsizeword(W);
     377        1841 :   w = (GEN)pari_malloc(s*sizeof(long));
     378        1841 :   av = (pari_sp)(w + s);
     379        1841 :   _prodprimes_addr = w;
     380        1841 :   _prodprimes = gcopy_avma(W, &av);
     381        1841 :   set_avma(ltop);
     382        1841 : }
     383             : 
     384             : static void
     385        1841 : initprimes0(ulong maxnum, long *lenp, pari_prime *p1)
     386             : {
     387        1841 :   pari_sp av = avma, bot = pari_mainstack->bot;
     388             :   long alloced, psize;
     389             :   byteptr q, end, p;
     390             :   ulong remains, curlow, rootnum, asize, prime_above, last;
     391             :   pari_prime *end1, *curdiff, *p_prime_above;
     392             : 
     393        1841 :   if (!odd(maxnum)) maxnum--; /* make it odd. */
     394             :   /* base case */
     395        1841 :   if (maxnum < 1ul<<17) { initprimes1(maxnum>>1, lenp, p1); return; }
     396             : 
     397             :   /* Checked to be enough up to 40e6, attained at 155893 */
     398        1841 :   rootnum = usqrt(maxnum) | 1;
     399        1841 :   initprimes1(rootnum>>1, &psize, p1);
     400        1841 :   last = rootnum;
     401        1841 :   end1 = p1 + psize - 1;
     402        1841 :   remains = (maxnum - last) >> 1; /* number of odd numbers to check */
     403             :   /* we access primes array of psize too; but we access it consecutively,
     404             :    * thus we do not include it in fixed_to_cache */
     405        1841 :   asize = good_arena_size((ulong)(rootnum * slow2_in_roots), remains+1, 0,
     406             :                           &cache_model) - 1;
     407             :   /* enough room on the stack ? */
     408        1841 :   alloced = (((byteptr)avma) <= ((byteptr)bot) + asize);
     409        1841 :   p = (byteptr)(alloced? pari_malloc(asize+1): stack_malloc(asize+1));
     410        1841 :   end = p + asize; /* the 0 sentinel goes at end. */
     411        1841 :   curlow = last + 2; /* First candidate: know primes up to last (odd). */
     412        1841 :   curdiff = end1;
     413             : 
     414             :   /* During each iteration p..end-1 represents a range of odd
     415             :      numbers.   */
     416        1841 :   p_prime_above = p1 + 2;
     417        1841 :   prime_above = 3;
     418        8894 :   while (remains)
     419             :   { /* cycle over arenas; performance not crucial */
     420             :     pari_prime was_delta;
     421        7053 :     if (asize > remains) { asize = remains; end = p + asize; }
     422             :     /* Fake the upper limit appropriate for the given arena */
     423      321864 :     while (prime_above*prime_above <= curlow + (asize << 1) && *p_prime_above)
     424      314811 :       prime_above = *p_prime_above++;
     425        7053 :     was_delta = *p_prime_above;
     426        7053 :     *p_prime_above = 0; /* sentinel for sieve_chunk */
     427        7053 :     sieve_chunk(p1, curlow, p, asize);
     428        7053 :     *p_prime_above = was_delta; /* restore */
     429             : 
     430        7053 :     p[asize] = 0; /* sentinel */
     431        7053 :     for (q = p; ; q++)
     432             :     { /* q runs over addresses corresponding to primes */
     433   964278669 :       while (*q) q++; /* use sentinel at end */
     434   150698426 :       if (q >= end) break;
     435   150691373 :       *curdiff++ = (pari_prime) 2*(q-p) + curlow;
     436             :     }
     437        7053 :     remains -= asize;
     438        7053 :     curlow += (asize<<1);
     439             :   }
     440        1841 :   *curdiff++ = 0; /* sentinel */
     441        1841 :   *lenp = curdiff - p1;
     442        1841 :   if (alloced) pari_free(p); else set_avma(av);
     443             : }
     444             : 
     445             : ulong
     446    44322027 : maxprime(void) { return pari_PRIMES? pari_PRIMES[pari_PRIMES[0]]: 0; }
     447             : ulong
     448    44282827 : maxprimelim(void) { return pari_PRIMES? _maxprimelim: 0; }
     449             : ulong
     450         196 : maxprimeN(void) { return pari_PRIMES? pari_PRIMES[0]: 0; }
     451             : GEN
     452     2686089 : prodprimes(void) { return pari_PRIMES? _prodprimes: NULL; }
     453             : void
     454           0 : maxprime_check(ulong c) { if (maxprime() < c) pari_err_MAXPRIME(c); }
     455             : 
     456             : static pari_prime*
     457        1841 : initprimes(ulong maxnum)
     458             : {
     459             :   pari_prime *t;
     460             :   long len;
     461             :   ulong N;
     462        1841 :   if (maxnum < 65537)
     463             :   {
     464           0 :     maxnum = 65537;
     465           0 :     N = 6543;
     466             :   }
     467             :   else
     468        1841 :     N = (long) ceil(primepi_upper_bound((double)maxnum));
     469        1841 :   t = (pari_prime*) pari_malloc(sizeof(*t) * (N+2));
     470        1841 :   initprimes0(maxnum, &len, t+1); t[0] = (pari_prime)(len-1);
     471        1841 :   _maxprimelim = maxnum;
     472        1841 :   return (pari_prime*) pari_realloc(t, sizeof(*t) * (len+1));
     473             : }
     474             : 
     475             : void
     476        1841 : initprimetable(ulong maxnum)
     477             : {
     478        1841 :   pari_prime *old = pari_PRIMES;
     479        1841 :   pari_PRIMES = initprimes(maxnum);
     480        1841 :   if (old) free(old);
     481        1841 :   set_prodprimes();
     482        1841 : }
     483             : 
     484             : /**********************************************************************/
     485             : /***                                                                ***/
     486             : /***                     forprime                                   ***/
     487             : /***                                                                ***/
     488             : /**********************************************************************/
     489             : 
     490             : /* return good chunk size for sieve, 16 | chunk + 2 */
     491             : static ulong
     492     7308200 : optimize_chunk(ulong a, ulong b)
     493             : {
     494             :   /* TODO: Optimize size (surely < 512k to stay in L2 cache, but not so large
     495             :    * as to force recalculating too often). */
     496     7308200 :   ulong chunk = 0x80000UL;
     497     7308200 :   ulong tmp = (b - a) / chunk + 1;
     498             : 
     499     7308200 :   if (tmp == 1)
     500           7 :     chunk = b - a + 16;
     501             :   else
     502     7308193 :     chunk = (b - a) / tmp + 15;
     503             :   /* ensure 16 | chunk + 2 */
     504     7308200 :   return (((chunk + 2)>>4)<<4) - 2;
     505             : }
     506             : static void
     507     7308187 : sieve_init(forprime_t *T, ulong a, ulong b)
     508             : {
     509     7308187 :   T->sieveb = b;
     510     7308187 :   T->chunk = optimize_chunk(a, b);
     511             :   /* >> 1 [only odds] + 3 [convert from bits to bytes] */
     512     7308209 :   T->isieve = (unsigned char*)stack_malloc(((T->chunk+2) >> 4) + 1);
     513     7308204 :   T->cache[0] = 0;
     514     7308204 :   T->a = a;
     515     7308204 :   T->end = minuu(a + T->chunk, b);
     516     7308197 :   T->pos = T->maxpos = 0;
     517     7308197 : }
     518             : 
     519             : enum {PRST_none, PRST_diffptr, PRST_sieve, PRST_unextprime, PRST_nextprime};
     520             : 
     521             : static void
     522    12844020 : u_forprime_set_prime_table(forprime_t *T, ulong a)
     523             : {
     524    12844020 :   T->strategy = PRST_diffptr;
     525    12844020 :   if (a < 3)
     526             :   {
     527     1517818 :     T->p = 0;
     528     1517818 :     T->n = 0;
     529             :   }
     530             :   else
     531             :   {
     532    11326202 :     long n = PRIMES_search(a - 1);
     533    11329079 :     if (n < 0) n = - n - 1;
     534    11329079 :     T->n = n;
     535    11329079 :     T->p = pari_PRIMES[n];
     536             :   }
     537    12846897 : }
     538             : 
     539             : /* Set p so that p + q the smallest integer = c (mod q) and > original p.
     540             :  * Assume 0 < c < q. */
     541             : static void
     542      101975 : arith_set(forprime_t *T)
     543             : {
     544      101975 :   ulong r = T->p % T->q; /* 0 <= r <= min(p, q-1) */
     545      101975 :   pari_sp av = avma;
     546      101975 :   GEN d = adduu(T->p - r, T->c); /* = c mod q */
     547      101975 :   if (T->c > r) d = subiu(d, T->q);
     548             :   /* d = c mod q,  d = c > r? p-r+c-q: p-r+c, so that
     549             :    *  d <= p  and  d+q = c>r? p-r+c  : p-r+c+q > p */
     550      101975 :   if (signe(d) <= 0)
     551             :   {
     552          20 :     T->p = 0;
     553          20 :     T->strategy = PRST_nextprime;
     554          20 :     affii(d, T->pp);
     555             :   }
     556             :   else
     557      101955 :     T->p = itou_or_0(d);
     558      101975 :   set_avma(av);
     559      101975 : }
     560             : 
     561             : /* run through primes in arithmetic progression = c (mod q) */
     562             : static int
     563    26355090 : u_forprime_sieve_arith_init(forprime_t *T, struct pari_sieve *psieve,
     564             :                             ulong a, ulong b, ulong c, ulong q)
     565             : {
     566             :   ulong maxp, maxp2;
     567    26355090 :   if (!odd(b) && b > 2) b--;
     568    26353419 :   if (a > b || b < 2)
     569             :   {
     570      882916 :     T->strategy = PRST_diffptr; /* paranoia */
     571      882916 :     T->p = 0; /* empty */
     572      882916 :     T->b = 0; /* empty */
     573      882916 :     T->n = 0;
     574      882916 :     return 0;
     575             :   }
     576    25470503 :   maxp = maxprime();
     577    25469673 :   if (q != 1)
     578             :   {
     579             :     ulong D;
     580      587449 :     c %= q; D = ugcd(c, q);
     581      587426 :     if (D != 1) { a = maxuu(a,D); b = minuu(b,D); }
     582      587426 :     if (odd(q) && (a > 2 || c != 2))
     583             :     { /* only *odd* primes. If a <= c = 2, then p = 2 must be included :-( */
     584      508757 :       if (!odd(c)) c += q;
     585      510074 :       q <<= 1;
     586             :     }
     587             :   }
     588    25470930 :   T->q = q;
     589    25470930 :   T->c = c;
     590    25470930 :   T->strategy = PRST_none; /* unknown */
     591    25470930 :   T->psieve = psieve; /* unused for now */
     592    25470930 :   T->isieve = NULL; /* unused for now */
     593    25470930 :   T->b = b;
     594    25470930 :   if (maxp >= b) { /* [a,b] \subset prime table */
     595    10124558 :     u_forprime_set_prime_table(T, a);
     596    10127562 :     return 1;
     597             :   }
     598             :   /* b > maxp */
     599    15346372 :   if (a >= maxp)
     600             :   {
     601    12627646 :     T->p = a - 1;
     602    12627646 :     if (T->q != 1) arith_set(T);
     603             :   }
     604             :   else
     605     2718726 :     u_forprime_set_prime_table(T, a);
     606             : 
     607    15346348 :   maxp2 = (maxp & HIGHMASK)? 0 : maxp*maxp;
     608             :   /* FIXME: should sieve as well if q != 1, adapt sieve code */
     609    15346348 :   if (q != 1 || (maxp2 && maxp2 <= a)
     610     7309390 :              || T->b - maxuu(a,maxp) < maxp / expu(b))
     611     8038264 :   { if (T->strategy==PRST_none) T->strategy = PRST_unextprime; }
     612             :   else
     613             :   { /* worth sieving */
     614             : #ifdef LONG_IS_64BIT
     615     5195445 :     const ulong UPRIME_MAX = 18446744073709551557UL;
     616             : #else
     617     2112750 :     const ulong UPRIME_MAX = 4294967291UL;
     618             : #endif
     619             :     ulong sieveb;
     620     7308195 :     if (b > UPRIME_MAX) b = UPRIME_MAX;
     621     7308195 :     sieveb = b;
     622     7308195 :     if (maxp2 && maxp2 < b) sieveb = maxp2;
     623     7308195 :     if (T->strategy==PRST_none) T->strategy = PRST_sieve;
     624     7308195 :     sieve_init(T, maxuu(maxp+2, a), sieveb);
     625             :   }
     626    15346335 :   return 1;
     627             : }
     628             : 
     629             : int
     630    21275153 : u_forprime_arith_init(forprime_t *T, ulong a, ulong b, ulong c, ulong q)
     631    21275153 : { return u_forprime_sieve_arith_init(T, NULL, a, b, c, q); }
     632             : 
     633             : /* will run through primes in [a,b] */
     634             : int
     635    20684272 : u_forprime_init(forprime_t *T, ulong a, ulong b)
     636    20684272 : { return u_forprime_arith_init(T, a,b, 0,1); }
     637             : 
     638             : /* will run through primes in [a,b] */
     639             : static int
     640     5072861 : u_forprime_sieve_init(forprime_t *T, struct pari_sieve *s, ulong b)
     641     5072861 : { return u_forprime_sieve_arith_init(T, s, s->start, b, s->c, s->q); }
     642             : 
     643             : /* now only run through primes <= c; assume c <= b above */
     644             : void
     645          63 : u_forprime_restrict(forprime_t *T, ulong c) { T->b = c; }
     646             : 
     647             : /* b = NULL: loop forever */
     648             : int
     649        2182 : forprimestep_init(forprime_t *T, GEN a, GEN b, GEN q)
     650             : {
     651        2182 :   GEN c = NULL;
     652             :   long lb;
     653             : 
     654        2182 :   a = gceil(a); if (typ(a) != t_INT) pari_err_TYPE("forprime_init",a);
     655        2182 :   T->qq = NULL; T->q = 1; T->c = 0;
     656        2182 :   if (q)
     657             :   {
     658         133 :     switch(typ(q))
     659             :     {
     660          56 :       case t_INT:
     661          56 :         c = a; break;
     662          77 :       case t_INTMOD:
     663          77 :         c = gel(q,2); q = gel(q,1);
     664             :         /* first int >= initial a which is = c (mod q) */
     665          77 :         a = addii(a, modii(subii(c,a), q)); break;
     666           0 :       default: pari_err_TYPE("forprimestep_init",q);
     667             :     }
     668         133 :     if (signe(q) <= 0) pari_err_TYPE("forprimestep_init (q <= 0)",q);
     669         133 :     if (equali1(q)) c = q = NULL;
     670             :     else
     671             :     {
     672         133 :       GEN D = gcdii(c, q);
     673         133 :       if (!is_pm1(D))
     674             :       { /* at most one prime: c */
     675          42 :         if (cmpii(a, D) < 0) a = D;
     676          42 :         if (gcmp(b, D) > 0) b = D;
     677             :       }
     678         133 :       if ((T->q = itou_or_0(q)))
     679         125 :         T->c = umodiu(c, T->q);
     680             :       else
     681           8 :         T->qq = q;
     682             :     }
     683             :   }
     684        2182 :   if (signe(a) <= 0) a = q? modii(a, q): gen_1;
     685        2182 :   if (b && typ(b) != t_INFINITY)
     686             :   {
     687         789 :     b = gfloor(b);
     688         789 :     if (typ(b) != t_INT) pari_err_TYPE("forprime_init",b);
     689         789 :     if (signe(b) < 0 || cmpii(a,b) > 0)
     690             :     {
     691          21 :       T->strategy = PRST_nextprime; /* paranoia */
     692          21 :       T->bb = T->pp = gen_0; return 0;
     693             :     }
     694         768 :     lb = lgefint(b);
     695         768 :     T->bb = b;
     696             :   }
     697        1393 :   else if (!b || inf_get_sign(b) > 0)
     698             :   {
     699        1393 :     lb = lgefint(a) + 4;
     700        1393 :     T->bb = NULL;
     701             :   }
     702             :   else /* b == -oo */
     703             :   {
     704           0 :     T->strategy = PRST_nextprime; /* paranoia */
     705           0 :     T->bb = T->pp = gen_0; return 0;
     706             :   }
     707        2161 :   T->pp = cgeti(T->qq? maxuu(lb, lgefint(T->qq)): lb);
     708             :   /* a, b are positive integers, a <= b */
     709        2161 :   if (!T->qq && lgefint(a) == 3) /* lb == 3 implies b != NULL */
     710        2018 :     return u_forprime_arith_init(T, uel(a,2), lb == 3? uel(b,2): ULONG_MAX,
     711             :                                     T->c, T->q);
     712         143 :   T->strategy = PRST_nextprime;
     713         143 :   affii(T->qq? subii(a,T->qq): subiu(a,T->q), T->pp); return 1;
     714             : }
     715             : int
     716        1301 : forprime_init(forprime_t *T, GEN a, GEN b)
     717        1301 : { return forprimestep_init(T,a,b,NULL); }
     718             : 
     719             : /* assume a <= b <= maxprime()^2, a,b odd, sieve[n] corresponds to
     720             :  *   a+16*n, a+16*n+2, ..., a+16*n+14 (bits 0 to 7)
     721             :  * maxpos = index of last sieve cell.
     722             :  * b-a+2 must be divisible by 16 for use by u_forprime_next */
     723             : static void
     724        8954 : sieve_block(ulong a, ulong b, ulong maxpos, unsigned char* sieve)
     725             : {
     726        8954 :   ulong i, lim = usqrt(b), sz = (b-a) >> 1;
     727        8954 :   (void)memset(sieve, 0, maxpos+1);
     728        8954 :   for (i = 2;; i++)
     729    24344387 :   { /* p is odd */
     730    24353341 :     ulong k, r, p = pari_PRIMES[i]; /* starts at p = 3 */
     731    24353341 :     if (p > lim) break;
     732             : 
     733             :     /* solve a + 2k = 0 (mod p) */
     734    24344387 :     r = a % p;
     735    24344387 :     if (r == 0)
     736       16114 :       k = 0;
     737             :     else
     738             :     {
     739    24328273 :       k = p - r;
     740    24328273 :       if (odd(k)) k += p;
     741    24328273 :       k >>= 1;
     742             :     }
     743             :     /* m = a + 2k is the smallest odd m >= a, p | m */
     744             :     /* position n (corresponds to a+2n) is sieve[n>>3], bit n&7 */
     745  5705915113 :     while (k <= sz) { sieve[k>>3] |= 1 << (k&7); k += p; /* 2k += 2p */ }
     746             :   }
     747        8954 : }
     748             : 
     749             : static void
     750        1841 : pari_sieve_init(struct pari_sieve *s, ulong a, ulong b)
     751             : {
     752        1841 :   ulong maxpos= (b - a) >> 4;
     753        1841 :   s->start = a; s->end = b;
     754        1841 :   s->sieve = (unsigned char*) pari_malloc(maxpos+1);
     755        1841 :   s->c = 0; s->q = 1;
     756        1841 :   sieve_block(a, b, maxpos, s->sieve);
     757        1841 :   s->maxpos = maxpos; /* must be last in case of SIGINT */
     758        1841 : }
     759             : 
     760             : static struct pari_sieve pari_sieve_modular;
     761             : 
     762             : #ifdef LONG_IS_64BIT
     763             : #define PARI_MODULAR_BASE ((1UL<<((BITS_IN_LONG-2)>>1))+1)
     764             : #else
     765             : #define PARI_MODULAR_BASE ((1UL<<(BITS_IN_LONG-1))+1)
     766             : #endif
     767             : 
     768             : void
     769        1841 : pari_init_primes(ulong maxprime)
     770             : {
     771        1841 :   ulong a = PARI_MODULAR_BASE, b = a + (1UL<<20)-2;
     772        1841 :   initprimetable(maxprime);
     773        1841 :   pari_sieve_init(&pari_sieve_modular, a, b);
     774        1841 : }
     775             : 
     776             : void
     777        1841 : pari_close_primes(void)
     778             : {
     779        1841 :   if (pari_PRIMES)
     780             :   {
     781        1841 :     pari_free(pari_PRIMES);
     782        1841 :     pari_free(_prodprimes_addr);
     783             :   }
     784        1841 :   pari_free(pari_sieve_modular.sieve);
     785        1841 : }
     786             : 
     787             : void
     788     4464464 : init_modular_small(forprime_t *S)
     789             : {
     790             : #ifdef LONG_IS_64BIT
     791     3826397 :   u_forprime_sieve_init(S, &pari_sieve_modular, ULONG_MAX);
     792             : #else
     793      638067 :   ulong a = (1UL<<((BITS_IN_LONG-2)>>1))+1;
     794      638067 :   u_forprime_init(S, a, ULONG_MAX);
     795             : #endif
     796     4464447 : }
     797             : 
     798             : void
     799     8696089 : init_modular_big(forprime_t *S)
     800             : {
     801             : #ifdef LONG_IS_64BIT
     802     7449621 :   u_forprime_init(S, HIGHBIT + 1, ULONG_MAX);
     803             : #else
     804     1246468 :   u_forprime_sieve_init(S, &pari_sieve_modular, ULONG_MAX);
     805             : #endif
     806     8696048 : }
     807             : 
     808             : /* T->cache is a 0-terminated list of primes, return the first one and
     809             :  * remove it from list. Most of the time the list contains a single prime */
     810             : static ulong
     811   129453693 : shift_cache(forprime_t *T)
     812             : {
     813             :   long i;
     814   129453693 :   T->p = T->cache[0];
     815   173119335 :   for (i = 1;; i++)  /* remove one prime from cache */
     816   173119335 :     if (! (T->cache[i-1] = T->cache[i]) ) break;
     817   129453693 :   return T->p;
     818             : }
     819             : 
     820             : ulong
     821   202917968 : u_forprime_next(forprime_t *T)
     822             : {
     823   202917968 :   if (T->strategy == PRST_diffptr)
     824             :   {
     825             :     for(;;)
     826             :     {
     827   218319639 :       if (++T->n <= pari_PRIMES[0])
     828             :       {
     829   218319485 :         T->p = pari_PRIMES[T->n];
     830   218319485 :         if (T->p > T->b) return 0;
     831   218132825 :         if (T->q == 1 || T->p % T->q == T->c) return T->p;
     832             :       }
     833             :       else
     834             :       { /* beyond the table */
     835         154 :         T->strategy = T->isieve? PRST_sieve: PRST_unextprime;
     836         154 :         if (T->q != 1) { arith_set(T); if (!T->p) return 0; }
     837             :         /* T->p possibly not a prime if q != 1 */
     838         154 :         break;
     839             :       }
     840             :     }
     841             :   }
     842   142732515 :   if (T->strategy == PRST_sieve)
     843             :   {
     844             :     ulong n;
     845   129453922 :     if (T->cache[0]) return shift_cache(T);
     846    92520242 : NEXT_CHUNK:
     847    92527356 :     if (T->psieve)
     848             :     {
     849     5072842 :       T->sieve = T->psieve->sieve;
     850     5072842 :       T->end = T->psieve->end;
     851     5072842 :       if (T->end > T->sieveb) T->end = T->sieveb;
     852     5072842 :       T->maxpos = T->psieve->maxpos;
     853     5072842 :       T->pos = 0;
     854     5072842 :       T->psieve = NULL;
     855             :     }
     856   140235121 :     for (n = T->pos; n < T->maxpos; n++)
     857   140225180 :       if (T->sieve[n] != 0xFF)
     858             :       {
     859    92517415 :         unsigned char mask = T->sieve[n];
     860    92517415 :         ulong p = T->a + (n<<4);
     861    92517415 :         long i = 0;
     862    92517415 :         T->pos = n;
     863    92517415 :         if (!(mask &  1)) T->cache[i++] = p;
     864    92517415 :         if (!(mask &  2)) T->cache[i++] = p+2;
     865    92517415 :         if (!(mask &  4)) T->cache[i++] = p+4;
     866    92517415 :         if (!(mask &  8)) T->cache[i++] = p+6;
     867    92517415 :         if (!(mask & 16)) T->cache[i++] = p+8;
     868    92517415 :         if (!(mask & 32)) T->cache[i++] = p+10;
     869    92517415 :         if (!(mask & 64)) T->cache[i++] = p+12;
     870    92517415 :         if (!(mask &128)) T->cache[i++] = p+14;
     871    92517415 :         T->cache[i] = 0;
     872    92517415 :         T->pos = n+1;
     873    92517415 :         return shift_cache(T);
     874             :       }
     875             :     /* n = T->maxpos, last cell: check p <= b */
     876        9941 :     if (T->maxpos && n == T->maxpos && T->sieve[n] != 0xFF)
     877             :     {
     878        2752 :       unsigned char mask = T->sieve[n];
     879        2752 :       ulong p = T->a + (n<<4);
     880        2752 :       long i = 0;
     881        2752 :       T->pos = n;
     882        2752 :       if (!(mask &  1) && p <= T->sieveb) T->cache[i++] = p;
     883        2752 :       if (!(mask &  2) && p <= T->sieveb-2) T->cache[i++] = p+2;
     884        2752 :       if (!(mask &  4) && p <= T->sieveb-4) T->cache[i++] = p+4;
     885        2752 :       if (!(mask &  8) && p <= T->sieveb-6) T->cache[i++] = p+6;
     886        2752 :       if (!(mask & 16) && p <= T->sieveb-8) T->cache[i++] = p+8;
     887        2752 :       if (!(mask & 32) && p <= T->sieveb-10) T->cache[i++] = p+10;
     888        2752 :       if (!(mask & 64) && p <= T->sieveb-12) T->cache[i++] = p+12;
     889        2752 :       if (!(mask &128) && p <= T->sieveb-14) T->cache[i++] = p+14;
     890        2752 :       if (i)
     891             :       {
     892        2598 :         T->cache[i] = 0;
     893        2598 :         T->pos = n+1;
     894        2598 :         return shift_cache(T);
     895             :       }
     896             :     }
     897             : 
     898        7343 :     if (T->maxpos && T->end >= T->sieveb) /* done with sieves ? */
     899             :     {
     900         230 :       if (T->sieveb == T->b && T->b != ULONG_MAX) return 0;
     901           1 :       T->strategy = PRST_unextprime;
     902             :     }
     903             :     else
     904             :     { /* initialize next chunk */
     905        7113 :       T->sieve = T->isieve;
     906        7113 :       if (T->maxpos == 0)
     907        3177 :         T->a |= 1; /* first time; ensure odd */
     908             :       else
     909        3936 :         T->a = (T->end + 2) | 1;
     910        7113 :       T->end = T->a + T->chunk; /* may overflow */
     911        7113 :       if (T->end < T->a || T->end > T->sieveb) T->end = T->sieveb;
     912             :       /* end and a are odd; sieve[k] contains the a + 8*2k + (0,2,...,14).
     913             :        * The largest k is (end-a) >> 4 */
     914        7113 :       T->pos = 0;
     915        7113 :       T->maxpos = (T->end - T->a) >> 4;
     916        7113 :       sieve_block(T->a, T->end, T->maxpos, T->sieve);
     917        7114 :       goto NEXT_CHUNK;
     918             :     }
     919             :   }
     920    13278594 :   if (T->strategy == PRST_unextprime)
     921             :   {
     922    13274710 :     if (T->q == 1)
     923             :     {
     924             : #ifdef LONG_IS_64BIT
     925    13121589 :       switch(T->p)
     926             :       {
     927             : #define retp(x) return T->p = (HIGHBIT+x <= T->b)? HIGHBIT+x: 0
     928     7449562 :         case HIGHBIT: retp(29);
     929     3158893 :         case HIGHBIT + 29: retp(99);
     930      334745 :         case HIGHBIT + 99: retp(123);
     931      186465 :         case HIGHBIT +123: retp(131);
     932      128184 :         case HIGHBIT +131: retp(155);
     933      107415 :         case HIGHBIT +155: retp(255);
     934       86173 :         case HIGHBIT +255: retp(269);
     935       76702 :         case HIGHBIT +269: retp(359);
     936       62394 :         case HIGHBIT +359: retp(435);
     937       54974 :         case HIGHBIT +435: retp(449);
     938       48276 :         case HIGHBIT +449: retp(453);
     939       45410 :         case HIGHBIT +453: retp(485);
     940       39417 :         case HIGHBIT +485: retp(491);
     941       36324 :         case HIGHBIT +491: retp(543);
     942       33943 :         case HIGHBIT +543: retp(585);
     943       31368 :         case HIGHBIT +585: retp(599);
     944       27401 :         case HIGHBIT +599: retp(753);
     945       26627 :         case HIGHBIT +753: retp(849);
     946       25673 :         case HIGHBIT +849: retp(879);
     947       24095 :         case HIGHBIT +879: retp(885);
     948       23405 :         case HIGHBIT +885: retp(903);
     949       22925 :         case HIGHBIT +903: retp(995);
     950             : #undef retp
     951             :       }
     952             : #endif
     953     1091286 :       T->p = unextprime(T->p + 1);
     954     1093222 :       if (T->p > T->b) return 0;
     955             :     }
     956             :     else do {
     957     2798464 :       T->p += T->q;
     958     2798464 :       if (T->p < T->q || T->p > T->b) { T->p = 0; break; } /* overflow */
     959     2798438 :     } while (!uisprime(T->p));
     960     1246043 :     if (T->p && T->p <= T->b) return T->p;
     961             :     /* overflow ulong, switch to GEN */
     962          66 :     T->strategy = PRST_nextprime;
     963             :   }
     964        3950 :   return 0; /* overflow */
     965             : }
     966             : 
     967             : GEN
     968    45113960 : forprime_next(forprime_t *T)
     969             : {
     970             :   pari_sp av;
     971             :   GEN p;
     972    45113960 :   if (T->strategy != PRST_nextprime)
     973             :   {
     974    45106081 :     ulong u = u_forprime_next(T);
     975    45106081 :     if (u) { affui(u, T->pp); return T->pp; }
     976             :     /* failure */
     977         590 :     if (T->strategy != PRST_nextprime) return NULL; /* we're done */
     978             :     /* overflow ulong, switch to GEN */
     979          48 :     u = ULONG_MAX;
     980          48 :     if (T->q > 1) u -= (ULONG_MAX-T->c) % T->q;
     981          48 :     affui(u, T->pp);
     982             :   }
     983        7927 :   av = avma; p = T->pp;
     984        7927 :   if (T->q == 1)
     985             :   {
     986        7749 :     p = nextprime(addiu(p, 1));
     987        7749 :     if (T->bb && abscmpii(p, T->bb) > 0) return gc_NULL(av);
     988             :   } else do {
     989        3341 :     p = T->qq? addii(p, T->qq): addiu(p, T->q);
     990        3341 :     if (T->bb && abscmpii(p, T->bb) > 0) return gc_NULL(av);
     991        3285 :   } while (!BPSW_psp(p));
     992        7744 :   affii(p, T->pp); return gc_const(av, T->pp);
     993             : }
     994             : 
     995             : void
     996         861 : forprimestep(GEN a, GEN b, GEN q, GEN code)
     997             : {
     998         861 :   pari_sp av = avma;
     999             :   forprime_t T;
    1000             : 
    1001         861 :   if (!forprimestep_init(&T, a,b,q)) { set_avma(av); return; }
    1002             : 
    1003         847 :   push_lex(T.pp,code);
    1004       37823 :   while(forprime_next(&T))
    1005             :   {
    1006       37403 :     closure_evalvoid(code); if (loop_break()) break;
    1007             :     /* p changed in 'code', complain */
    1008       36983 :     if (get_lex(-1) != T.pp)
    1009           7 :       pari_err(e_MISC,"prime index read-only: was changed to %Ps", get_lex(-1));
    1010             :   }
    1011         840 :   pop_lex(1); set_avma(av);
    1012             : }
    1013             : void
    1014         735 : forprime(GEN a, GEN b, GEN code) { return forprimestep(a,b,NULL,code); }
    1015             : 
    1016             : int
    1017          70 : forcomposite_init(forcomposite_t *C, GEN a, GEN b)
    1018             : {
    1019          70 :   pari_sp av = avma;
    1020          70 :   a = gceil(a);
    1021          70 :   if (typ(a)!=t_INT) pari_err_TYPE("forcomposite",a);
    1022          70 :   if (b) {
    1023          63 :     if (typ(b) == t_INFINITY) b = NULL;
    1024             :     else
    1025             :     {
    1026          56 :       b = gfloor(b);
    1027          56 :       if (typ(b)!=t_INT) pari_err_TYPE("forcomposite",b);
    1028             :     }
    1029             :   }
    1030          70 :   if (signe(a) < 0) pari_err_DOMAIN("forcomposite", "a", "<", gen_0, a);
    1031          70 :   if (abscmpiu(a, 4) < 0) a = utoipos(4);
    1032          70 :   C->first = 1;
    1033          70 :   if (!forprime_init(&C->T, a,b) && cmpii(a,b) > 0)
    1034             :   {
    1035           7 :     C->n = gen_1; /* in case caller forgets to check the return value */
    1036           7 :     C->b = gen_0; return gc_bool(av,0);
    1037             :   }
    1038          63 :   C->n = setloop(a);
    1039          63 :   C->b = b;
    1040          63 :   C->p = NULL; return 1;
    1041             : }
    1042             : 
    1043             : GEN
    1044         238 : forcomposite_next(forcomposite_t *C)
    1045             : {
    1046         238 :   if (C->first) /* first call ever */
    1047             :   {
    1048          63 :     C->first = 0;
    1049          63 :     C->p = forprime_next(&C->T);
    1050             :   }
    1051             :   else
    1052         175 :     C->n = incloop(C->n);
    1053         238 :   if (C->p)
    1054             :   {
    1055         161 :     if (cmpii(C->n, C->p) < 0) return C->n;
    1056          77 :     C->n = incloop(C->n);
    1057             :     /* n = p+1 */
    1058          77 :     C->p = forprime_next(&C->T); /* nextprime(p) > n */
    1059          77 :     if (C->p) return C->n;
    1060             :   }
    1061         105 :   if (!C->b || cmpii(C->n, C->b) <= 0) return C->n;
    1062          42 :   return NULL;
    1063             : }
    1064             : 
    1065             : void
    1066          70 : forcomposite(GEN a, GEN b, GEN code)
    1067             : {
    1068          70 :   pari_sp av = avma;
    1069             :   forcomposite_t T;
    1070             :   GEN n;
    1071          70 :   if (!forcomposite_init(&T,a,b)) return;
    1072          63 :   push_lex(T.n,code);
    1073         238 :   while((n = forcomposite_next(&T)))
    1074             :   {
    1075         196 :     closure_evalvoid(code); if (loop_break()) break;
    1076             :     /* n changed in 'code', complain */
    1077         182 :     if (get_lex(-1) != n)
    1078           7 :       pari_err(e_MISC,"index read-only: was changed to %Ps", get_lex(-1));
    1079             :   }
    1080          56 :   pop_lex(1); set_avma(av);
    1081             : }

Generated by: LCOV version 1.16