Code coverage tests

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

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

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

LCOV - code coverage report
Current view: top level - language - init.c (source / functions) Hit Total Coverage
Test: PARI/GP v2.10.0 lcov report (development 21703-2e13792) Lines: 832 1185 70.2 %
Date: 2018-01-17 06:18:31 Functions: 101 134 75.4 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /* Copyright (C) 2000-2003  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             : /*******************************************************************/
      15             : /*                                                                 */
      16             : /*        INITIALIZING THE SYSTEM, ERRORS, STACK MANAGEMENT        */
      17             : /*                                                                 */
      18             : /*******************************************************************/
      19             : /* _GNU_SOURCE is needed before first include to get RUSAGE_THREAD */
      20             : #undef _GNU_SOURCE /* avoid warning */
      21             : #define _GNU_SOURCE
      22             : #include <string.h>
      23             : #if defined(_WIN32) || defined(__CYGWIN32__)
      24             : #  include "../systems/mingw/mingw.h"
      25             : #  include <process.h>
      26             : #endif
      27             : #include "paricfg.h"
      28             : #if defined(STACK_CHECK) && !defined(__EMX__)
      29             : #  include <sys/types.h>
      30             : #  include <sys/time.h>
      31             : #  include <sys/resource.h>
      32             : #endif
      33             : #if defined(HAS_WAITPID) && defined(HAS_SETSID)
      34             : #  include <sys/wait.h>
      35             : #endif
      36             : #ifdef HAS_MMAP
      37             : #  include <sys/mman.h>
      38             : #endif
      39             : #if defined(USE_GETTIMEOFDAY) || defined(USE_GETRUSAGE) || defined(USE_TIMES)
      40             : #  include <sys/time.h>
      41             : #endif
      42             : #if defined(USE_GETRUSAGE)
      43             : #  include <sys/resource.h>
      44             : #endif
      45             : #if defined(USE_FTIME) || defined(USE_FTIMEFORWALLTIME)
      46             : #  include <sys/timeb.h>
      47             : #endif
      48             : #if defined(USE_CLOCK_GETTIME) || defined(USE_TIMES)
      49             : #  include <time.h>
      50             : #endif
      51             : #if defined(USE_TIMES)
      52             : #  include <sys/times.h>
      53             : #endif
      54             : #include "pari.h"
      55             : #include "paripriv.h"
      56             : #include "anal.h"
      57             : 
      58             : const double LOG2    = 0.6931471805599453; /* log(2) */
      59             : const double LOG10_2 = 0.3010299956639812; /* log_10(2) */
      60             : const double LOG2_10 = 3.321928094887362;  /* log_2(10) */
      61             : 
      62             : GEN gnil, gen_0, gen_1, gen_m1, gen_2, gen_m2, ghalf, err_e_STACK;
      63             : 
      64             : static const ulong readonly_constants[] = {
      65             :   evaltyp(t_INT) | _evallg(2),  /* gen_0 */
      66             :   evallgefint(2),
      67             :   evaltyp(t_INT) | _evallg(2),  /* gnil */
      68             :   evallgefint(2),
      69             :   evaltyp(t_INT) | _evallg(3),  /* gen_1 */
      70             :   evalsigne(1) | evallgefint(3),
      71             :   1,
      72             :   evaltyp(t_INT) | _evallg(3),  /* gen_2 */
      73             :   evalsigne(1) | evallgefint(3),
      74             :   2,
      75             :   evaltyp(t_INT) | _evallg(3),  /* gen_m1 */
      76             :   evalsigne(-1) | evallgefint(3),
      77             :   1,
      78             :   evaltyp(t_INT) | _evallg(3),  /* gen_m2 */
      79             :   evalsigne(-1) | evallgefint(3),
      80             :   2,
      81             :   evaltyp(t_FRAC) | _evallg(3), /* ghalf */
      82             :   (ulong)(readonly_constants+4),
      83             :   (ulong)(readonly_constants+7)
      84             : };
      85             : static const ulong readonly_err_STACK[] = {
      86             :   evaltyp(t_ERROR) | _evallg(2),
      87             :   e_STACK
      88             : };
      89             : THREAD GEN bernzone, primetab;
      90             : byteptr diffptr;
      91             : FILE    *pari_outfile, *pari_errfile, *pari_logfile, *pari_infile;
      92             : char    *current_logfile, *current_psfile, *pari_datadir;
      93             : long    gp_colors[c_LAST];
      94             : int     disable_color;
      95             : ulong   DEBUGFILES, DEBUGLEVEL, DEBUGMEM;
      96             : long    DEBUGVAR;
      97             : ulong   pari_mt_nbthreads;
      98             : long    precreal;
      99             : ulong   precdl, logstyle;
     100             : gp_data *GP_DATA;
     101             : 
     102             : GEN colormap, pari_graphcolors;
     103             : 
     104             : entree  **varentries;
     105             : THREAD long *varpriority;
     106             : 
     107             : THREAD pari_sp avma;
     108             : THREAD struct pari_mainstack *pari_mainstack;
     109             : 
     110             : static void ** MODULES;
     111             : static pari_stack s_MODULES;
     112             : const long functions_tblsz = 135; /* size of functions_hash */
     113             : entree **functions_hash, **defaults_hash;
     114             : 
     115             : char *(*cb_pari_fgets_interactive)(char *s, int n, FILE *f);
     116             : int (*cb_pari_get_line_interactive)(const char*, const char*, filtre_t *F);
     117             : void (*cb_pari_quit)(long);
     118             : void (*cb_pari_init_histfile)(void);
     119             : void (*cb_pari_ask_confirm)(const char *);
     120             : int  (*cb_pari_handle_exception)(long);
     121             : int  (*cb_pari_err_handle)(GEN);
     122             : int  (*cb_pari_whatnow)(PariOUT *out, const char *, int);
     123             : void (*cb_pari_sigint)(void);
     124             : void (*cb_pari_pre_recover)(long);
     125             : void (*cb_pari_err_recover)(long);
     126             : int (*cb_pari_break_loop)(int);
     127             : int (*cb_pari_is_interactive)(void);
     128             : void (*cb_pari_start_output)();
     129             : 
     130             : const char * pari_library_path = NULL;
     131             : 
     132             : static THREAD GEN global_err_data;
     133             : THREAD jmp_buf *iferr_env;
     134             : const long CATCH_ALL = -1;
     135             : 
     136             : static void pari_init_timer(void);
     137             : 
     138             : /*********************************************************************/
     139             : /*                                                                   */
     140             : /*                       BLOCKS & CLONES                             */
     141             : /*                                                                   */
     142             : /*********************************************************************/
     143             : /*#define DEBUG*/
     144             : static THREAD long next_block;
     145             : static THREAD GEN cur_block; /* current block in block list */
     146             : #ifdef DEBUG
     147             : static THREAD long NUM;
     148             : #endif
     149             : 
     150             : static void
     151       87699 : pari_init_blocks(void)
     152             : {
     153       87699 :   next_block = 0; cur_block = NULL;
     154             : #ifdef DEBUG
     155             :   NUM = 0;
     156             : #endif
     157       87699 : }
     158             : 
     159             : static void
     160       88481 : pari_close_blocks(void)
     161             : {
     162       88481 :   while (cur_block) killblock(cur_block);
     163       88783 : }
     164             : 
     165             : /* Return x, where:
     166             :  * x[-4]: reference count
     167             :  * x[-3]: adress of next block
     168             :  * x[-2]: adress of preceding block.
     169             :  * x[-1]: number of allocated blocs.
     170             :  * x[0..n-1]: malloc-ed memory. */
     171             : GEN
     172   208976063 : newblock(size_t n)
     173             : {
     174   208976063 :   long *x = (long *) pari_malloc((n + BL_HEAD)*sizeof(long)) + BL_HEAD;
     175             : 
     176   208978252 :   bl_refc(x) = 1;
     177   208978252 :   bl_next(x) = NULL;
     178   208978252 :   bl_prev(x) = cur_block;
     179   208978252 :   bl_num(x)  = next_block++;
     180   208978252 :   if (cur_block) bl_next(cur_block) = x;
     181             : #ifdef DEBUG
     182             :   err_printf("+ %ld\n", ++NUM);
     183             : #endif
     184   208978252 :   if (DEBUGMEM > 2)
     185           0 :     err_printf("new block, size %6lu (no %ld): %08lx\n", n, next_block-1, x);
     186   208977367 :   return cur_block = x;
     187             : }
     188             : 
     189             : GEN
     190       50402 : gcloneref(GEN x)
     191             : {
     192       50402 :   if (isclone(x)) { ++bl_refc(x); return x; }
     193       50402 :   else return gclone(x);
     194             : }
     195             : 
     196             : void
     197           0 : gclone_refc(GEN x) { ++bl_refc(x); }
     198             : 
     199             : void
     200   251854392 : gunclone(GEN x)
     201             : {
     202   503709322 :   if (--bl_refc(x) > 0) return;
     203   208978272 :   BLOCK_SIGINT_START;
     204   208979116 :   if (bl_next(x)) bl_prev(bl_next(x)) = bl_prev(x);
     205             :   else
     206             :   {
     207    25760777 :     cur_block = bl_prev(x);
     208    25760777 :     next_block = bl_num(x);
     209             :   }
     210   208979116 :   if (bl_prev(x)) bl_next(bl_prev(x)) = bl_next(x);
     211   208979116 :   if (DEBUGMEM > 2)
     212           0 :     err_printf("killing block (no %ld): %08lx\n", bl_num(x), x);
     213   208979278 :   free((void*)bl_base(x)); /* pari_free not needed: we already block */
     214   208979278 :   BLOCK_SIGINT_END;
     215             : #ifdef DEBUG
     216             :   err_printf("- %ld\n", NUM--);
     217             : #endif
     218             : }
     219             : 
     220             : /* Recursively look for clones in the container and kill them. Then kill
     221             :  * container if clone. SIGINT could be blocked until it returns */
     222             : void
     223  2813270216 : gunclone_deep(GEN x)
     224             : {
     225             :   long i, lx;
     226             :   GEN v;
     227  5626540432 :   if (isclone(x) && bl_refc(x) > 1) { --bl_refc(x); return; }
     228  2768929059 :   BLOCK_SIGINT_START;
     229  2768929059 :   switch(typ(x))
     230             :   {
     231             :     case t_VEC: case t_COL: case t_MAT:
     232   106647704 :       lx = lg(x);
     233   106647704 :       for (i=1;i<lx;i++) gunclone_deep(gel(x,i));
     234   106647704 :       break;
     235             :     case t_LIST:
     236         270 :       v = list_data(x); lx = v? lg(v): 1;
     237         270 :       for (i=1;i<lx;i++) gunclone_deep(gel(v,i));
     238         270 :       if (v) killblock(v);
     239         270 :       break;
     240             :   }
     241  2768929059 :   if (isclone(x)) gunclone(x);
     242  2768929059 :   BLOCK_SIGINT_END;
     243             : }
     244             : 
     245             : int
     246      150345 : pop_entree_block(entree *ep, long loc)
     247             : {
     248      150345 :   GEN x = (GEN)ep->value;
     249      150345 :   if (bl_num(x) < loc) return 0; /* older */
     250         168 :   if (DEBUGMEM>2)
     251           0 :     err_printf("popping %s (block no %ld)\n", ep->name, bl_num(x));
     252         168 :   gunclone_deep(x); return 1;
     253             : }
     254             : 
     255             : /*********************************************************************/
     256             : /*                                                                   */
     257             : /*                       C STACK SIZE CONTROL                        */
     258             : /*                                                                   */
     259             : /*********************************************************************/
     260             : /* Avoid core dump on deep recursion. Adapted Perl code by Dominic Dunlop */
     261             : THREAD void *PARI_stack_limit = NULL;
     262             : 
     263             : #ifdef STACK_CHECK
     264             : 
     265             : #  ifdef __EMX__                                /* Emulate */
     266             : void
     267             : pari_stackcheck_init(void *pari_stack_base)
     268             : {
     269             :   (void) pari_stack_base;
     270             :   if (!pari_stack_base) { PARI_stack_limit = NULL; return; }
     271             :   PARI_stack_limit = get_stack(1./16, 32*1024);
     272             : }
     273             : #  else /* !__EMX__ */
     274             : /* Set PARI_stack_limit to (a little above) the lowest safe address that can be
     275             :  * used on the stack. Leave PARI_stack_limit at its initial value (NULL) to
     276             :  * show no check should be made [init failed]. Assume stack grows downward. */
     277             : void
     278        1519 : pari_stackcheck_init(void *pari_stack_base)
     279             : {
     280             :   struct rlimit rip;
     281             :   ulong size;
     282        1519 :   if (!pari_stack_base) { PARI_stack_limit = NULL; return; }
     283        1519 :   if (getrlimit(RLIMIT_STACK, &rip)) return;
     284        1519 :   size = rip.rlim_cur;
     285        1519 :   if (size == (ulong)RLIM_INFINITY || size > (ulong)pari_stack_base)
     286           0 :     PARI_stack_limit = (void*)(((ulong)pari_stack_base) / 16);
     287             :   else
     288        1519 :     PARI_stack_limit = (void*)((ulong)pari_stack_base - (size/16)*15);
     289             : }
     290             : #  endif /* !__EMX__ */
     291             : 
     292             : #else
     293             : void
     294             : pari_stackcheck_init(void *pari_stack_base)
     295             : {
     296             :   (void) pari_stack_base; PARI_stack_limit = NULL;
     297             : }
     298             : #endif /* STACK_CHECK */
     299             : 
     300             : /*******************************************************************/
     301             : /*                         HEAP TRAVERSAL                          */
     302             : /*******************************************************************/
     303             : struct getheap_t { long n, l; };
     304             : static void
     305        6699 : f_getheap(GEN x, void *D)
     306             : {
     307        6699 :   struct getheap_t *T = (struct getheap_t*)D;
     308        6699 :   T->n++;
     309        6699 :   T->l += gsizeword(x);
     310        6699 : }
     311             : GEN
     312          84 : getheap(void)
     313             : {
     314          84 :   struct getheap_t T = { 0, 0 };
     315          84 :   traverseheap(&f_getheap, &T);
     316          84 :   return mkvec2s(T.n, T.l + BL_HEAD * T.n);
     317             : }
     318             : 
     319             : void
     320          84 : traverseheap( void(*f)(GEN, void *), void *data )
     321             : {
     322             :   GEN x;
     323          84 :   for (x = cur_block; x; x = bl_prev(x)) f(x, data);
     324          84 : }
     325             : 
     326             : /*********************************************************************/
     327             : /*                          DAEMON / FORK                            */
     328             : /*********************************************************************/
     329             : #if defined(HAS_WAITPID) && defined(HAS_SETSID)
     330             : /* Properly fork a process, detaching from main process group without creating
     331             :  * zombies on exit. Parent returns 1, son returns 0 */
     332             : int
     333          60 : pari_daemon(void)
     334             : {
     335          60 :   pid_t pid = fork();
     336         120 :   switch(pid) {
     337           0 :       case -1: return 1; /* father, fork failed */
     338             :       case 0:
     339          60 :         (void)setsid(); /* son becomes process group leader */
     340          60 :         if (fork()) _exit(0); /* now son exits, also when fork fails */
     341           0 :         break; /* grandson: its father is the son, which exited,
     342             :                 * hence father becomes 'init', that'll take care of it */
     343             :       default: /* father, fork succeeded */
     344          60 :         (void)waitpid(pid,NULL,0); /* wait for son to exit, immediate */
     345          60 :         return 1;
     346             :   }
     347             :   /* grandson */
     348           0 :   return 0;
     349             : }
     350             : #else
     351             : int
     352             : pari_daemon(void)
     353             : {
     354             :   pari_err_IMPL("pari_daemon without waitpid & setsid");
     355             :   return 0;
     356             : }
     357             : #endif
     358             : 
     359             : /*********************************************************************/
     360             : /*                                                                   */
     361             : /*                       SYSTEM INITIALIZATION                       */
     362             : /*                                                                   */
     363             : /*********************************************************************/
     364             : static int try_to_recover = 0;
     365             : THREAD VOLATILE int PARI_SIGINT_block = 0, PARI_SIGINT_pending = 0;
     366             : 
     367             : /*********************************************************************/
     368             : /*                         SIGNAL HANDLERS                           */
     369             : /*********************************************************************/
     370             : static void
     371           0 : dflt_sigint_fun(void) { pari_err(e_MISC, "user interrupt"); }
     372             : 
     373             : #if defined(_WIN32) || defined(__CYGWIN32__)
     374             : int win32ctrlc = 0, win32alrm = 0;
     375             : void
     376             : dowin32ctrlc(void)
     377             : {
     378             :   win32ctrlc = 0;
     379             :   cb_pari_sigint();
     380             : }
     381             : #endif
     382             : 
     383             : static void
     384           0 : pari_handle_SIGINT(void)
     385             : {
     386             : #ifdef _WIN32
     387             :   if (++win32ctrlc >= 5) _exit(3);
     388             : #else
     389           0 :   cb_pari_sigint();
     390             : #endif
     391           0 : }
     392             : 
     393             : typedef void (*pari_sighandler_t)(int);
     394             : 
     395             : pari_sighandler_t
     396       16697 : os_signal(int sig, pari_sighandler_t f)
     397             : {
     398             : #ifdef HAS_SIGACTION
     399             :   struct sigaction sa, oldsa;
     400             : 
     401       16697 :   sa.sa_handler = f;
     402       16697 :   sigemptyset(&sa.sa_mask);
     403       16697 :   sa.sa_flags = SA_NODEFER;
     404             : 
     405       16697 :   if (sigaction(sig, &sa, &oldsa)) return NULL;
     406       16697 :   return oldsa.sa_handler;
     407             : #else
     408             :   return signal(sig,f);
     409             : #endif
     410             : }
     411             : 
     412             : void
     413           0 : pari_sighandler(int sig)
     414             : {
     415             :   const char *msg;
     416             : #ifndef HAS_SIGACTION
     417             :   /*SYSV reset the signal handler in the handler*/
     418             :   (void)os_signal(sig,pari_sighandler);
     419             : #endif
     420           0 :   switch(sig)
     421             :   {
     422             : #ifdef SIGBREAK
     423             :     case SIGBREAK:
     424             :       if (PARI_SIGINT_block==1)
     425             :       {
     426             :         PARI_SIGINT_pending=SIGBREAK;
     427             :         mt_sigint();
     428             :       }
     429             :       else pari_handle_SIGINT();
     430             :       return;
     431             : #endif
     432             : 
     433             : #ifdef SIGINT
     434             :     case SIGINT:
     435           0 :       if (PARI_SIGINT_block==1)
     436             :       {
     437           0 :         PARI_SIGINT_pending=SIGINT;
     438           0 :         mt_sigint();
     439             :       }
     440           0 :       else pari_handle_SIGINT();
     441           0 :       return;
     442             : #endif
     443             : 
     444             : #ifdef SIGSEGV
     445             :     case SIGSEGV:
     446           0 :       msg="PARI/GP (Segmentation Fault)"; break;
     447             : #endif
     448             : #ifdef SIGBUS
     449             :     case SIGBUS:
     450           0 :       msg="PARI/GP (Bus Error)"; break;
     451             : #endif
     452             : #ifdef SIGFPE
     453             :     case SIGFPE:
     454           0 :       msg="PARI/GP (Floating Point Exception)"; break;
     455             : #endif
     456             : 
     457             : #ifdef SIGPIPE
     458             :     case SIGPIPE:
     459             :     {
     460           0 :       pariFILE *f = GP_DATA->pp->file;
     461           0 :       if (f && pari_outfile == f->file)
     462             :       {
     463           0 :         GP_DATA->pp->file = NULL; /* to avoid oo recursion on error */
     464           0 :         pari_outfile = stdout; pari_fclose(f);
     465           0 :         pari_err(e_MISC, "Broken Pipe, resetting file stack...");
     466             :       }
     467             :       return; /* LCOV_EXCL_LINE */
     468             :     }
     469             : #endif
     470             : 
     471           0 :     default: msg="signal handling"; break;
     472             :   }
     473           0 :   pari_err_BUG(msg);
     474             : }
     475             : 
     476             : void
     477        3036 : pari_sig_init(void (*f)(int))
     478             : {
     479             : #ifdef SIGBUS
     480        3036 :   (void)os_signal(SIGBUS,f);
     481             : #endif
     482             : #ifdef SIGFPE
     483        3036 :   (void)os_signal(SIGFPE,f);
     484             : #endif
     485             : #ifdef SIGINT
     486        3036 :   (void)os_signal(SIGINT,f);
     487             : #endif
     488             : #ifdef SIGBREAK
     489             :   (void)os_signal(SIGBREAK,f);
     490             : #endif
     491             : #ifdef SIGPIPE
     492        3036 :   (void)os_signal(SIGPIPE,f);
     493             : #endif
     494             : #ifdef SIGSEGV
     495        3036 :   (void)os_signal(SIGSEGV,f);
     496             : #endif
     497        3036 : }
     498             : 
     499             : /*********************************************************************/
     500             : /*                      STACK AND UNIVERSAL CONSTANTS                */
     501             : /*********************************************************************/
     502             : static void
     503        1519 : init_universal_constants(void)
     504             : {
     505        1519 :   gen_0  = (GEN)readonly_constants;
     506        1519 :   gnil   = (GEN)readonly_constants+2;
     507        1519 :   gen_1  = (GEN)readonly_constants+4;
     508        1519 :   gen_2  = (GEN)readonly_constants+7;
     509        1519 :   gen_m1 = (GEN)readonly_constants+10;
     510        1519 :   gen_m2 = (GEN)readonly_constants+13;
     511        1519 :   ghalf  = (GEN)readonly_constants+16;
     512        1519 :   err_e_STACK = (GEN)readonly_err_STACK;
     513        1519 : }
     514             : 
     515             : static void
     516       87811 : pari_init_errcatch(void)
     517             : {
     518       87811 :   iferr_env = NULL;
     519       87811 :   global_err_data = NULL;
     520       87811 : }
     521             : 
     522             : /*********************************************************************/
     523             : /*                           INIT DEFAULTS                           */
     524             : /*********************************************************************/
     525             : void
     526        1519 : pari_init_defaults(void)
     527             : {
     528             :   long i;
     529        1519 :   initout(1);
     530             : 
     531             : #ifdef LONG_IS_64BIT
     532        1304 :   precreal = 128;
     533             : #else
     534         215 :   precreal = 96;
     535             : #endif
     536             : 
     537        1519 :   precdl = 16;
     538        1519 :   DEBUGFILES = DEBUGLEVEL = 0;
     539        1519 :   DEBUGMEM = 1;
     540        1519 :   disable_color = 1;
     541        1519 :   logstyle = logstyle_none;
     542             : 
     543        1519 :   current_psfile = pari_strdup("pari.ps");
     544        1519 :   current_logfile= pari_strdup("pari.log");
     545        1519 :   pari_logfile = NULL;
     546             : 
     547        1519 :   pari_datadir = os_getenv("GP_DATA_DIR");
     548        1519 :   if (!pari_datadir)
     549             :   {
     550             : #if defined(_WIN32) || defined(__CYGWIN32__)
     551             :     if (paricfg_datadir[0]=='@' && paricfg_datadir[1]==0)
     552             :       pari_datadir = win32_datadir();
     553             :     else
     554             : #endif
     555        1519 :       pari_datadir = pari_strdup(paricfg_datadir);
     556             :   }
     557           0 :   else pari_datadir= pari_strdup(pari_datadir);
     558        1519 :   for (i=0; i<c_LAST; i++) gp_colors[i] = c_NONE;
     559        1519 :   colormap = NULL; pari_graphcolors = NULL;
     560        1519 : }
     561             : 
     562             : /*********************************************************************/
     563             : /*                   FUNCTION HASHTABLES, MODULES                    */
     564             : /*********************************************************************/
     565             : extern entree functions_basic[], functions_default[];
     566             : static void
     567        1519 : pari_init_functions(void)
     568             : {
     569        1519 :   pari_stack_init(&s_MODULES, sizeof(*MODULES),(void**)&MODULES);
     570        1519 :   pari_stack_pushp(&s_MODULES,functions_basic);
     571        1519 :   functions_hash = (entree**) pari_calloc(sizeof(entree*)*functions_tblsz);
     572        1519 :   pari_fill_hashtable(functions_hash, functions_basic);
     573        1519 :   defaults_hash = (entree**) pari_calloc(sizeof(entree*)*functions_tblsz);
     574        1519 :   pari_add_defaults_module(functions_default);
     575        1519 : }
     576             : 
     577             : void
     578        1517 : pari_add_module(entree *ep)
     579             : {
     580        1517 :   pari_fill_hashtable(functions_hash, ep);
     581        1517 :   pari_stack_pushp(&s_MODULES, ep);
     582        1517 : }
     583             : 
     584             : void
     585        1519 : pari_add_defaults_module(entree *ep)
     586        1519 : { pari_fill_hashtable(defaults_hash, ep); }
     587             : 
     588             : /*********************************************************************/
     589             : /*                       PARI MAIN STACK                             */
     590             : /*********************************************************************/
     591             : 
     592             : #ifdef HAS_MMAP
     593             : #define PARI_STACK_ALIGN (sysconf(_SC_PAGE_SIZE))
     594             : #ifndef MAP_ANONYMOUS
     595             : #define MAP_ANONYMOUS MAP_ANON
     596             : #endif
     597             : #ifndef MAP_NORESERVE
     598             : #define MAP_NORESERVE 0
     599             : #endif
     600             : static void *
     601       89264 : pari_mainstack_malloc(size_t size)
     602             : {
     603             :   /* Check that the system allows reserving "size" bytes. This is just
     604             :    * a check, we immediately free the memory. */
     605       89264 :   void *b = mmap(NULL, size, PROT_READ|PROT_WRITE,
     606             :                              MAP_PRIVATE|MAP_ANONYMOUS, -1, 0);
     607       89264 :   if (b == MAP_FAILED) return NULL;
     608       89264 :   munmap(b, size);
     609             : 
     610             :   /* Map again, this time with MAP_NORESERVE. On some operating systems
     611             :    * like Cygwin, this is needed because remapping with PROT_NONE and
     612             :    * MAP_NORESERVE does not work as expected. */
     613       89264 :   b = mmap(NULL, size, PROT_READ|PROT_WRITE,
     614             :                        MAP_PRIVATE|MAP_ANONYMOUS|MAP_NORESERVE, -1, 0);
     615       89264 :   if (b == MAP_FAILED) return NULL;
     616       89264 :   return b;
     617             : }
     618             : 
     619             : static void
     620       89262 : pari_mainstack_mfree(void *s, size_t size)
     621             : {
     622       89262 :   munmap(s, size);
     623       89262 : }
     624             : 
     625             : /* Completely discard the memory mapped between the addresses "from"
     626             :  * and "to" (which must be page-aligned).
     627             :  *
     628             :  * We use mmap() with PROT_NONE, which means that the underlying memory
     629             :  * is freed and that the kernel should not commit memory for it. We
     630             :  * still keep the mapping such that we can change the flags to
     631             :  * PROT_READ|PROT_WRITE later.
     632             :  *
     633             :  * NOTE: remapping with MAP_FIXED and PROT_NONE is not the same as
     634             :  * calling mprotect(..., PROT_NONE) because the latter will keep the
     635             :  * memory committed (this is in particular relevant on Linux with
     636             :  * vm.overcommit = 2). This remains true even when calling
     637             :  * madvise(..., MADV_DONTNEED). */
     638             : static void
     639      161013 : pari_mainstack_mreset(pari_sp from, pari_sp to)
     640             : {
     641      161013 :   size_t s = to - from;
     642             :   void *addr, *res;
     643      322026 :   if (!s) return;
     644             : 
     645           0 :   addr = (void*)from;
     646           0 :   res = mmap(addr, s, PROT_NONE,
     647             :              MAP_FIXED|MAP_PRIVATE|MAP_ANONYMOUS|MAP_NORESERVE, -1, 0);
     648           0 :   if (res != addr) pari_err(e_MEM);
     649             : }
     650             : 
     651             : /* Commit (make available) the virtual memory mapped between the
     652             :  * addresses "from" and "to" (which must be page-aligned).
     653             :  * Return 0 if successful, -1 if failed. */
     654             : static int
     655      161013 : pari_mainstack_mextend(pari_sp from, pari_sp to)
     656             : {
     657      161013 :   size_t s = to - from;
     658      161013 :   return mprotect((void*)from, s, PROT_READ|PROT_WRITE);
     659             : }
     660             : 
     661             : /* Set actual stack size to the given size. This sets st->size and
     662             :  * st->bot. If not enough system memory is available, this can fail.
     663             :  * Return 1 if successful, 0 if failed (in that case, st->size is not
     664             :  * changed) */
     665             : static int
     666      161013 : pari_mainstack_setsize(struct pari_mainstack *st, size_t size)
     667             : {
     668      161013 :   pari_sp newbot = st->top - size;
     669             :   /* Align newbot to pagesize */
     670      161013 :   pari_sp alignbot = newbot & ~(pari_sp)(PARI_STACK_ALIGN - 1);
     671      161013 :   if (pari_mainstack_mextend(alignbot, st->top))
     672             :   {
     673             :     /* Making the memory available did not work: limit vsize to the
     674             :      * current actual stack size. */
     675           0 :     st->vsize = st->size;
     676           0 :     pari_warn(warnstack, st->vsize);
     677           0 :     return 0;
     678             :   }
     679      161013 :   pari_mainstack_mreset(st->vbot, alignbot);
     680      161013 :   st->bot = newbot;
     681      161013 :   st->size = size;
     682      161013 :   return 1;
     683             : }
     684             : 
     685             : #else
     686             : #define PARI_STACK_ALIGN (0x40UL)
     687             : static void *
     688             : pari_mainstack_malloc(size_t s)
     689             : {
     690             :   return malloc(s); /* NOT pari_malloc, e_MEM would be deadly */
     691             : }
     692             : 
     693             : static void
     694             : pari_mainstack_mfree(void *s, size_t size) { (void) size; free(s); }
     695             : 
     696             : static int
     697             : pari_mainstack_setsize(struct pari_mainstack *st, size_t size)
     698             : {
     699             :   st->bot = st->top - size;
     700             :   st->size = size;
     701             :   return 1;
     702             : }
     703             : 
     704             : #endif
     705             : 
     706             : static const size_t MIN_STACK = 500032UL;
     707             : static size_t
     708      178526 : fix_size(size_t a)
     709             : {
     710      178526 :   size_t ps = PARI_STACK_ALIGN;
     711      178526 :   size_t b = a & ~(ps - 1); /* Align */
     712      178526 :   if (b < a && b < ~(ps - 1)) b += ps;
     713      178526 :   if (b < MIN_STACK) b = MIN_STACK;
     714      178526 :   return b;
     715             : }
     716             : 
     717             : static void
     718       89264 : pari_mainstack_alloc(struct pari_mainstack *st, size_t rsize, size_t vsize)
     719             : {
     720       89264 :   size_t sizemax = vsize ? vsize: rsize, s = fix_size(sizemax);
     721             :   for (;;)
     722             :   {
     723       89264 :     st->vbot = (pari_sp)pari_mainstack_malloc(s);
     724       89264 :     if (st->vbot) break;
     725           0 :     if (s == MIN_STACK) pari_err(e_MEM); /* no way out. Die */
     726           0 :     s = fix_size(s >> 1);
     727           0 :     pari_warn(warnstack, s);
     728           0 :   }
     729       89264 :   st->vsize = vsize ? s: 0;
     730       89264 :   st->rsize = minuu(rsize, s);
     731       89264 :   st->top = st->vbot+s;
     732       89264 :   if (!pari_mainstack_setsize(st, st->rsize))
     733             :   {
     734             :     /* This should never happen since we only decrease the allocated space */
     735           0 :     pari_err(e_MEM);
     736             :   }
     737       89264 :   st->memused = 0;
     738       89264 : }
     739             : 
     740             : static void
     741       89262 : pari_mainstack_free(struct pari_mainstack *st)
     742             : {
     743       89262 :   pari_mainstack_mfree((void*)st->vbot, st->vsize ? st->vsize : fix_size(st->rsize));
     744       89262 :   st->top = st->bot = st->vbot = 0;
     745       89262 :   st->size = st->vsize = 0;
     746       89262 : }
     747             : 
     748             : static void
     749         264 : pari_mainstack_resize(struct pari_mainstack *st, size_t rsize, size_t vsize)
     750             : {
     751         264 :   BLOCK_SIGINT_START;
     752         264 :   pari_mainstack_free(st);
     753         264 :   pari_mainstack_alloc(st, rsize, vsize);
     754         264 :   BLOCK_SIGINT_END;
     755         264 : }
     756             : 
     757             : static void
     758       87712 : pari_mainstack_use(struct pari_mainstack *st)
     759             : {
     760       87712 :   pari_mainstack = st;
     761       87712 :   avma = st->top;
     762       87712 : }
     763             : 
     764             : static void
     765        1519 : paristack_alloc(size_t rsize, size_t vsize)
     766             : {
     767        1519 :   pari_mainstack_alloc(pari_mainstack, rsize, vsize);
     768        1519 :   pari_mainstack_use(pari_mainstack);
     769        1519 : }
     770             : 
     771             : void
     772           0 : paristack_setsize(size_t rsize, size_t vsize)
     773             : {
     774           0 :   pari_mainstack_resize(pari_mainstack, rsize, vsize);
     775           0 :   pari_mainstack_use(pari_mainstack);
     776           0 : }
     777             : 
     778             : void
     779           0 : parivstack_resize(ulong newsize)
     780             : {
     781             :   size_t s;
     782           0 :   if (newsize && newsize < pari_mainstack->rsize)
     783           0 :     pari_err_DIM("stack sizes [parisizemax < parisize]");
     784           0 :   if (newsize == pari_mainstack->vsize) return;
     785           0 :   evalstate_reset();
     786           0 :   paristack_setsize(pari_mainstack->rsize, newsize);
     787           0 :   s = pari_mainstack->vsize ? pari_mainstack->vsize : pari_mainstack->rsize;
     788           0 :   if (DEBUGMEM)
     789           0 :     pari_warn(warner,"new maximum stack size = %lu (%.3f Mbytes)",
     790             :               s, s/1048576.);
     791           0 :   pari_init_errcatch();
     792           0 :   cb_pari_err_recover(-1);
     793             : }
     794             : 
     795             : void
     796         270 : paristack_newrsize(ulong newsize)
     797             : {
     798         270 :   size_t s, vsize = pari_mainstack->vsize;
     799         270 :   if (!newsize) newsize = pari_mainstack->rsize << 1;
     800         270 :   if (newsize != pari_mainstack->rsize)
     801         264 :     pari_mainstack_resize(pari_mainstack, newsize, vsize);
     802         270 :   evalstate_reset();
     803         270 :   s = pari_mainstack->rsize;
     804         270 :   if (DEBUGMEM)
     805         270 :     pari_warn(warner,"new stack size = %lu (%.3f Mbytes)", s, s/1048576.);
     806         270 :   pari_init_errcatch();
     807         270 :   cb_pari_err_recover(-1);
     808           0 : }
     809             : 
     810             : void
     811           0 : paristack_resize(ulong newsize)
     812             : {
     813           0 :   long size = pari_mainstack->size;
     814           0 :   if (!newsize)
     815           0 :     newsize = 2 * size;
     816           0 :   newsize = minuu(newsize, pari_mainstack->vsize);
     817           0 :   if (newsize <= pari_mainstack->size) return;
     818           0 :   if (pari_mainstack_setsize(pari_mainstack, newsize))
     819             :   {
     820           0 :     if (DEBUGMEM)
     821           0 :       pari_warn(warner, "increasing stack size to %lu", pari_mainstack->size);
     822             :   }
     823             :   else
     824             :   {
     825           0 :     pari_mainstack_setsize(pari_mainstack, size);
     826           0 :     pari_err(e_STACK);
     827             :   }
     828             : }
     829             : 
     830             : void
     831       71749 : parivstack_reset(void)
     832             : {
     833       71749 :   pari_mainstack_setsize(pari_mainstack, pari_mainstack->rsize);
     834       71749 :   if (avma < pari_mainstack->bot)
     835           0 :     pari_err_BUG("parivstack_reset [avma < bot]");
     836       71749 : }
     837             : 
     838             : /* Enlarge the stack if needed such that the unused portion of the stack
     839             :  * (between bot and avma) is large enough to contain x longs. */
     840             : void
     841           7 : new_chunk_resize(size_t x)
     842             : {
     843           7 :   if (pari_mainstack->vsize==0
     844           0 :     || x > (avma-pari_mainstack->vbot) / sizeof(long)) pari_err(e_STACK);
     845           0 :   while (x > (avma-pari_mainstack->bot) / sizeof(long))
     846           0 :     paristack_resize(0);
     847           0 : }
     848             : 
     849             : /*********************************************************************/
     850             : /*                       PARI THREAD                                 */
     851             : /*********************************************************************/
     852             : 
     853             : /* Initial PARI thread structure t with a stack of size s and virtual size v
     854             :  * and argument arg */
     855             : 
     856             : void
     857           0 : pari_thread_valloc(struct pari_thread *t, size_t s, size_t v, GEN arg)
     858             : {
     859           0 :   pari_mainstack_alloc(&t->st,s,v);
     860           0 :   t->data = arg;
     861           0 : }
     862             : 
     863             : /* Initial PARI thread structure t with a stack of size s and
     864             :  * argument arg */
     865             : 
     866             : void
     867       87481 : pari_thread_alloc(struct pari_thread *t, size_t s, GEN arg)
     868             : {
     869       87481 :   pari_mainstack_alloc(&t->st,s,0);
     870       87481 :   t->data = arg;
     871       87481 : }
     872             : 
     873             : void
     874       87481 : pari_thread_free(struct pari_thread *t)
     875             : {
     876       87481 :   pari_mainstack_free(&t->st);
     877       87481 : }
     878             : 
     879             : void
     880       88483 : pari_thread_init(void)
     881             : {
     882       88483 :   pari_init_blocks();
     883       87906 :   pari_init_errcatch();
     884       87799 :   pari_init_rand();
     885       87976 :   pari_init_floats();
     886       87711 :   pari_init_parser();
     887       88303 :   pari_init_compiler();
     888       88474 :   pari_init_evaluator();
     889       88065 :   pari_init_files();
     890       88133 :   pari_thread_init_primetab();
     891       87836 :   pari_thread_init_seadata();
     892       87629 : }
     893             : 
     894             : void
     895       10278 : pari_thread_sync(void)
     896             : {
     897       10278 :   pari_pthread_init_primetab();
     898       10278 :   pari_pthread_init_seadata();
     899       10278 :   pari_pthread_init_varstate();
     900       10278 : }
     901             : 
     902             : void
     903       88791 : pari_thread_close(void)
     904             : {
     905       88791 :   pari_thread_close_files();
     906       88817 :   pari_close_evaluator();
     907       88728 :   pari_close_compiler();
     908       88664 :   pari_close_parser();
     909       88596 :   pari_close_floats();
     910       88548 :   pari_close_blocks();
     911       88186 : }
     912             : 
     913             : GEN
     914       86170 : pari_thread_start(struct pari_thread *t)
     915             : {
     916       86170 :   pari_mainstack_use(&t->st);
     917       87084 :   pari_thread_init();
     918       86115 :   pari_thread_init_varstate();
     919       87412 :   return t->data;
     920             : }
     921             : 
     922             : /*********************************************************************/
     923             : /*                       LIBPARI INIT / CLOSE                        */
     924             : /*********************************************************************/
     925             : 
     926             : static void
     927           0 : pari_exit(void)
     928             : {
     929           0 :   err_printf("  ***   Error in the PARI system. End of program.\n");
     930           0 :   exit(1);
     931             : }
     932             : 
     933             : static void
     934           0 : dflt_err_recover(long errnum) { (void) errnum; pari_exit(); }
     935             : 
     936             : static void
     937           0 : dflt_pari_quit(long err) { (void)err; /*do nothing*/; }
     938             : 
     939             : static int pari_err_display(GEN err);
     940             : 
     941             : /* initialize PARI data. Initialize [new|old]fun to NULL for default set. */
     942             : void
     943        1519 : pari_init_opts(size_t parisize, ulong maxprime, ulong init_opts)
     944             : {
     945             :   ulong u;
     946             : 
     947        1519 :   pari_mt_nbthreads = 0;
     948        1519 :   cb_pari_quit = dflt_pari_quit;
     949        1519 :   cb_pari_init_histfile = NULL;
     950        1519 :   cb_pari_get_line_interactive = NULL;
     951        1519 :   cb_pari_fgets_interactive = NULL;
     952        1519 :   cb_pari_whatnow = NULL;
     953        1519 :   cb_pari_handle_exception = NULL;
     954        1519 :   cb_pari_err_handle = pari_err_display;
     955        1519 :   cb_pari_pre_recover = NULL;
     956        1519 :   cb_pari_break_loop = NULL;
     957        1519 :   cb_pari_is_interactive = NULL;
     958        1519 :   cb_pari_start_output = NULL;
     959        1519 :   cb_pari_sigint = dflt_sigint_fun;
     960        1519 :   if (init_opts&INIT_JMPm) cb_pari_err_recover = dflt_err_recover;
     961             : 
     962        1519 :   pari_stackcheck_init(&u);
     963        1519 :   pari_init_homedir();
     964        1519 :   if (init_opts&INIT_DFTm) {
     965           0 :     pari_init_defaults();
     966           0 :     GP_DATA = default_gp_data();
     967           0 :     pari_init_paths();
     968             :   }
     969             : 
     970        1519 :   pari_mainstack = (struct pari_mainstack *) malloc(sizeof(*pari_mainstack));
     971        1519 :   paristack_alloc(parisize, 0);
     972        1519 :   init_universal_constants();
     973        1519 :   diffptr = NULL;
     974        1519 :   if (!(init_opts&INIT_noPRIMEm))  pari_init_primes(maxprime);
     975        1519 :   if (!(init_opts&INIT_noINTGMPm)) pari_kernel_init();
     976        1519 :   pari_init_graphics();
     977        1519 :   pari_init_primetab();
     978        1519 :   pari_init_seadata();
     979        1519 :   pari_thread_init();
     980        1519 :   pari_init_functions();
     981        1519 :   pari_var_init();
     982        1519 :   pari_init_timer();
     983        1519 :   pari_init_buffers();
     984        1519 :   (void)getabstime();
     985        1519 :   try_to_recover = 1;
     986        1519 :   if (!(init_opts&INIT_noIMTm)) pari_mt_init();
     987        1519 :   if ((init_opts&INIT_SIGm)) pari_sig_init(pari_sighandler);
     988        1519 : }
     989             : 
     990             : void
     991           0 : pari_init(size_t parisize, ulong maxprime)
     992           0 : { pari_init_opts(parisize, maxprime, INIT_JMPm | INIT_SIGm | INIT_DFTm); }
     993             : 
     994             : void
     995        1517 : pari_close_opts(ulong init_opts)
     996             : {
     997             :   long i;
     998             : 
     999        1517 :   BLOCK_SIGINT_START;
    1000        1517 :   if ((init_opts&INIT_SIGm)) pari_sig_init(SIG_DFL);
    1001        1517 :   if (!(init_opts&INIT_noIMTm)) pari_mt_close();
    1002             : 
    1003      206312 :   for (i = 0; i < functions_tblsz; i++)
    1004             :   {
    1005      204795 :     entree *ep = functions_hash[i];
    1006     2161659 :     while (ep) {
    1007     1752069 :       entree *EP = ep->next;
    1008     1752069 :       if (!EpSTATIC(ep)) { freeep(ep); free(ep); }
    1009     1752069 :       ep = EP;
    1010             :     }
    1011             :   }
    1012        1517 :   pari_var_close();
    1013        1517 :   pari_close_mf();
    1014        1517 :   pari_thread_close();
    1015        1517 :   pari_close_files();
    1016        1517 :   pari_close_homedir();
    1017        1517 :   if (!(init_opts&INIT_noINTGMPm)) pari_kernel_close();
    1018             : 
    1019        1517 :   free((void*)functions_hash);
    1020        1517 :   free((void*)defaults_hash);
    1021        1517 :   if (diffptr) pari_close_primes();
    1022        1517 :   free(current_logfile);
    1023        1517 :   free(current_psfile);
    1024        1517 :   pari_mainstack_free(pari_mainstack);
    1025        1517 :   free((void*)pari_mainstack);
    1026        1517 :   pari_stack_delete(&s_MODULES);
    1027        1517 :   if (pari_datadir) free(pari_datadir);
    1028        1517 :   if (init_opts&INIT_DFTm)
    1029             :   { /* delete GP_DATA */
    1030        1517 :     pari_close_paths();
    1031        1517 :     if (GP_DATA->hist->v) free((void*)GP_DATA->hist->v);
    1032        1517 :     if (GP_DATA->pp->cmd) free((void*)GP_DATA->pp->cmd);
    1033        1517 :     if (GP_DATA->help) free((void*)GP_DATA->help);
    1034        1517 :     free((void*)GP_DATA->prompt);
    1035        1517 :     free((void*)GP_DATA->prompt_cont);
    1036        1517 :     free((void*)GP_DATA->histfile);
    1037             :   }
    1038        1517 :   BLOCK_SIGINT_END;
    1039        1517 : }
    1040             : 
    1041             : void
    1042        1517 : pari_close(void)
    1043        1517 : { pari_close_opts(INIT_JMPm | INIT_SIGm | INIT_DFTm); }
    1044             : 
    1045             : /*******************************************************************/
    1046             : /*                                                                 */
    1047             : /*                         ERROR RECOVERY                          */
    1048             : /*                                                                 */
    1049             : /*******************************************************************/
    1050             : void
    1051       86884 : gp_context_save(struct gp_context* rec)
    1052             : {
    1053       86884 :   rec->file = pari_last_tmp_file();
    1054       86884 :   if (DEBUGFILES>1)
    1055           0 :     err_printf("gp_context_save: %s\n", rec->file ? rec->file->name: "NULL");
    1056       86884 :   rec->prettyp = GP_DATA->fmt->prettyp;
    1057       86884 :   rec->listloc = next_block;
    1058       86884 :   rec->iferr_env = iferr_env;
    1059       86884 :   rec->err_data  = global_err_data;
    1060       86884 :   varstate_save(&rec->var);
    1061       86884 :   evalstate_save(&rec->eval);
    1062       86884 :   parsestate_save(&rec->parse);
    1063       86884 : }
    1064             : 
    1065             : void
    1066        8308 : gp_context_restore(struct gp_context* rec)
    1067             : {
    1068             :   long i;
    1069             : 
    1070       16616 :   if (!try_to_recover) return;
    1071             :   /* disable gp_context_restore() and SIGINT */
    1072        8308 :   try_to_recover = 0;
    1073        8308 :   BLOCK_SIGINT_START
    1074        8308 :   if (DEBUGMEM>2) err_printf("entering recover(), loc = %ld\n", rec->listloc);
    1075        8308 :   evalstate_restore(&rec->eval);
    1076        8308 :   parsestate_restore(&rec->parse);
    1077        8308 :   filestate_restore(rec->file);
    1078        8308 :   global_err_data = rec->err_data;
    1079        8308 :   iferr_env = rec->iferr_env;
    1080        8308 :   GP_DATA->fmt->prettyp = rec->prettyp;
    1081             : 
    1082     1129888 :   for (i = 0; i < functions_tblsz; i++)
    1083             :   {
    1084     1121580 :     entree *ep = functions_hash[i];
    1085    12481373 :     while (ep)
    1086             :     {
    1087    10238213 :       entree *EP = ep->next;
    1088    10238213 :       switch(EpVALENCE(ep))
    1089             :       {
    1090             :         case EpVAR:
    1091      172757 :           while (pop_val_if_newer(ep,rec->listloc)) /* empty */;
    1092      172757 :           break;
    1093      610952 :         case EpNEW: break;
    1094             :       }
    1095    10238213 :       ep = EP;
    1096             :     }
    1097             :   }
    1098        8308 :   varstate_restore(&rec->var);
    1099        8308 :   if (DEBUGMEM>2) err_printf("leaving recover()\n");
    1100        8308 :   BLOCK_SIGINT_END
    1101        8308 :   try_to_recover = 1;
    1102             : }
    1103             : 
    1104             : static void
    1105        8238 : err_recover(long numerr)
    1106             : {
    1107        8238 :   if (cb_pari_pre_recover)
    1108        8238 :     cb_pari_pre_recover(numerr);
    1109           0 :   evalstate_reset();
    1110           0 :   killallfiles();
    1111           0 :   pari_init_errcatch();
    1112           0 :   cb_pari_err_recover(numerr);
    1113           0 : }
    1114             : 
    1115             : static void
    1116        8753 : err_init(void)
    1117             : {
    1118             :   /* make sure pari_err msg starts at the beginning of line */
    1119        8753 :   if (!pari_last_was_newline()) pari_putc('\n');
    1120        8753 :   pariOut->flush();
    1121        8753 :   pariErr->flush();
    1122        8753 :   out_term_color(pariErr, c_ERR);
    1123        8753 : }
    1124             : 
    1125             : static void
    1126        8718 : err_init_msg(int user)
    1127             : {
    1128             :   const char *gp_function_name;
    1129        8718 :   out_puts(pariErr, "  *** ");
    1130        8718 :   if (!user && (gp_function_name = closure_func_err()))
    1131        5268 :     out_printf(pariErr, "%s: ", gp_function_name);
    1132             :   else
    1133        3450 :     out_puts(pariErr, "  ");
    1134        8718 : }
    1135             : 
    1136             : void
    1137         494 : pari_warn(int numerr, ...)
    1138             : {
    1139             :   char *ch1;
    1140             :   va_list ap;
    1141             : 
    1142         494 :   va_start(ap,numerr);
    1143             : 
    1144         494 :   err_init();
    1145         494 :   err_init_msg(numerr==warnuser || numerr==warnstack);
    1146         494 :   switch (numerr)
    1147             :   {
    1148             :     case warnuser:
    1149           7 :       out_puts(pariErr, "user warning: ");
    1150           7 :       out_print0(pariErr, NULL, va_arg(ap, GEN), f_RAW);
    1151           7 :       break;
    1152             : 
    1153             :     case warnmem:
    1154           0 :       out_puts(pariErr, "collecting garbage in "); ch1=va_arg(ap, char*);
    1155           0 :       out_vprintf(pariErr, ch1,ap); out_putc(pariErr, '.');
    1156           0 :       break;
    1157             : 
    1158             :     case warner:
    1159         487 :       out_puts(pariErr, "Warning: "); ch1=va_arg(ap, char*);
    1160         487 :       out_vprintf(pariErr, ch1,ap); out_putc(pariErr, '.');
    1161         487 :       break;
    1162             : 
    1163             :     case warnprec:
    1164           0 :       out_vprintf(pariErr, "Warning: increasing prec in %s; new prec = %ld",
    1165             :                       ap);
    1166           0 :       break;
    1167             : 
    1168             :     case warnfile:
    1169           0 :       out_puts(pariErr, "Warning: failed to "),
    1170           0 :       ch1 = va_arg(ap, char*);
    1171           0 :       out_printf(pariErr, "%s: %s", ch1, va_arg(ap, char*));
    1172           0 :       break;
    1173             : 
    1174             :     case warnstack:
    1175             :     {
    1176           0 :       ulong  s = va_arg(ap, ulong);
    1177             :       char buf[128];
    1178           0 :       sprintf(buf,"Warning: not enough memory, new stack %lu", (ulong)s);
    1179           0 :       out_puts(pariErr,buf);
    1180           0 :       break;
    1181             :     }
    1182             : 
    1183             :   }
    1184         494 :   va_end(ap);
    1185         494 :   out_term_color(pariErr, c_NONE);
    1186         494 :   out_putc(pariErr, '\n');
    1187         494 :   pariErr->flush();
    1188         494 : }
    1189             : void
    1190           0 : pari_sigint(const char *time_s)
    1191             : {
    1192           0 :   int recover=0;
    1193           0 :   BLOCK_SIGALRM_START
    1194           0 :   err_init();
    1195           0 :   closure_err(0);
    1196           0 :   err_init_msg(0);
    1197           0 :   out_puts(pariErr, "user interrupt after ");
    1198           0 :   out_puts(pariErr, time_s);
    1199           0 :   out_term_color(pariErr, c_NONE);
    1200           0 :   pariErr->flush();
    1201           0 :   if (cb_pari_handle_exception)
    1202           0 :     recover = cb_pari_handle_exception(-1);
    1203           0 :   if (!recover && !block)
    1204           0 :     PARI_SIGINT_pending = 0;
    1205           0 :   BLOCK_SIGINT_END
    1206           0 :   if (!recover) err_recover(e_MISC);
    1207           0 : }
    1208             : 
    1209             : #define retmkerr2(x,y)\
    1210             :   do { GEN _v = cgetg(3, t_ERROR);\
    1211             :        _v[1] = (x);\
    1212             :        gel(_v,2) = (y); return _v; } while(0)
    1213             : #define retmkerr3(x,y,z)\
    1214             :   do { GEN _v = cgetg(4, t_ERROR);\
    1215             :        _v[1] = (x);\
    1216             :        gel(_v,2) = (y);\
    1217             :        gel(_v,3) = (z); return _v; } while(0)
    1218             : #define retmkerr4(x,y,z,t)\
    1219             :   do { GEN _v = cgetg(5, t_ERROR);\
    1220             :        _v[1] = (x);\
    1221             :        gel(_v,2) = (y);\
    1222             :        gel(_v,3) = (z);\
    1223             :        gel(_v,4) = (t); return _v; } while(0)
    1224             : #define retmkerr5(x,y,z,t,u)\
    1225             :   do { GEN _v = cgetg(6, t_ERROR);\
    1226             :        _v[1] = (x);\
    1227             :        gel(_v,2) = (y);\
    1228             :        gel(_v,3) = (z);\
    1229             :        gel(_v,4) = (t);\
    1230             :        gel(_v,5) = (u); return _v; } while(0)
    1231             : #define retmkerr6(x,y,z,t,u,v)\
    1232             :   do { GEN _v = cgetg(7, t_ERROR);\
    1233             :        _v[1] = (x);\
    1234             :        gel(_v,2) = (y);\
    1235             :        gel(_v,3) = (z);\
    1236             :        gel(_v,4) = (t);\
    1237             :        gel(_v,5) = (u);\
    1238             :        gel(_v,6) = (v); return _v; } while(0)
    1239             : 
    1240             : static GEN
    1241       42917 : pari_err2GEN(long numerr, va_list ap)
    1242             : {
    1243       42917 :   switch ((enum err_list) numerr)
    1244             :   {
    1245             :   case e_SYNTAX:
    1246             :     {
    1247          35 :       const char *msg = va_arg(ap, char*);
    1248          35 :       const char *s = va_arg(ap,char *);
    1249          35 :       const char *entry = va_arg(ap,char *);
    1250          35 :       retmkerr3(numerr,strtoGENstr(msg), mkvecsmall2((long)s,(long)entry));
    1251             :     }
    1252             :   case e_MISC: case e_ALARM:
    1253             :     {
    1254         114 :       const char *ch1 = va_arg(ap, char*);
    1255         114 :       retmkerr2(numerr, gvsprintf(ch1,ap));
    1256             :     }
    1257             :   case e_NOTFUNC:
    1258             :   case e_USER:
    1259        2737 :     retmkerr2(numerr,va_arg(ap, GEN));
    1260             :   case e_FILE:
    1261             :     {
    1262           0 :       const char *f = va_arg(ap, const char*);
    1263           0 :       retmkerr3(numerr, strtoGENstr(f), strtoGENstr(va_arg(ap, char*)));
    1264             :     }
    1265             :   case e_OVERFLOW:
    1266             :   case e_IMPL:
    1267             :   case e_DIM:
    1268             :   case e_CONSTPOL:
    1269             :   case e_ROOTS0:
    1270             :   case e_FLAG:
    1271             :   case e_PREC:
    1272             :   case e_BUG:
    1273             :   case e_ARCH:
    1274             :   case e_PACKAGE:
    1275        1238 :     retmkerr2(numerr, strtoGENstr(va_arg(ap, char*)));
    1276             :   case e_MODULUS:
    1277             :   case e_VAR:
    1278             :     {
    1279         938 :       const char *f = va_arg(ap, const char*);
    1280         938 :       GEN x = va_arg(ap, GEN);
    1281         938 :       GEN y = va_arg(ap, GEN);
    1282         938 :       retmkerr4(numerr, strtoGENstr(f), x,y);
    1283             :     }
    1284             :   case e_INV:
    1285             :   case e_IRREDPOL:
    1286             :   case e_PRIME:
    1287             :   case e_SQRTN:
    1288             :   case e_TYPE:
    1289             :     {
    1290       32086 :       const char *f = va_arg(ap, const char*);
    1291       32086 :       GEN x = va_arg(ap, GEN);
    1292       32086 :       retmkerr3(numerr, strtoGENstr(f), x);
    1293             :     }
    1294             :   case e_COPRIME: case e_OP: case e_TYPE2:
    1295             :     {
    1296        3423 :       const char *f = va_arg(ap, const char*);
    1297        3423 :       GEN x = va_arg(ap, GEN);
    1298        3423 :       GEN y = va_arg(ap, GEN);
    1299        3423 :       retmkerr4(numerr,strtoGENstr(f),x,y);
    1300             :     }
    1301             :   case e_COMPONENT:
    1302             :     {
    1303         245 :       const char *f= va_arg(ap, const char *);
    1304         245 :       const char *op = va_arg(ap, const char *);
    1305         245 :       GEN l = va_arg(ap, GEN);
    1306         245 :       GEN x = va_arg(ap, GEN);
    1307         245 :       retmkerr5(numerr,strtoGENstr(f),strtoGENstr(op),l,x);
    1308             :     }
    1309             :   case e_DOMAIN:
    1310             :     {
    1311        1947 :       const char *f = va_arg(ap, const char*);
    1312        1947 :       const char *v = va_arg(ap, const char *);
    1313        1947 :       const char *op = va_arg(ap, const char *);
    1314        1947 :       GEN l = va_arg(ap, GEN);
    1315        1947 :       GEN x = va_arg(ap, GEN);
    1316        1947 :       retmkerr6(numerr,strtoGENstr(f),strtoGENstr(v),strtoGENstr(op),l,x);
    1317             :     }
    1318             :   case e_PRIORITY:
    1319             :     {
    1320         147 :       const char *f = va_arg(ap, const char*);
    1321         147 :       GEN x = va_arg(ap, GEN);
    1322         147 :       const char *op = va_arg(ap, const char *);
    1323         147 :       long v = va_arg(ap, long);
    1324         147 :       retmkerr5(numerr,strtoGENstr(f),x,strtoGENstr(op),stoi(v));
    1325             :     }
    1326             :   case e_MAXPRIME:
    1327           0 :     retmkerr2(numerr, utoi(va_arg(ap, ulong)));
    1328             :   case e_STACK:
    1329           7 :     return err_e_STACK;
    1330             :   case e_STACKTHREAD:
    1331           0 :     retmkerr3(numerr, utoi(va_arg(ap, ulong)), utoi(va_arg(ap, ulong)));
    1332             :   default:
    1333           0 :     return mkerr(numerr);
    1334             :   }
    1335             : }
    1336             : 
    1337             : static char *
    1338        6258 : type_dim(GEN x)
    1339             : {
    1340        6258 :   char *v = stack_malloc(64);
    1341        6258 :   switch(typ(x))
    1342             :   {
    1343             :     case t_MAT:
    1344             :     {
    1345          63 :       long l = lg(x), r = (l == 1)? 1: lgcols(x);
    1346          63 :       sprintf(v, "t_MAT (%ldx%ld)", r-1,l-1);
    1347          63 :       break;
    1348             :     }
    1349             :     case t_COL:
    1350          91 :       sprintf(v, "t_COL (%ld elts)", lg(x)-1);
    1351          91 :       break;
    1352             :     case t_VEC:
    1353         133 :       sprintf(v, "t_VEC (%ld elts)", lg(x)-1);
    1354         133 :       break;
    1355             :     default:
    1356        5971 :       v = (char*)type_name(typ(x));
    1357             :   }
    1358        6258 :   return v;
    1359             : }
    1360             : 
    1361             : static char *
    1362        2282 : gdisplay(GEN x)
    1363             : {
    1364        2282 :   char *s = GENtostr_raw(x);
    1365        2282 :   if (strlen(s) < 1600) return s;
    1366          21 :   if (! GP_DATA->breakloop) return (char*)"(...)";
    1367           0 :   return stack_sprintf("\n  ***  (...) Huge %s omitted; you can access it via dbg_err()", type_name(typ(x)));
    1368             : }
    1369             : 
    1370             : char *
    1371       15539 : pari_err2str(GEN e)
    1372             : {
    1373       15539 :   long numerr = err_get_num(e);
    1374       15539 :   switch ((enum err_list) numerr)
    1375             :   {
    1376             :   case e_ALARM:
    1377           0 :     return pari_sprintf("alarm interrupt after %Ps.",gel(e,2));
    1378             :   case e_MISC:
    1379         112 :     return pari_sprintf("%Ps.",gel(e,2));
    1380             : 
    1381             :   case e_ARCH:
    1382           0 :     return pari_sprintf("sorry, '%Ps' not available on this system.",gel(e,2));
    1383             :   case e_BUG:
    1384          14 :     return pari_sprintf("bug in %Ps, please report.",gel(e,2));
    1385             :   case e_CONSTPOL:
    1386           7 :     return pari_sprintf("constant polynomial in %Ps.", gel(e,2));
    1387             :   case e_COPRIME:
    1388         189 :     return pari_sprintf("elements not coprime in %Ps:\n    %s\n    %s",
    1389         189 :                         gel(e,2), gdisplay(gel(e,3)), gdisplay(gel(e,4)));
    1390             :   case e_DIM:
    1391         459 :     return pari_sprintf("inconsistent dimensions in %Ps.", gel(e,2));
    1392             :   case e_FILE:
    1393           0 :     return pari_sprintf("error opening %Ps: `%Ps'.", gel(e,2), gel(e,3));
    1394             :   case e_FLAG:
    1395          35 :     return pari_sprintf("invalid flag in %Ps.", gel(e,2));
    1396             :   case e_IMPL:
    1397         336 :     return pari_sprintf("sorry, %Ps is not yet implemented.", gel(e,2));
    1398             :   case e_PACKAGE:
    1399           0 :     return pari_sprintf("package %Ps is required, please install it.", gel(e,2));
    1400             :   case e_INV:
    1401         539 :     return pari_sprintf("impossible inverse in %Ps: %s.", gel(e,2),
    1402         539 :                         gdisplay(gel(e,3)));
    1403             :   case e_IRREDPOL:
    1404          42 :     return pari_sprintf("not an irreducible polynomial in %Ps: %s.",
    1405          42 :                         gel(e,2), gdisplay(gel(e,3)));
    1406             :   case e_MAXPRIME:
    1407             :     {
    1408           0 :       const char * msg = "not enough precomputed primes";
    1409           0 :       ulong c = itou(gel(e,2));
    1410           0 :       if (c) return pari_sprintf("%s, need primelimit ~ %lu.",msg, c);
    1411           0 :       else   return pari_strdup(msg);
    1412             :     }
    1413             :   case e_MEM:
    1414           0 :     return pari_strdup("not enough memory");
    1415             :   case e_MODULUS:
    1416             :     {
    1417         735 :       GEN x = gel(e,3), y = gel(e,4);
    1418        1470 :       return pari_sprintf("inconsistent moduli in %Ps: %s != %s",
    1419         735 :                           gel(e,2), gdisplay(x), gdisplay(y));
    1420             :     }
    1421           0 :   case e_NONE: return NULL;
    1422             :   case e_NOTFUNC:
    1423        2723 :     return pari_strdup("not a function in function call");
    1424             :   case e_OP: case e_TYPE2:
    1425             :     {
    1426        3129 :       pari_sp av = avma;
    1427             :       char *v;
    1428        3129 :       const char *f, *op = GSTR(gel(e,2));
    1429        3129 :       const char *what = numerr == e_OP? "inconsistent": "forbidden";
    1430        3129 :       GEN x = gel(e,3);
    1431        3129 :       GEN y = gel(e,4);
    1432        3129 :       switch(*op)
    1433             :       {
    1434          14 :       case '+': f = "addition"; break;
    1435         105 :       case '*': f = "multiplication"; break;
    1436        2387 :       case '/': case '%': case '\\': f = "division"; break;
    1437           0 :       case '=': op = "-->"; f = "assignment"; break;
    1438         623 :       default:  f = op; op = ","; break;
    1439             :       }
    1440        3129 :       v = pari_sprintf("%s %s %s %s %s.", what,f,type_dim(x),op,type_dim(y));
    1441        3129 :       avma = av; return v;
    1442             :     }
    1443             :   case e_COMPONENT:
    1444             :     {
    1445         245 :       const char *f= GSTR(gel(e,2));
    1446         245 :       const char *op= GSTR(gel(e,3));
    1447         245 :       GEN l = gel(e,4);
    1448         245 :       if (!*f)
    1449         147 :         return pari_sprintf("non-existent component: index %s %Ps",op,l);
    1450          98 :       return pari_sprintf("non-existent component in %s: index %s %Ps",f,op,l);
    1451             :     }
    1452             :   case e_DOMAIN:
    1453             :     {
    1454        1849 :       const char *f = GSTR(gel(e,2));
    1455        1849 :       const char *v = GSTR(gel(e,3));
    1456        1849 :       const char *op= GSTR(gel(e,4));
    1457        1849 :       GEN l = gel(e,5);
    1458        1849 :       if (!*op)
    1459          28 :         return pari_sprintf("domain error in %s: %s out of range",f,v);
    1460        1821 :       return pari_sprintf("domain error in %s: %s %s %Ps",f,v,op,l);
    1461             :     }
    1462             :   case e_PRIORITY:
    1463             :     {
    1464          98 :       const char *f = GSTR(gel(e,2));
    1465          98 :       long vx = gvar(gel(e,3));
    1466          98 :       const char *op= GSTR(gel(e,4));
    1467          98 :       long v = itos(gel(e,5));
    1468          98 :       return pari_sprintf("incorrect priority in %s: variable %Ps %s %Ps",f,
    1469             :              pol_x(vx), op, pol_x(v));
    1470             :     }
    1471             :   case e_OVERFLOW:
    1472          91 :     return pari_sprintf("overflow in %Ps.", gel(e,2));
    1473             :   case e_PREC:
    1474         210 :     return pari_sprintf("precision too low in %Ps.", gel(e,2));
    1475             :   case e_PRIME:
    1476          98 :     return pari_sprintf("not a prime number in %Ps: %s.",
    1477          98 :                         gel(e,2), gdisplay(gel(e,3)));
    1478             :   case e_ROOTS0:
    1479          42 :     return pari_sprintf("zero polynomial in %Ps.", gel(e,2));
    1480             :   case e_SQRTN:
    1481         154 :     return pari_sprintf("not an n-th power residue in %Ps: %s.",
    1482         154 :                         gel(e,2), gdisplay(gel(e,3)));
    1483             :   case e_STACK:
    1484             :   case e_STACKTHREAD:
    1485             :     {
    1486           7 :       const char *stack = numerr == e_STACK? "PARI": "thread";
    1487           7 :       const char *var = numerr == e_STACK? "parisizemax": "threadsizemax";
    1488           8 :       size_t rsize = numerr == e_STACKTHREAD && GP_DATA->threadsize ?
    1489           7 :                                 GP_DATA->threadsize: pari_mainstack->rsize;
    1490           7 :       size_t vsize = numerr == e_STACK? pari_mainstack->vsize:
    1491           0 :                                         GP_DATA->threadsizemax;
    1492           7 :       char *buf = (char *) pari_malloc(512*sizeof(char));
    1493           7 :       if (vsize)
    1494             :       {
    1495           0 :         sprintf(buf, "the %s stack overflows !\n"
    1496             :             "  current stack size: %lu (%.3f Mbytes)\n"
    1497             :             "  [hint] you can increase '%s' using default()\n",
    1498           0 :             stack, (ulong)vsize, (double)vsize/1048576., var);
    1499             :       }
    1500             :       else
    1501             :       {
    1502           7 :         sprintf(buf, "the %s stack overflows !\n"
    1503             :             "  current stack size: %lu (%.3f Mbytes)\n"
    1504             :             "  [hint] set '%s' to a non-zero value in your GPRC\n",
    1505           7 :             stack, (ulong)rsize, (double)rsize/1048576., var);
    1506             :       }
    1507           7 :       return buf;
    1508             :     }
    1509             :   case e_SYNTAX:
    1510           0 :     return pari_strdup(GSTR(gel(e,2)));
    1511             :   case e_TYPE:
    1512        8962 :     return pari_sprintf("incorrect type in %Ps (%s).",
    1513        8962 :                         gel(e,2), type_name(typ(gel(e,3))));
    1514             :   case e_USER:
    1515          14 :     return pari_sprint0("user error: ", gel(e,2), f_RAW);
    1516             :   case e_VAR:
    1517             :     {
    1518         203 :       GEN x = gel(e,3), y = gel(e,4);
    1519         609 :       return pari_sprintf("inconsistent variables in %Ps, %Ps != %Ps.",
    1520         609 :                           gel(e,2), pol_x(varn(x)), pol_x(varn(y)));
    1521             :     }
    1522             :   }
    1523             :   return NULL; /*LCOV_EXCL_LINE*/
    1524             : }
    1525             : 
    1526             : static int
    1527        8259 : pari_err_display(GEN err)
    1528             : {
    1529        8259 :   long numerr=err_get_num(err);
    1530        8259 :   err_init();
    1531        8259 :   if (numerr==e_SYNTAX)
    1532             :   {
    1533          35 :     const char *msg = GSTR(gel(err,2));
    1534          35 :     const char *s     = (const char *) gmael(err,3,1);
    1535          35 :     const char *entry = (const char *) gmael(err,3,2);
    1536          35 :     print_errcontext(pariErr, msg, s, entry);
    1537             :   }
    1538             :   else
    1539             :   {
    1540        8224 :     char *s = pari_err2str(err);
    1541        8224 :     closure_err(0);
    1542        8224 :     err_init_msg(numerr==e_USER);
    1543        8224 :     pariErr->puts(s);
    1544        8224 :     if (numerr==e_NOTFUNC)
    1545             :     {
    1546        2723 :       GEN fun = gel(err,2);
    1547        2723 :       if (gequalX(fun))
    1548             :       {
    1549        2723 :         entree *ep = varentries[varn(fun)];
    1550        2723 :         const char *s = ep->name;
    1551        2723 :         if (cb_pari_whatnow) cb_pari_whatnow(pariErr,s,1);
    1552             :       }
    1553             :     }
    1554        8210 :     pari_free(s);
    1555             :   }
    1556        8245 :   out_term_color(pariErr, c_NONE);
    1557        8245 :   pariErr->flush(); return 0;
    1558             : }
    1559             : 
    1560             : void
    1561       42933 : pari_err(int numerr, ...)
    1562             : {
    1563             :   va_list ap;
    1564             :   GEN E;
    1565             : 
    1566       42933 :   va_start(ap,numerr);
    1567             : 
    1568       42933 :   if (numerr)
    1569       42917 :     E = pari_err2GEN(numerr,ap);
    1570             :   else
    1571             :   {
    1572          16 :     E = va_arg(ap,GEN);
    1573          16 :     numerr = err_get_num(E);
    1574             :   }
    1575       42933 :   global_err_data = E;
    1576       42933 :   if (*iferr_env) longjmp(*iferr_env, numerr);
    1577        8261 :   mt_err_recover(numerr);
    1578        8259 :   va_end(ap);
    1579       16504 :   if (cb_pari_err_handle &&
    1580        8259 :       cb_pari_err_handle(E)) return;
    1581       16483 :   if (cb_pari_handle_exception &&
    1582        8245 :       cb_pari_handle_exception(numerr)) return;
    1583        8238 :   err_recover(numerr);
    1584             : }
    1585             : 
    1586             : GEN
    1587       69332 : pari_err_last(void) { return global_err_data; }
    1588             : 
    1589             : const char *
    1590       26947 : numerr_name(long numerr)
    1591             : {
    1592       26947 :   switch ((enum err_list) numerr)
    1593             :   {
    1594           0 :   case e_ALARM:    return "e_ALARM";
    1595           0 :   case e_ARCH:     return "e_ARCH";
    1596           0 :   case e_BUG:      return "e_BUG";
    1597           0 :   case e_COMPONENT: return "e_COMPONENT";
    1598           0 :   case e_CONSTPOL: return "e_CONSTPOL";
    1599           0 :   case e_COPRIME:  return "e_COPRIME";
    1600           0 :   case e_DIM:      return "e_DIM";
    1601          56 :   case e_DOMAIN:   return "e_DOMAIN";
    1602           0 :   case e_FILE:     return "e_FILE";
    1603           7 :   case e_FLAG:     return "e_FLAG";
    1604          28 :   case e_IMPL:     return "e_IMPL";
    1605       19093 :   case e_INV:      return "e_INV";
    1606           0 :   case e_IRREDPOL: return "e_IRREDPOL";
    1607           0 :   case e_MAXPRIME: return "e_MAXPRIME";
    1608           0 :   case e_MEM:      return "e_MEM";
    1609           0 :   case e_MISC:     return "e_MISC";
    1610           0 :   case e_MODULUS:  return "e_MODULUS";
    1611           0 :   case e_NONE:     return "e_NONE";
    1612           0 :   case e_NOTFUNC:  return "e_NOTFUNC";
    1613           0 :   case e_OP:       return "e_OP";
    1614           0 :   case e_OVERFLOW: return "e_OVERFLOW";
    1615           0 :   case e_PACKAGE:  return "e_PACKAGE";
    1616           0 :   case e_PREC:     return "e_PREC";
    1617           0 :   case e_PRIME:    return "e_PRIME";
    1618          49 :   case e_PRIORITY: return "e_PRIORITY";
    1619           0 :   case e_ROOTS0:   return "e_ROOTS0";
    1620           0 :   case e_SQRTN:    return "e_SQRTN";
    1621           0 :   case e_STACK:    return "e_STACK";
    1622           0 :   case e_SYNTAX:   return "e_SYNTAX";
    1623           0 :   case e_STACKTHREAD:   return "e_STACKTHREAD";
    1624           0 :   case e_TYPE2:    return "e_TYPE2";
    1625        7714 :   case e_TYPE:     return "e_TYPE";
    1626           0 :   case e_USER:     return "e_USER";
    1627           0 :   case e_VAR:      return "e_VAR";
    1628             :   }
    1629           0 :   return "invalid error number";
    1630             : }
    1631             : 
    1632             : long
    1633          21 : name_numerr(const char *s)
    1634             : {
    1635          21 :   if (!strcmp(s,"e_ALARM"))    return e_ALARM;
    1636          21 :   if (!strcmp(s,"e_ARCH"))     return e_ARCH;
    1637          21 :   if (!strcmp(s,"e_BUG"))      return e_BUG;
    1638          21 :   if (!strcmp(s,"e_COMPONENT")) return e_COMPONENT;
    1639          21 :   if (!strcmp(s,"e_CONSTPOL")) return e_CONSTPOL;
    1640          21 :   if (!strcmp(s,"e_COPRIME"))  return e_COPRIME;
    1641          21 :   if (!strcmp(s,"e_DIM"))      return e_DIM;
    1642          21 :   if (!strcmp(s,"e_DOMAIN"))   return e_DOMAIN;
    1643          21 :   if (!strcmp(s,"e_FILE"))     return e_FILE;
    1644          21 :   if (!strcmp(s,"e_FLAG"))     return e_FLAG;
    1645          21 :   if (!strcmp(s,"e_IMPL"))     return e_IMPL;
    1646          21 :   if (!strcmp(s,"e_INV"))      return e_INV;
    1647           0 :   if (!strcmp(s,"e_IRREDPOL")) return e_IRREDPOL;
    1648           0 :   if (!strcmp(s,"e_MAXPRIME")) return e_MAXPRIME;
    1649           0 :   if (!strcmp(s,"e_MEM"))      return e_MEM;
    1650           0 :   if (!strcmp(s,"e_MISC"))     return e_MISC;
    1651           0 :   if (!strcmp(s,"e_MODULUS"))  return e_MODULUS;
    1652           0 :   if (!strcmp(s,"e_NONE"))     return e_NONE;
    1653           0 :   if (!strcmp(s,"e_NOTFUNC"))  return e_NOTFUNC;
    1654           0 :   if (!strcmp(s,"e_OP"))       return e_OP;
    1655           0 :   if (!strcmp(s,"e_OVERFLOW")) return e_OVERFLOW;
    1656           0 :   if (!strcmp(s,"e_PACKAGE"))  return e_PACKAGE;
    1657           0 :   if (!strcmp(s,"e_PREC"))     return e_PREC;
    1658           0 :   if (!strcmp(s,"e_PRIME"))    return e_PRIME;
    1659           0 :   if (!strcmp(s,"e_PRIORITY")) return e_PRIORITY;
    1660           0 :   if (!strcmp(s,"e_ROOTS0"))   return e_ROOTS0;
    1661           0 :   if (!strcmp(s,"e_SQRTN"))    return e_SQRTN;
    1662           0 :   if (!strcmp(s,"e_STACK"))    return e_STACK;
    1663           0 :   if (!strcmp(s,"e_SYNTAX"))   return e_SYNTAX;
    1664           0 :   if (!strcmp(s,"e_TYPE"))     return e_TYPE;
    1665           0 :   if (!strcmp(s,"e_TYPE2"))    return e_TYPE2;
    1666           0 :   if (!strcmp(s,"e_USER"))     return e_USER;
    1667           0 :   if (!strcmp(s,"e_VAR"))      return e_VAR;
    1668           0 :   pari_err(e_MISC,"unknown error name");
    1669             :   return -1; /* LCOV_EXCL_LINE */
    1670             : }
    1671             : 
    1672             : GEN
    1673       26947 : errname(GEN err)
    1674             : {
    1675       26947 :   if (typ(err)!=t_ERROR) pari_err_TYPE("errname",err);
    1676       26947 :   return strtoGENstr(numerr_name(err_get_num(err)));
    1677             : }
    1678             : 
    1679             : /* Try f (trapping error e), recover using r (break_loop, if NULL) */
    1680             : GEN
    1681          21 : trap0(const char *e, GEN r, GEN f)
    1682             : {
    1683          21 :   long numerr = CATCH_ALL;
    1684             :   GEN x;
    1685          21 :   if (!e || !*e) numerr = CATCH_ALL;
    1686          21 :   else numerr = name_numerr(e);
    1687          21 :   if (!f) {
    1688           0 :     pari_warn(warner,"default handlers are no longer supported --> ignored");
    1689           0 :     return gnil;
    1690             :   }
    1691          21 :   x = closure_trapgen(f, numerr);
    1692          14 :   if (x == (GEN)1L) x = r? closure_evalgen(r): gnil;
    1693          14 :   return x;
    1694             : }
    1695             : 
    1696             : /*******************************************************************/
    1697             : /*                                                                */
    1698             : /*                       CLONING & COPY                            */
    1699             : /*                  Replicate an existing GEN                      */
    1700             : /*                                                                 */
    1701             : /*******************************************************************/
    1702             : /* lontyp[tx] = 0 (non recursive type) or number of codewords for type tx */
    1703             : const  long lontyp[] = { 0,0,0,1,1,2,1,2,1,1, 2,2,0,1,1,1,1,1,1,1, 2,0,0,2,2,1 };
    1704             : 
    1705             : static GEN
    1706         690 : list_internal_copy(GEN z, long nmax)
    1707             : {
    1708             :   long i, l;
    1709             :   GEN a;
    1710         690 :   if (!z) return NULL;
    1711         560 :   l = lg(z);
    1712         560 :   a = newblock(nmax+1);
    1713         560 :   for (i = 1; i < l; i++) gel(a,i) = gel(z,i)? gclone(gel(z,i)): gen_0;
    1714         560 :   a[0] = z[0]; return a;
    1715             : }
    1716             : 
    1717             : static void
    1718         690 : listassign(GEN x, GEN y)
    1719             : {
    1720         690 :   long nmax = list_nmax(x);
    1721         690 :   GEN L = list_data(x);
    1722         690 :   if (!nmax && L) nmax = lg(L) + 32; /* not malloc'ed yet */
    1723         690 :   y[1] = evaltyp(list_typ(x))|evallg(nmax);
    1724         690 :   list_data(y) = list_internal_copy(L, nmax);
    1725         690 : }
    1726             : 
    1727             : /* copy list on the PARI stack */
    1728             : GEN
    1729         207 : listcopy(GEN x)
    1730             : {
    1731         207 :   GEN y = mklist(), L = list_data(x);
    1732         207 :   if (L) list_data(y) = gcopy(L);
    1733         207 :   y[1] = evaltyp(list_typ(x));
    1734         207 :   return y;
    1735             : }
    1736             : 
    1737             : GEN
    1738  3283114323 : gcopy(GEN x)
    1739             : {
    1740  3283114323 :   long tx = typ(x), lx, i;
    1741             :   GEN y;
    1742  3283114323 :   switch(tx)
    1743             :   { /* non recursive types */
    1744  2982293008 :     case t_INT: return signe(x)? icopy(x): gen_0;
    1745             :     case t_REAL:
    1746             :     case t_STR:
    1747   168247005 :     case t_VECSMALL: return leafcopy(x);
    1748             :     /* one more special case */
    1749         207 :     case t_LIST: return listcopy(x);
    1750             :   }
    1751   132574103 :   y = cgetg_copy(x, &lx);
    1752   132574104 :   if (lontyp[tx] == 1) i = 1; else { y[1] = x[1]; i = 2; }
    1753   132574104 :   for (; i<lx; i++) gel(y,i) = gcopy(gel(x,i));
    1754   132574104 :   return y;
    1755             : }
    1756             : 
    1757             : /* as gcopy, but truncate to the first lx components if recursive type
    1758             :  * [ leaves use their own lg ]. No checks. */
    1759             : GEN
    1760         588 : gcopy_lg(GEN x, long lx)
    1761             : {
    1762         588 :   long tx = typ(x), i;
    1763             :   GEN y;
    1764         588 :   switch(tx)
    1765             :   { /* non recursive types */
    1766           0 :     case t_INT: return signe(x)? icopy(x): gen_0;
    1767             :     case t_REAL:
    1768             :     case t_STR:
    1769           0 :     case t_VECSMALL: return leafcopy(x);
    1770             :     /* one more special case */
    1771           0 :     case t_LIST: return listcopy(x);
    1772             :   }
    1773         588 :   y = cgetg(lx, tx);
    1774         588 :   if (lontyp[tx] == 1) i = 1; else { y[1] = x[1]; i = 2; }
    1775         588 :   for (; i<lx; i++) gel(y,i) = gcopy(gel(x,i));
    1776         588 :   return y;
    1777             : }
    1778             : 
    1779             : /* cf cgetg_copy: "allocate" (by updating first codeword only) for subsequent
    1780             :  * copy of x, as if avma = *AVMA */
    1781             : INLINE GEN
    1782   322058698 : cgetg_copy_avma(GEN x, long *plx, pari_sp *AVMA) {
    1783             :   GEN z;
    1784   322058698 :   *plx = lg(x);
    1785   322058698 :   z = ((GEN)*AVMA) - *plx;
    1786   322058698 :   z[0] = x[0] & (TYPBITS|LGBITS);
    1787   322058698 :   *AVMA = (pari_sp)z; return z;
    1788             : }
    1789             : INLINE GEN
    1790         133 : cgetlist_avma(pari_sp *AVMA)
    1791             : {
    1792         133 :   GEN y = ((GEN)*AVMA) - 3;
    1793         133 :   y[0] = _evallg(3) | evaltyp(t_LIST);
    1794         133 :   *AVMA = (pari_sp)y; return y;
    1795             : }
    1796             : 
    1797             : /* copy x as if avma = *AVMA, update *AVMA */
    1798             : GEN
    1799  2762284029 : gcopy_avma(GEN x, pari_sp *AVMA)
    1800             : {
    1801  2762284029 :   long i, lx, tx = typ(x);
    1802             :   GEN y;
    1803             : 
    1804  2762284029 :   switch(typ(x))
    1805             :   { /* non recursive types */
    1806             :     case t_INT:
    1807  2635539938 :       if (lgefint(x) == 2) return gen_0;
    1808  2178430635 :       *AVMA = (pari_sp)icopy_avma(x, *AVMA);
    1809  2178430640 :       return (GEN)*AVMA;
    1810             :     case t_REAL: case t_STR: case t_VECSMALL:
    1811    27319266 :       *AVMA = (pari_sp)leafcopy_avma(x, *AVMA);
    1812    27319266 :       return (GEN)*AVMA;
    1813             : 
    1814             :     /* one more special case */
    1815             :     case t_LIST:
    1816         133 :       y = cgetlist_avma(AVMA);
    1817         133 :       listassign(x, y); return y;
    1818             : 
    1819             :   }
    1820    99424692 :   y = cgetg_copy_avma(x, &lx, AVMA);
    1821    99424679 :   if (lontyp[tx] == 1) i = 1; else { y[1] = x[1]; i = 2; }
    1822    99424679 :   for (; i<lx; i++) gel(y,i) = gcopy_avma(gel(x,i), AVMA);
    1823    99424703 :   return y;
    1824             : }
    1825             : 
    1826             : /* [copy_bin/bin_copy:] same as gcopy_avma but use NULL to code an exact 0, and
    1827             :  * make shallow copies of finalized t_LISTs */
    1828             : static GEN
    1829  1059787283 : gcopy_av0(GEN x, pari_sp *AVMA)
    1830             : {
    1831  1059787283 :   long i, lx, tx = typ(x);
    1832             :   GEN y;
    1833             : 
    1834  1059787283 :   switch(tx)
    1835             :   { /* non recursive types */
    1836             :     case t_INT:
    1837   712039111 :       if (!signe(x)) return NULL; /* special marker */
    1838   401182538 :       *AVMA = (pari_sp)icopy_avma(x, *AVMA);
    1839   401182541 :       return (GEN)*AVMA;
    1840             :     case t_LIST:
    1841          49 :       if (list_data(x) && !list_nmax(x)) break; /* not finalized, need copy */
    1842             :       /* else finalized: shallow copy */
    1843             :     case t_REAL: case t_STR: case t_VECSMALL:
    1844   125114398 :       *AVMA = (pari_sp)leafcopy_avma(x, *AVMA);
    1845   125114398 :       return (GEN)*AVMA;
    1846             :   }
    1847   222633774 :   y = cgetg_copy_avma(x, &lx, AVMA);
    1848   222633882 :   if (lontyp[tx] == 1) i = 1; else { y[1] = x[1]; i = 2; }
    1849   222633882 :   for (; i<lx; i++) gel(y,i) = gcopy_av0(gel(x,i), AVMA);
    1850   222634049 :   return y;
    1851             : }
    1852             : 
    1853             : INLINE GEN
    1854           0 : icopy_avma_canon(GEN x, pari_sp AVMA)
    1855             : {
    1856           0 :   long i, lx = lgefint(x);
    1857           0 :   GEN y = ((GEN)AVMA) - lx;
    1858           0 :   y[0] = evaltyp(t_INT)|evallg(lx); /* kills isclone */
    1859           0 :   y[1] = x[1]; x = int_MSW(x);
    1860           0 :   for (i=2; i<lx; i++, x = int_precW(x)) y[i] = *x;
    1861           0 :   return y;
    1862             : }
    1863             : 
    1864             : /* [copy_bin_canon:] same as gcopy_av0, but copy integers in
    1865             :  * canonical (native kernel) form and make a full copy of t_LISTs */
    1866             : static GEN
    1867           0 : gcopy_av0_canon(GEN x, pari_sp *AVMA)
    1868             : {
    1869           0 :   long i, lx, tx = typ(x);
    1870             :   GEN y;
    1871             : 
    1872           0 :   switch(tx)
    1873             :   { /* non recursive types */
    1874             :     case t_INT:
    1875           0 :       if (!signe(x)) return NULL; /* special marker */
    1876           0 :       *AVMA = (pari_sp)icopy_avma_canon(x, *AVMA);
    1877           0 :       return (GEN)*AVMA;
    1878             :     case t_REAL: case t_STR: case t_VECSMALL:
    1879           0 :       *AVMA = (pari_sp)leafcopy_avma(x, *AVMA);
    1880           0 :       return (GEN)*AVMA;
    1881             : 
    1882             :     /* one more special case */
    1883             :     case t_LIST:
    1884             :     {
    1885           0 :       long t = list_typ(x);
    1886           0 :       GEN y = cgetlist_avma(AVMA), z = list_data(x);
    1887           0 :       if (z) {
    1888           0 :         list_data(y) = gcopy_av0_canon(z, AVMA);
    1889           0 :         y[1] = evaltyp(t)|evallg(lg(z)-1);
    1890             :       } else {
    1891           0 :         list_data(y) = NULL;
    1892           0 :         y[1] = evaltyp(t);
    1893             :       }
    1894           0 :       return y;
    1895             :     }
    1896             :   }
    1897           0 :   y = cgetg_copy_avma(x, &lx, AVMA);
    1898           0 :   if (lontyp[tx] == 1) i = 1; else { y[1] = x[1]; i = 2; }
    1899           0 :   for (; i<lx; i++) gel(y,i) = gcopy_av0_canon(gel(x,i), AVMA);
    1900           0 :   return y;
    1901             : }
    1902             : 
    1903             : /* [copy_bin/bin_copy:] size (number of words) required for
    1904             :  * gcopy_av0_canon(x) */
    1905             : static long
    1906           0 : taille0_canon(GEN x)
    1907             : {
    1908           0 :   long i,n,lx, tx = typ(x);
    1909           0 :   switch(tx)
    1910             :   { /* non recursive types */
    1911           0 :     case t_INT: return signe(x)? lgefint(x): 0;
    1912             :     case t_REAL:
    1913             :     case t_STR:
    1914           0 :     case t_VECSMALL: return lg(x);
    1915             : 
    1916             :     /* one more special case */
    1917             :     case t_LIST:
    1918             :     {
    1919           0 :       GEN L = list_data(x);
    1920           0 :       return L? 3 + taille0_canon(L): 3;
    1921             :     }
    1922             :   }
    1923           0 :   n = lx = lg(x);
    1924           0 :   for (i=lontyp[tx]; i<lx; i++) n += taille0_canon(gel(x,i));
    1925           0 :   return n;
    1926             : }
    1927             : 
    1928             : /* [copy_bin/bin_copy:] size (number of words) required for gcopy_av0(x) */
    1929             : static long
    1930  1059788274 : taille0(GEN x)
    1931             : {
    1932  1059788274 :   long i,n,lx, tx = typ(x);
    1933  1059788274 :   switch(tx)
    1934             :   { /* non recursive types */
    1935             :     case t_INT:
    1936   712039800 :       lx = lgefint(x);
    1937   712039800 :       return lx == 2? 0: lx;
    1938             :     case t_LIST:
    1939             :     {
    1940          49 :       GEN L = list_data(x);
    1941          49 :       if (L && !list_nmax(x)) break; /* not finalized, deep copy */
    1942             :     }
    1943             :     /* else finalized: shallow */
    1944             :     case t_REAL:
    1945             :     case t_STR:
    1946             :     case t_VECSMALL:
    1947   125114432 :       return lg(x);
    1948             :   }
    1949   222634042 :   n = lx = lg(x);
    1950   222634042 :   for (i=lontyp[tx]; i<lx; i++) n += taille0(gel(x,i));
    1951   222633996 :   return n;
    1952             : }
    1953             : 
    1954             : static long
    1955  2831337670 : gsizeclone_i(GEN x)
    1956             : {
    1957  2831337670 :   long i,n,lx, tx = typ(x);
    1958  2831337670 :   switch(tx)
    1959             :   { /* non recursive types */
    1960  2636120581 :     case t_INT: lx = lgefint(x); return lx == 2? 0: lx;;
    1961             :     case t_REAL:
    1962             :     case t_STR:
    1963    34520774 :     case t_VECSMALL: return lg(x);
    1964             : 
    1965         697 :     case t_LIST: return 3;
    1966             :     default:
    1967   160695618 :       n = lx = lg(x);
    1968   160695618 :       for (i=lontyp[tx]; i<lx; i++) n += gsizeclone_i(gel(x,i));
    1969   160695607 :       return n;
    1970             :   }
    1971             : }
    1972             : 
    1973             : /* #words needed to clone x; t_LIST is a special case since list_data() is
    1974             :  * malloc'ed later, in list_internal_copy() */
    1975             : static long
    1976   208764201 : gsizeclone(GEN x) { return (typ(x) == t_INT)? lgefint(x): gsizeclone_i(x); }
    1977             : 
    1978             : long
    1979        6846 : gsizeword(GEN x)
    1980             : {
    1981             :   GEN L;
    1982        6846 :   if (typ(x) != t_LIST) return gsizeclone(x);
    1983             :   /* For t_LIST, return the actual list size, gsizeclone() is always 3 */
    1984         147 :   L = list_data(x);
    1985         147 :   return L? 3 + gsizeclone(L): 3;
    1986             : }
    1987             : long
    1988         147 : gsizebyte(GEN x) { return gsizeword(x) * sizeof(long); }
    1989             : 
    1990             : /* return a clone of x structured as a gcopy */
    1991             : GENbin*
    1992    51357020 : copy_bin(GEN x)
    1993             : {
    1994    51357020 :   long t = taille0(x);
    1995    51356990 :   GENbin *p = (GENbin*)pari_malloc(sizeof(GENbin) + t*sizeof(long));
    1996    51357096 :   pari_sp AVMA = (pari_sp)(GENbinbase(p) + t);
    1997    51357052 :   p->rebase = &shiftaddress;
    1998    51357052 :   p->len = t;
    1999    51357052 :   p->x   = gcopy_av0(x, &AVMA);
    2000    51357025 :   p->base= (GEN)AVMA; return p;
    2001             : }
    2002             : 
    2003             : /* same, writing t_INT in canonical native form */
    2004             : GENbin*
    2005           0 : copy_bin_canon(GEN x)
    2006             : {
    2007           0 :   long t = taille0_canon(x);
    2008           0 :   GENbin *p = (GENbin*)pari_malloc(sizeof(GENbin) + t*sizeof(long));
    2009           0 :   pari_sp AVMA = (pari_sp)(GENbinbase(p) + t);
    2010           0 :   p->rebase = &shiftaddress_canon;
    2011           0 :   p->len = t;
    2012           0 :   p->x   = gcopy_av0_canon(x, &AVMA);
    2013           0 :   p->base= (GEN)AVMA; return p;
    2014             : }
    2015             : 
    2016             : GEN
    2017   208757376 : gclone(GEN x)
    2018             : {
    2019   208757376 :   long i,lx,tx = typ(x), t = gsizeclone(x);
    2020   208757376 :   GEN y = newblock(t);
    2021   208757376 :   switch(tx)
    2022             :   { /* non recursive types */
    2023             :     case t_INT:
    2024   142030162 :       lx = lgefint(x);
    2025   142030162 :       y[0] = evaltyp(t_INT)|evallg(lx);
    2026   142030162 :       for (i=1; i<lx; i++) y[i] = x[i];
    2027   142030162 :       break;
    2028             :     case t_REAL:
    2029             :     case t_STR:
    2030             :     case t_VECSMALL:
    2031     6090513 :       lx = lg(x);
    2032     6090513 :       for (i=0; i<lx; i++) y[i] = x[i];
    2033     6090513 :       break;
    2034             : 
    2035             :     /* one more special case */
    2036             :     case t_LIST:
    2037         557 :       y[0] = evaltyp(t_LIST)|_evallg(3);
    2038         557 :       listassign(x, y);
    2039         557 :       break;
    2040             :     default: {
    2041    60636144 :       pari_sp AVMA = (pari_sp)(y + t);
    2042    60636144 :       lx = lg(x);
    2043    60636144 :       y[0] = x[0];
    2044    60636144 :       if (lontyp[tx] == 1) i = 1; else { y[1] = x[1]; i = 2; }
    2045    60636144 :       for (; i<lx; i++) gel(y,i) = gcopy_avma(gel(x,i), &AVMA);
    2046             :     }
    2047             :   }
    2048   208757376 :   setisclone(y); return y;
    2049             : }
    2050             : 
    2051             : void
    2052   748931265 : shiftaddress(GEN x, long dec)
    2053             : {
    2054   748931265 :   long i, lx, tx = typ(x);
    2055   748931265 :   if (is_recursive_t(tx))
    2056             :   {
    2057   222634102 :     if (tx == t_LIST)
    2058             :     {
    2059   748931263 :       if (!list_data(x) || list_nmax(x)) return; /* empty or finalized */
    2060             :       /* not finalized, update pointers  */
    2061             :     }
    2062   222634060 :     lx = lg(x);
    2063  1231065167 :     for (i=lontyp[tx]; i<lx; i++) {
    2064  1008431078 :       if (!x[i]) gel(x,i) = gen_0;
    2065             :       else
    2066             :       {
    2067   697581064 :         x[i] += dec;
    2068   697581064 :         shiftaddress(gel(x,i), dec);
    2069             :       }
    2070             :     }
    2071             :   }
    2072             : }
    2073             : 
    2074             : void
    2075           0 : shiftaddress_canon(GEN x, long dec)
    2076             : {
    2077           0 :   long i, lx, tx = typ(x);
    2078           0 :   switch(tx)
    2079             :   { /* non recursive types */
    2080             :     case t_INT: {
    2081             :       GEN y;
    2082           0 :       lx = lgefint(x); if (lx <= 3) return;
    2083           0 :       y = x + 2;
    2084           0 :       x = int_MSW(x);  if (x == y) return;
    2085           0 :       while (x > y) { lswap(*x, *y); x = int_precW(x); y++; }
    2086           0 :       break;
    2087             :     }
    2088             :     case t_REAL:
    2089             :     case t_STR:
    2090             :     case t_VECSMALL:
    2091           0 :       break;
    2092             : 
    2093             :     /* one more special case */
    2094             :     case t_LIST: {
    2095           0 :       GEN Lx = list_data(x);
    2096           0 :       if (Lx) {
    2097           0 :         pari_sp av = avma;
    2098           0 :         GEN L = (GEN)((long)Lx+dec);
    2099           0 :         shiftaddress_canon(L, dec);
    2100           0 :         list_data(x) = list_internal_copy(L, lg(L)); avma = av;
    2101             :       }
    2102             :     }
    2103             :     default:
    2104           0 :       lx = lg(x);
    2105           0 :       for (i=lontyp[tx]; i<lx; i++) {
    2106           0 :         if (!x[i]) gel(x,i) = gen_0;
    2107             :         else
    2108             :         {
    2109           0 :           x[i] += dec;
    2110           0 :           shiftaddress_canon(gel(x,i), dec);
    2111             :         }
    2112             :       }
    2113             :   }
    2114             : }
    2115             : 
    2116             : /********************************************************************/
    2117             : /**                                                                **/
    2118             : /**                INSERT DYNAMIC OBJECT IN STRUCTURE              **/
    2119             : /**                                                                **/
    2120             : /********************************************************************/
    2121             : GEN
    2122          28 : obj_reinit(GEN S)
    2123             : {
    2124          28 :   GEN s, T = leafcopy(S);
    2125          28 :   long a = lg(T)-1;
    2126          28 :   s = gel(T,a);
    2127          28 :   gel(T,a) = zerovec(lg(s)-1);
    2128          28 :   return T;
    2129             : }
    2130             : 
    2131             : GEN
    2132     1858115 : obj_init(long d, long n)
    2133             : {
    2134     1858115 :   GEN S = cgetg(d+2, t_VEC);
    2135     1858113 :   gel(S, d+1) = zerovec(n);
    2136     1858115 :   return S;
    2137             : }
    2138             : /* insert O in S [last position] at position K, return it */
    2139             : GEN
    2140     1099935 : obj_insert(GEN S, long K, GEN O)
    2141     1099935 : { return obj_insert_shallow(S, K, gclone(O)); }
    2142             : /* as obj_insert. WITHOUT cloning (for libpari, when creating a *new* obj
    2143             :  * from an existing one) */
    2144             : GEN
    2145     1104289 : obj_insert_shallow(GEN S, long K, GEN O)
    2146             : {
    2147     1104289 :   GEN o, v = gel(S, lg(S)-1);
    2148     1104289 :   if (typ(v) != t_VEC) pari_err_TYPE("obj_insert", S);
    2149     1104289 :   o = gel(v,K);
    2150     1104289 :   gel(v,K) = O; /*SIGINT: before unclone(o)*/
    2151     1104289 :   if (isclone(o)) gunclone(o);
    2152     1104289 :   return gel(v,K);
    2153             : }
    2154             : 
    2155             : /* Does S [last position] contain data at position K ? Return it, or NULL */
    2156             : GEN
    2157     2427837 : obj_check(GEN S, long K)
    2158             : {
    2159     2427837 :   GEN O, v = gel(S, lg(S)-1);
    2160     2427837 :   if (typ(v) != t_VEC || K >= lg(v)) pari_err_TYPE("obj_check", S);
    2161     2427837 :   O = gel(v,K); return isintzero(O)? NULL: O;
    2162             : }
    2163             : 
    2164             : GEN
    2165      741444 : obj_checkbuild(GEN S, long tag, GEN (*build)(GEN))
    2166             : {
    2167      741444 :   GEN O = obj_check(S, tag);
    2168      741444 :   if (!O)
    2169      634117 :   { pari_sp av = avma; O = obj_insert(S, tag, build(S)); avma = av; }
    2170      741437 :   return O;
    2171             : }
    2172             : 
    2173             : GEN
    2174       35155 : obj_checkbuild_prec(GEN S, long tag, GEN (*build)(GEN,long),
    2175             :   long (*pr)(GEN), long prec)
    2176             : {
    2177       35155 :   pari_sp av = avma;
    2178       35155 :   GEN w = obj_check(S, tag);
    2179       35155 :   if (!w || pr(w) < prec) w = obj_insert(S, tag, build(S, prec));
    2180       35155 :   avma = av; return gcopy(w);
    2181             : }
    2182             : GEN
    2183       11481 : obj_checkbuild_realprec(GEN S, long tag, GEN (*build)(GEN,long), long prec)
    2184       11481 : { return obj_checkbuild_prec(S,tag,build,gprecision,prec); }
    2185             : GEN
    2186         497 : obj_checkbuild_padicprec(GEN S, long tag, GEN (*build)(GEN,long), long prec)
    2187         497 : { return obj_checkbuild_prec(S,tag,build,padicprec_relative,prec); }
    2188             : 
    2189             : /* Reset S [last position], freeing all clones */
    2190             : void
    2191        9772 : obj_free(GEN S)
    2192             : {
    2193        9772 :   GEN v = gel(S, lg(S)-1);
    2194             :   long i;
    2195        9772 :   if (typ(v) != t_VEC) pari_err_TYPE("obj_free", S);
    2196       53844 :   for (i = 1; i < lg(v); i++)
    2197             :   {
    2198       44072 :     GEN o = gel(v,i);
    2199       44072 :     gel(v,i) = gen_0;
    2200       44072 :     gunclone_deep(o);
    2201             :   }
    2202        9772 : }
    2203             : 
    2204             : /*******************************************************************/
    2205             : /*                                                                 */
    2206             : /*                         STACK MANAGEMENT                        */
    2207             : /*                                                                 */
    2208             : /*******************************************************************/
    2209             : INLINE void
    2210  1886679791 : dec_gerepile(pari_sp *x, pari_sp av0, pari_sp av, pari_sp tetpil, size_t dec)
    2211             : {
    2212  1886679791 :   if (*x < av && *x >= av0)
    2213             :   { /* update address if in stack */
    2214  1563163315 :     if (*x < tetpil) *x += dec;
    2215           0 :     else pari_err_BUG("gerepile, significant pointers lost");
    2216             :   }
    2217  1886679791 : }
    2218             : 
    2219             : void
    2220      189367 : gerepileallsp(pari_sp av, pari_sp tetpil, int n, ...)
    2221             : {
    2222      189367 :   const pari_sp av0 = avma;
    2223      189367 :   const size_t dec = av-tetpil;
    2224             :   int i;
    2225      189367 :   va_list a; va_start(a, n);
    2226      189367 :   (void)gerepile(av,tetpil,NULL);
    2227      189367 :   for (i=0; i<n; i++) dec_gerepile((pari_sp*)va_arg(a,GEN*), av0,av,tetpil,dec);
    2228      189367 :   va_end(a);
    2229      189367 : }
    2230             : 
    2231             : /* Takes an array of pointers to GENs, of length n.
    2232             :  * Cleans up the stack between av and tetpil, updating those GENs. */
    2233             : void
    2234     5480764 : gerepilemanysp(pari_sp av, pari_sp tetpil, GEN* gptr[], int n)
    2235             : {
    2236     5480764 :   const pari_sp av0 = avma;
    2237     5480764 :   const size_t dec = av-tetpil;
    2238             :   int i;
    2239     5480764 :   (void)gerepile(av,tetpil,NULL);
    2240     5480764 :   for (i=0; i<n; i++) dec_gerepile((pari_sp*)gptr[i], av0, av, tetpil, dec);
    2241     5480764 : }
    2242             : 
    2243             : /* Takes an array of GENs (cast to longs), of length n.
    2244             :  * Cleans up the stack between av and tetpil, updating those GENs. */
    2245             : void
    2246    95866570 : gerepilecoeffssp(pari_sp av, pari_sp tetpil, long *g, int n)
    2247             : {
    2248    95866570 :   const pari_sp av0 = avma;
    2249    95866570 :   const size_t dec = av-tetpil;
    2250             :   int i;
    2251    95866570 :   (void)gerepile(av,tetpil,NULL);
    2252    95866570 :   for (i=0; i<n; i++,g++) dec_gerepile((pari_sp*)g, av0, av, tetpil, dec);
    2253    95866570 : }
    2254             : 
    2255             : static int
    2256           0 : dochk_gerepileupto(GEN av, GEN x)
    2257             : {
    2258             :   long i,lx,tx;
    2259           0 :   if (!isonstack(x)) return 1;
    2260           0 :   if (x > av)
    2261             :   {
    2262           0 :     pari_warn(warner,"bad object %Ps",x);
    2263           0 :     return 0;
    2264             :   }
    2265           0 :   tx = typ(x);
    2266           0 :   if (! is_recursive_t(tx)) return 1;
    2267             : 
    2268           0 :   lx = lg(x);
    2269           0 :   for (i=lontyp[tx]; i<lx; i++)
    2270           0 :     if (!dochk_gerepileupto(av, gel(x,i)))
    2271             :     {
    2272           0 :       pari_warn(warner,"bad component %ld in object %Ps",i,x);
    2273           0 :       return 0;
    2274             :     }
    2275           0 :   return 1;
    2276             : }
    2277             : /* check that x and all its components are out of stack, or have been
    2278             :  * created after av */
    2279             : int
    2280           0 : chk_gerepileupto(GEN x) { return dochk_gerepileupto(x, x); }
    2281             : 
    2282             : /* print stack between avma & av */
    2283             : void
    2284           0 : dbg_gerepile(pari_sp av)
    2285             : {
    2286           0 :   GEN x = (GEN)avma;
    2287           0 :   while (x < (GEN)av)
    2288             :   {
    2289           0 :     const long tx = typ(x), lx = lg(x);
    2290             :     GEN *a;
    2291             : 
    2292           0 :     pari_printf(" [%ld] %Ps:", x - (GEN)avma, x);
    2293           0 :     if (! is_recursive_t(tx)) { pari_putc('\n'); x += lx; continue; }
    2294           0 :     a = (GEN*)x + lontyp[tx]; x += lx;
    2295           0 :     for (  ; a < (GEN*)x; a++)
    2296             :     {
    2297           0 :       if (*a == gen_0)
    2298           0 :         pari_puts("  gen_0");
    2299           0 :       else if (*a == gen_1)
    2300           0 :         pari_puts("  gen_1");
    2301           0 :       else if (*a == gen_m1)
    2302           0 :         pari_puts("  gen_m1");
    2303           0 :       else if (*a == gen_2)
    2304           0 :         pari_puts("  gen_2");
    2305           0 :       else if (*a == gen_m2)
    2306           0 :         pari_puts("  gen_m2");
    2307           0 :       else if (*a == ghalf)
    2308           0 :         pari_puts("  ghalf");
    2309           0 :       else if (isclone(*a))
    2310           0 :         pari_printf("  %Ps (clone)", *a);
    2311             :       else
    2312           0 :         pari_printf("  %Ps [%ld]", *a, *a - (GEN)avma);
    2313           0 :       if (a+1 < (GEN*)x) pari_putc(',');
    2314             :     }
    2315           0 :     pari_printf("\n");
    2316             :   }
    2317           0 : }
    2318             : void
    2319           0 : dbg_gerepileupto(GEN q)
    2320             : {
    2321           0 :   err_printf("%Ps:\n", q);
    2322           0 :   dbg_gerepile((pari_sp) (q+lg(q)));
    2323           0 : }
    2324             : 
    2325             : GEN
    2326   404263653 : gerepile(pari_sp av, pari_sp tetpil, GEN q)
    2327             : {
    2328   404263653 :   const size_t dec = av - tetpil;
    2329   404263653 :   const pari_sp av0 = avma;
    2330             :   GEN x, a;
    2331             : 
    2332   404263653 :   if (dec == 0) return q;
    2333   334391946 :   if ((long)dec < 0) pari_err(e_MISC,"lbot>ltop in gerepile");
    2334             : 
    2335             :   /* dec_gerepile(&q, av0, av, tetpil, dec), saving 1 comparison */
    2336   334384719 :   if (q >= (GEN)av0 && q < (GEN)tetpil)
    2337   236568220 :     q = (GEN) (((pari_sp)q) + dec);
    2338             : 
    2339   334384719 :   for (x = (GEN)av, a = (GEN)tetpil; a > (GEN)av0; ) *--x = *--a;
    2340   334384719 :   avma = (pari_sp)x;
    2341  2538598956 :   while (x < (GEN)av)
    2342             :   {
    2343  1869821287 :     const long tx = typ(x), lx = lg(x);
    2344             : 
    2345  1869821287 :     if (! is_recursive_t(tx)) { x += lx; continue; }
    2346   439451874 :     a = x + lontyp[tx]; x += lx;
    2347   439451874 :     for (  ; a < x; a++) dec_gerepile((pari_sp*)a, av0, av, tetpil, dec);
    2348             :   }
    2349   334392950 :   return q;
    2350             : }
    2351             : 
    2352             : void
    2353           0 : fill_stack(void)
    2354             : {
    2355           0 :   GEN x = ((GEN)pari_mainstack->bot);
    2356           0 :   while (x < (GEN)avma) *x++ = 0xfefefefeUL;
    2357           0 : }
    2358             : 
    2359             : void
    2360           0 : debug_stack(void)
    2361             : {
    2362           0 :   pari_sp top = pari_mainstack->top, bot = pari_mainstack->bot;
    2363             :   GEN z;
    2364           0 :   err_printf("bot=0x%lx\ttop=0x%lx\tavma=0x%lx\n", bot, top, avma);
    2365           0 :   for (z = ((GEN)top)-1; z >= (GEN)avma; z--)
    2366           0 :     err_printf("%p:\t0x%lx\t%lu\n",z,*z,*z);
    2367           0 : }
    2368             : 
    2369             : void
    2370           0 : setdebugvar(long n) { DEBUGVAR=n; }
    2371             : 
    2372             : long
    2373           0 : getdebugvar(void) { return DEBUGVAR; }
    2374             : 
    2375             : long
    2376           7 : getstack(void) { return pari_mainstack->top-avma; }
    2377             : 
    2378             : /*******************************************************************/
    2379             : /*                                                                 */
    2380             : /*                               timer_delay                             */
    2381             : /*                                                                 */
    2382             : /*******************************************************************/
    2383             : 
    2384             : #if defined(USE_CLOCK_GETTIME)
    2385             : #if defined(_POSIX_THREAD_CPUTIME)
    2386             : static THREAD clockid_t time_type = CLOCK_THREAD_CPUTIME_ID;
    2387             : #else
    2388             : static const THREAD clockid_t time_type = CLOCK_PROCESS_CPUTIME_ID;
    2389             : #endif
    2390             : static void
    2391             : pari_init_timer(void)
    2392             : {
    2393             : #if defined(_POSIX_THREAD_CPUTIME)
    2394             :   time_type = CLOCK_PROCESS_CPUTIME_ID;
    2395             : #endif
    2396             : }
    2397             : 
    2398             : void
    2399             : timer_start(pari_timer *T)
    2400             : {
    2401             :   struct timespec t;
    2402             :   clock_gettime(time_type,&t);
    2403             :   T->us = t.tv_nsec / 1000;
    2404             :   T->s  = t.tv_sec;
    2405             : }
    2406             : #elif defined(USE_GETRUSAGE)
    2407             : #ifdef RUSAGE_THREAD
    2408             : static THREAD int rusage_type = RUSAGE_THREAD;
    2409             : #else
    2410             : static const THREAD int rusage_type = RUSAGE_SELF;
    2411             : #endif /*RUSAGE_THREAD*/
    2412             : static void
    2413        1519 : pari_init_timer(void)
    2414             : {
    2415             : #ifdef RUSAGE_THREAD
    2416        1519 :   rusage_type = RUSAGE_SELF;
    2417             : #endif
    2418        1519 : }
    2419             : 
    2420             : void
    2421      502215 : timer_start(pari_timer *T)
    2422             : {
    2423             :   struct rusage r;
    2424      502215 :   getrusage(rusage_type,&r);
    2425      502215 :   T->us = r.ru_utime.tv_usec;
    2426      502215 :   T->s  = r.ru_utime.tv_sec;
    2427      502215 : }
    2428             : #elif defined(USE_FTIME)
    2429             : 
    2430             : static void
    2431             : pari_init_timer(void) { }
    2432             : 
    2433             : void
    2434             : timer_start(pari_timer *T)
    2435             : {
    2436             :   struct timeb t;
    2437             :   ftime(&t);
    2438             :   T->us = ((long)t.millitm) * 1000;
    2439             :   T->s  = t.time;
    2440             : }
    2441             : 
    2442             : #else
    2443             : 
    2444             : static void
    2445             : _get_time(pari_timer *T, long Ticks, long TickPerSecond)
    2446             : {
    2447             :   T->us = (long) ((Ticks % TickPerSecond) * (1000000. / TickPerSecond));
    2448             :   T->s  = Ticks / TickPerSecond;
    2449             : }
    2450             : 
    2451             : # ifdef USE_TIMES
    2452             : static void
    2453             : pari_init_timer(void) { }
    2454             : 
    2455             : void
    2456             : timer_start(pari_timer *T)
    2457             : {
    2458             : # ifdef _SC_CLK_TCK
    2459             :   long tck = sysconf(_SC_CLK_TCK);
    2460             : # else
    2461             :   long tck = CLK_TCK;
    2462             : # endif
    2463             :   struct tms t; times(&t);
    2464             :   _get_time(T, t.tms_utime, tck);
    2465             : }
    2466             : # elif defined(_WIN32)
    2467             : static void
    2468             : pari_init_timer(void) { }
    2469             : 
    2470             : void
    2471             : timer_start(pari_timer *T)
    2472             : { _get_time(T, win32_timer(), 1000); }
    2473             : # else
    2474             : #  include <time.h>
    2475             : #  ifndef CLOCKS_PER_SEC
    2476             : #   define CLOCKS_PER_SEC 1000000 /* may be false on YOUR system */
    2477             : #  endif
    2478             : static void
    2479             : pari_init_timer(void) { }
    2480             : 
    2481             : void
    2482             : timer_start(pari_timer *T)
    2483             : { _get_time(T, clock(), CLOCKS_PER_SEC); }
    2484             : # endif
    2485             : #endif
    2486             : 
    2487             : static long
    2488       67777 : timer_aux(pari_timer *T, pari_timer *U)
    2489             : {
    2490       67777 :   long s = T->s, us = T->us; timer_start(U);
    2491       67777 :   return 1000 * (U->s - s) + (U->us - us + 500) / 1000;
    2492             : }
    2493             : /* return delay, reset timer */
    2494             : long
    2495       66256 : timer_delay(pari_timer *T) { return timer_aux(T, T); }
    2496             : /* return delay, don't reset timer */
    2497             : long
    2498        1521 : timer_get(pari_timer *T) { pari_timer t; return timer_aux(T, &t); }
    2499             : 
    2500             : static void
    2501           0 : timer_vprintf(pari_timer *T, const char *format, va_list args)
    2502             : {
    2503           0 :   out_puts(pariErr, "Time ");
    2504           0 :   out_vprintf(pariErr, format,args);
    2505           0 :   out_printf(pariErr, ": %ld\n", timer_delay(T));
    2506           0 :   pariErr->flush();
    2507           0 : }
    2508             : void
    2509           0 : timer_printf(pari_timer *T, const char *format, ...)
    2510             : {
    2511           0 :   va_list args; va_start(args, format);
    2512           0 :   timer_vprintf(T, format, args);
    2513           0 :   va_end(args);
    2514           0 : }
    2515             : 
    2516             : long
    2517           0 : timer(void)  { static THREAD pari_timer T; return timer_delay(&T);}
    2518             : long
    2519        3008 : gettime(void)  { static THREAD pari_timer T; return timer_delay(&T);}
    2520             : 
    2521             : static THREAD pari_timer timer2_T, abstimer_T;
    2522             : long
    2523           0 : timer2(void) {  return timer_delay(&timer2_T);}
    2524             : void
    2525           0 : msgtimer(const char *format, ...)
    2526             : {
    2527           0 :   va_list args; va_start(args, format);
    2528           0 :   timer_vprintf(&timer2_T, format, args);
    2529           0 :   va_end(args);
    2530           0 : }
    2531             : long
    2532        1519 : getabstime(void)  { return timer_get(&abstimer_T);}
    2533             : #if defined(USE_CLOCK_GETTIME) || defined(USE_GETTIMEOFDAY) \
    2534             :  || defined(USE_FTIMEFORWALLTIME)
    2535             : static GEN
    2536           0 : timetoi(ulong s, ulong m)
    2537             : {
    2538           0 :   pari_sp av = avma;
    2539           0 :   GEN r = addiu(muliu(utoi(s), 1000), m);
    2540           0 :   return gerepileuptoint(av, r);
    2541             : }
    2542             : #endif
    2543             : GEN
    2544           0 : getwalltime(void)
    2545             : {
    2546             : #if defined(USE_CLOCK_GETTIME)
    2547             :   struct timespec t;
    2548             :   if (!clock_gettime(CLOCK_REALTIME,&t))
    2549             :     return timetoi(t.tv_sec, (t.tv_nsec + 500000)/1000000);
    2550             : #elif defined(USE_GETTIMEOFDAY)
    2551             :   struct timeval tv;
    2552           0 :   if (!gettimeofday(&tv, NULL))
    2553           0 :     return timetoi(tv.tv_sec, (tv.tv_usec + 500)/1000);
    2554             : #elif defined(USE_FTIMEFORWALLTIME)
    2555             :   struct timeb tp;
    2556             :   ftime(&tp); return timetoi(tp.time, tp.millitm);
    2557             : #endif
    2558           0 :   return utoi(getabstime());
    2559             : }
    2560             : 
    2561             : /*******************************************************************/
    2562             : /*                                                                 */
    2563             : /*                   FUNCTIONS KNOWN TO THE ANALYZER               */
    2564             : /*                                                                 */
    2565             : /*******************************************************************/
    2566             : GEN
    2567           7 : pari_version(void)
    2568             : {
    2569           7 :   const ulong mask = (1UL<<PARI_VERSION_SHIFT) - 1;
    2570           7 :   ulong major, minor, patch, n = paricfg_version_code;
    2571           7 :   patch = n & mask; n >>= PARI_VERSION_SHIFT;
    2572           7 :   minor = n & mask; n >>= PARI_VERSION_SHIFT;
    2573           7 :   major = n;
    2574           7 :   if (*paricfg_vcsversion) {
    2575           7 :     const char *ver = paricfg_vcsversion;
    2576           7 :     const char *s = strchr(ver, '-');
    2577             :     char t[8];
    2578           7 :     const long len = s-ver;
    2579             :     GEN v;
    2580           7 :     if (!s || len > 6) pari_err_BUG("pari_version()"); /* paranoia */
    2581           7 :     memcpy(t, ver, len); t[len] = 0;
    2582           7 :     v = cgetg(6, t_VEC);
    2583           7 :     gel(v,1) = utoi(major);
    2584           7 :     gel(v,2) = utoi(minor);
    2585           7 :     gel(v,3) = utoi(patch);
    2586           7 :     gel(v,4) = stoi( atoi(t) );
    2587           7 :     gel(v,5) = strtoGENstr(s+1);
    2588           7 :     return v;
    2589             :   } else {
    2590           0 :     GEN v = cgetg(4, t_VEC);
    2591           0 :     gel(v,1) = utoi(major);
    2592           0 :     gel(v,2) = utoi(minor);
    2593           0 :     gel(v,3) = utoi(patch);
    2594           0 :     return v;
    2595             :   }
    2596             : }
    2597             : 
    2598             : /* List of GP functions: generated from the description system.
    2599             :  * Format (struct entree) :
    2600             :  *   char *name   : name (under GP).
    2601             :  *   ulong valence: (EpNEW, EpALIAS,EpVAR, EpINSTALL)|EpSTATIC
    2602             :  *   void *value  : For PREDEFINED FUNCTIONS: C function to call.
    2603             :  *                  For USER FUNCTIONS: pointer to defining data (block) =
    2604             :  *                   entree*: NULL, list of entree (arguments), NULL
    2605             :  *                   char*  : function text
    2606             :  *   long menu    : which help section do we belong to
    2607             :  *                   1: Standard monadic or dyadic OPERATORS
    2608             :  *                   2: CONVERSIONS and similar elementary functions
    2609             :  *                   3: functions related to COMBINATORICS
    2610             :  *                   4: TRANSCENDENTAL functions, etc.
    2611             :  *   char *code   : GP prototype, aka Parser Code (see libpari's manual)
    2612             :  *                  if NULL, use valence instead.
    2613             :  *   char *help   : short help text (init to NULL).
    2614             :  *   void *pvalue : push_val history.
    2615             :  *   long arity   : maximum number of arguments.
    2616             :  *   entree *next : next entree (init to NULL, used in hashing code). */
    2617             : #include "init.h"
    2618             : #include "default.h"

Generated by: LCOV version 1.11