Testing latest pari + WASM + node.js... and it works?! Wow.
License: GPL3
ubuntu2004
\def\TITLE{Developer's Guide to the PARI library}1\input parimacro.tex23% START TYPESET4\begintitle5\vskip 2.5truecm6\centerline{\mine Developer's Guide}7\vskip 1.truecm8\centerline{\mine to}9\vskip 1.truecm10\centerline{\mine the PARI library}11\vskip 1.truecm12\centerline{\sectiontitlebf (version \vers)}13\vskip 1.truecm14\authors15\endtitle1617\copyrightpage18\tableofcontents19\openin\std=develop.aux20\ifeof\std21\else22\input develop.aux23\fi24\chapno=02526\chapter{Work in progress}2728This draft documents private internal functions and structures for hard-core29PARI developers. Anything in here is liable to change on short notice. Don't30use anything in the present document, unless you are implementing new31features for the PARI library. Try to fix the interfaces before using them,32or document them in a better way.33If you find an undocumented hack somewhere, add it here.3435Hopefully, this will eventually document everything that we buried in36\kbd{paripriv.h} or even more private header files like \kbd{anal.h}.37Possibly, even implementation choices! Way to go.3839\section{The type \typ{CLOSURE}}\kbdsidx{t_CLOSURE}\sidx{closure}40This type holds closures and functions in compiled form, so is deeply41linked to the internals of the GP compiler and evaluator.42The length of this type can be $6$, $7$ or $8$ depending whether the43object is an ``inline closure'', a ``function'' or a ``true closure''.4445A function is a regular GP function. The GP input line is treated as a46function of arity $0$.4748A true closure is a GP function defined in a nonempty lexical context.4950An inline closure is a closure that appears in the code without51the preceding \kbd{->} token. They are generally attached to the prototype52code 'E' and 'I'. Inline closures can only exist as data of other closures,53see below.5455In the following example,56\bprog57f(a=Euler)=x->sin(x+a);58g=f(Pi/2);59plot(x=0,2*Pi,g(x))60@eprog\noindent61\kbd{f} is a function, \kbd{g} is a true closure and both \kbd{Euler} and62\kbd{g(x)} are inline closures.6364This type has a second codeword \kbd{z[1]}, which is the arity of the65function or closure. This is zero for inline closures. To access it, use6667\fun{long}{closure_arity}{GEN C}6869\item \kbd{z[2]} points to a \typ{STR} which holds the opcodes. To access it, use7071\fun{GEN}{closure_get_code}{GEN C}.7273\fun{const char *}{closure_codestr}{GEN C} returns as an array of \kbd{char}74starting at $1$.7576\item \kbd{z[3]} points to a \typ{VECSMALL} which holds the operands of the opcodes.77To access it, use7879\fun{GEN}{closure_get_oper}{GEN C}8081\item \kbd{z[4]} points to a \typ{VEC} which hold the data referenced by the82\kbd{pushgen} opcodes, which can be \typ{CLOSURE}, and in particular83inline closures. To access it, use8485\fun{GEN}{closure_get_data}{GEN C}8687\item \kbd{z[5]} points to a \typ{VEC} which hold extra data needed for88error-reporting and debugging. See \secref{se:dbgclosure} for details.89To access it, use9091\fun{GEN}{closure_get_dbg}{GEN C}9293Additionally, for functions and true closures,9495\item \kbd{z[6]} usually points to a \typ{VEC} with two components which are \typ{STR}.96The first one displays the list of arguments of the closure without the97enclosing parentheses, the second one the GP code of the function at the98right of the \kbd{->} token. They are used to display the closure, either in99implicit or explicit form. However for closures that were not generated from GP100code, \kbd{z[6]} can point to a \typ{STR} instead. To access it, use101102\fun{GEN}{closure_get_text}{GEN C}103104Additionally, for true closure,105106\item \kbd{z[7]} points to a \typ{VEC} which holds the values of all lexical107variables defined in the scope the closure was defined. To access it, use108109\fun{GEN}{closure_get_frame}{GEN C}110111\subsec{Debugging information in closure}\label{se:dbgclosure}112113Every \typ{CLOSURE} object \kbd{z} has a component \kbd{dbg=z[5]}114which hold extra data needed for error-reporting and debugging.115The object \kbd{dbg} is a \typ{VEC} with $3$ components:116117\kbd{dbg[1]} is a \typ{VECSMALL} of the same length than \kbd{z[3]}. For each118opcode, it holds the position of the corresponding GP source code in the119strings stored in \kbd{z[6]} for function or true closures, positive indices120referring to the second strings, and negative indices referring to the first121strings, the last element being indexed as $-1$. For inline closures, the122string of the parent function or true closure is used instead.123124\kbd{dbg[2]} is a \typ{VECSMALL} that lists opcodes index where new lexical125local variables are created. The value $0$ denotes the position before the126first offset and variables created by the prototype code 'V'.127128\kbd{dbg[3]} is a \typ{VEC} of \typ{VECSMALL}s that give the list of129\kbd{entree*} of the lexical local variables created at a given index in130\kbd{dbg[2]}.131132\section{The type \typ{LIST}}\kbdsidx{t_LIST}\sidx{list} This type needs to go133through various hoops to support GP's inconvenient memory model. Don't134use \typ{LIST}s in pure library mode, reimplement ordinary lists! This135dynamic type is implemented by a \kbd{GEN} of length 3: two codewords and a136vector containing the actual entries. In a normal setup (a finished list,137ready to be used),138139\item the vector is malloc'ed, so that it can be realloc'ated without moving140the parent \kbd{GEN}.141142\item all the entries are clones, possibly with cloned subcomponents; they143must be deleted with \tet{gunclone_deep}, not \tet{gunclone}.144145The following macros are proper lvalues and access the components146147\fun{long}{list_nmax}{GEN L}: current maximal number of elements. This grows148as needed.149150\fun{GEN}{list_data}{GEN L}: the elements. If \kbd{v = list\_data(L)}, then151either \kbd{v} is \kbd{NULL} (empty list) or \kbd{l = lg(v)} is defined, and152the elements are \kbd{v[1]}, \dots, \kbd{v[l-1]}.153154In most \kbd{gerepile} scenarios, the list components are not inspected155and a shallow copy of the malloc'ed vector is made. The functions156\kbd{gclone}, \kbd{copy\_bin\_canon} are exceptions, and make a full copy of157the list.158159The main problem with lists is to avoid memory leaks; in the above setup,160a statement like \kbd{a = List(1)} would already leak memory, since161\kbd{List(1)} allocates memory, which is cloned (second allocation) when162assigned to \kbd{a}; and the original list is lost. The solution we163implemented is164165\item to create anonymous lists (from \kbd{List}, \kbd{gtolist},166\kbd{concat} or \kbd{vecsort}) entirely on the stack, \emph{not} as described167above, and to set \kbd{list\_nmax} to $0$. Such a list is not yet proper and168trying to append elements to it fails:169\bprog170? listput(List(),1)171*** variable name expected: listput(List(),1)172*** ^----------------173@eprog\noindent174If we had been malloc'ing memory for the175\kbd{List([1,2,3])}, it would have leaked already.176177\item as soon as a list is assigned to a variable (or a component thereof)178by the GP evaluator, the assigned list is converted to the proper format179(with \kbd{list\_nmax} set) previously described.180181\fun{GEN}{listcopy}{GEN L} return a full copy of the \typ{LIST}~\kbd{L},182allocated on the \emph{stack} (hence \kbd{list\_nmax} is $0$). Shortcut for183\kbd{gcopy}.184185\fun{GEN}{mklistcopy}{GEN x} returns a list with a single element $x$,186allocated on the stack. Used to implement most cases of \kbd{gtolist}187(except vectors and lists).188189A typical low-level construct:190\bprog191long l;192/* assume L is a t_LIST */193L = list_data(L); /* discard t_LIST wrapper */194l = L? lg(L): 1;195for (i = 1; i < l; i++) output( gel(L, i) );196for (i = 1; i < l; i++) gel(L, i) = gclone( ... );197@eprog198199\subsec{Maps as Lists}200201GP's maps are implemented on top of \typ{LIST}s so as to benefit from202their peculiar memory models. Lists thus come in two subtypes: \typ{LIST\_RAW}203(actual lists) and \typ{LIST\_MAP} (a map).204205\fun{GEN}{mklist_typ}{long t} create a list of subtype $t$.206\fun{GEN}{mklist}{void} is an alias for207\bprog208mklist_typ(t_LIST_RAW);209@eprog210and211\fun{GEN}{mkmap}{void} is an alias for212\bprog213mklist_typ(t_LIST_MAP);214@eprog215216\fun{long}{list_typ}{GEN L} return the list subtype, either \typ{LIST\_RAW} or217\typ{LIST\_MAP}.218219\fun{void}{listpop}{GEN L, long index} as \kbd{listpop0},220assuming that $L$ is a \typ{LIST\_RAW}.221222\fun{GEN}{listput}{GEN list, GEN object, long index} as \kbd{listput0},223assuming that $L$ is a \typ{LIST\_RAW}.224225\fun{GEN}{mapdomain}{GEN T} vector of keys of the map $T$.226227\fun{GEN}{mapdomain_shallow}{GEN T} shallow version of \kbd{mapdomain}.228229\fun{GEN}{maptomat}{GEN T} convert a map to a factorization matrix.230231\fun{GEN}{maptomat_shallow}{GEN T} shallow version of \kbd{maptomat}.232233\section{Protection of noninterruptible code}234235GP allows the user to interrupt a computation by issuing SIGINT236(usually by entering control-C) or SIGALRM (usually using alarm()).237To avoid such interruption to occurs in section of code which are not238reentrant (in particular \kbd{malloc} and \kbd{free})239the following mechanism is provided:240241\fun{}{BLOCK_SIGINT_START}{}242Start a noninterruptible block code. Block both \kbd{SIGINT} and \kbd{SIGARLM}.243244\fun{}{BLOCK_SIGALRM_START}{}245Start a noninterruptible block code. Block only \kbd{SIGARLM}.246This is used in the \kbd{SIGINT} handler itself to delay an eventual pending247alarm.248249\fun{}{BLOCK_SIGINT_END}{}250End a noninterruptible block code251252The above macros make use of the following global variables:253254\tet{PARI_SIGINT_block}: set to $1$ (resp. $2$) by \kbd{BLOCK\_SIGINT\_START}255(resp. \kbd{BLOCK\_SIGALRM\_START}).256257\tet{PARI_SIGINT_pending}: Either $0$ (no signal was blocked), \kbd{SIGINT}258(\kbd{SIGINT} was blocked) or \kbd{SIGALRM} (\kbd{SIGALRM} was blocked).259This need to be set by the signal handler.260261Within a block, an automatic variable \kbd{int block} is defined which262records the value of \kbd{PARI\_SIGINT\_block} when entering the block.263264\subsec{Multithread interruptions}265266To support multithreaded programs, \kbd{BLOCK\_SIGINT\_START} and267\kbd{BLOCK\_SIGALRM\_START} call \kbd{MT\_SIGINT\_BLOCK(block)}, and268\kbd{BLOCK\_SIGINT\_END} calls \kbd{MT\_SIGINT\_UNBLOCK(block)}.269270\tet{MT_SIGINT_BLOCK} and \tet{MT_SIGINT_UNBLOCK} are defined by the271multithread engine. They can calls the following public functions defined by272the multithread engine.273274\fun{void}{mt_sigint_block}{void}275276\fun{void}{mt_sigint_unblock}{void}277278In practice this mechanism is used by the POSIX thread engine to protect against279asychronous cancellation.280281\section{$\F_{l^2}$ field for small primes $l$}282Let $l>2$ be a prime \kbd{ulong}. A \kbd{Fl2} is an element of the finite283field $\F_{l^2}$ represented (currently) by a \kbd{Flx} of degree at most $1$284modulo a polynomial of the form $x^2-D$ for some non square $0\leq D<p$.285Below \kbd{pi} denotes the pseudo inverse of \kbd{p}, see \kbd{Fl\_mul\_pre}286287\fun{int}{Fl2_equal1}{GEN x} return $1$ if $x=1$, else return $0$.288289\fun{GEN}{Fl2_mul_pre}{GEN x, GEN y, ulong D, ulong p, ulong pi} return $x\*y$.290291\fun{GEN}{Fl2_sqr_pre}{GEN x, ulong D, ulong p, ulong pi} return $x^2$.292293\fun{GEN}{Fl2_inv_pre}{GEN x, ulong D, ulong p, ulong pi} return $x^{-1}$.294295\fun{GEN}{Fl2_pow_pre}{GEN x, GEN n, ulong D, ulong p, ulong pi} return296$x^n$.297298\fun{GEN}{Fl2_sqrtn_pre}{GEN a, GEN n, ulong D, ulong p, ulong pi, GEN *zeta}299$n$-th root, as \kbd{Flxq\_sqrtn}.300301\fun{GEN}{Fl2_norm_pre}{GEN x, GEN n, ulong D, ulong p, ulong pi} return the302norm of $x$.303304\fun{GEN}{Flx_Fl2_eval_pre}{GEN P, GEN x, ulong D, ulong p, ulong pi}305return $P(x)$.306307\section{Public functions useless outside of GP context}308309These functions implement GP functionality for which the C language or310other libpari routines provide a better equivalent; or which are so tied311to the \kbd{gp} interpreter as to be virtually useless in \kbd{libpari}. Some312may be generated by \kbd{gp2c}. We document them here for completeness.313314\subsec{Conversions}315316\fun{GEN}{toser_i}{GEN x} internal shallow function, used to implement317automatic conversions to power series in GP (as in \kbd{cos(x)}).318Converts a \typ{POL} or a \typ{RFRAC} to a \typ{SER} in the same variable and319precision \kbd{precdl} (the global variable corresponding to320\kbd{seriesprecision}). Returns $x$ itself for a \typ{SER}, and \kbd{NULL}321for other argument types. The fact that it uses a global variable makes it322awkward whenever you're not implementing a new transcendental function in GP.323Use \tet{RgX_to_ser} or \tet{rfrac_to_ser} for a fast clean alternative to324\kbd{gtoser}.325326\fun{GEN}{listinit}{GEN x} a \typ{LIST} (from \kbd{List} or \kbd{Map}) may327exist in two different forms due to GP memory model:328329\item an ordinary \emph{read-only} copy on the PARI stack (as produced by330\kbd{gtolist} or \kbd{gtomap}) to which one may not assign elements331(\kbd{listput} will fail) unless the list is empty.332333\item a feature-complete GP list using (malloc'ed) \kbd{block}s to334allow dynamic insertions. An empty list is automaticaly promoted to this335status on first insertion.336337The \kbd{listinit} function creates a copy of existing \typ{SER} $x$ and338makes sure it is of the second kind. Variants of this are automatically339called by \kbd{gp} when assigning a \typ{LIST} to a GP variable; the340mecanism avoid memory leaks when creating a constant list, e.g.341\kbd{List([1,2,3])} (read-only), without assigning it to a variable. Whereas342after \kbd{L = List([1,2,3])} (GP list), we keep a pointer to the object and343may delete it when $L$ goes out of scope.344345This \kbd{libpari} function allows \kbd{gp2c} to simulate this process by346generating \kbd{listinit} calls at appropriate places.347348\subsec{Output}349350\fun{void}{print0}{GEN g, long flag} internal function underlying the351\kbd{print} GP function. Prints the entries of the \typ{VEC} $g$, one by one,352without any separator; entries of type \typ{STR} are printed without enclosing353quotes. \fl is one of \tet{f_RAW}, \tet{f_PRETTYMAT} or \tet{f_TEX}, using the354current default output context.355356\fun{void}{out_print0}{PariOUT *out, const char *sep, GEN g, long flag} as357\tet{print0}, using output context \kbd{out} and separator \kbd{sep} between358successive entries (no separator if \kbd{NULL}).359360\fun{void}{printsep}{const char *s, GEN g, long flag} \tet{out_print0} on361\tet{pariOut} followed by a newline.362363\fun{void}{printsep1}{const char *s, GEN g, long flag} \tet{out_print0} on364\tet{pariOut}.365366\fun{char*}{pari_sprint0}{const char *s, GEN g, long flag} displays $s$,367then \kbd{print0(g, flag)}.368369\fun{void}{print}{GEN g} equivalent to \kbd{print0(g, f\_RAW)}, followed370by a \kbd{\bs n} then an \kbd{fflush}.371372\fun{void}{printp}{GEN g} equivalent to \kbd{print0(g, f\_PRETTYMAT)},373followed by a \kbd{\bs n} then an \kbd{fflush}.374375\fun{void}{print1}{GEN g} as above, without the \kbd{\bs n}. Use376\tet{pari_printf} or \tet{output} instead.377378\fun{void}{printtex}{GEN g} equivalent to \kbd{print0(g, t\_TEX)}, followed379by a \kbd{\bs n} then an \kbd{fflush}. Use \tet{GENtoTeXstr} and380\tet{pari_printf} instead.381382\fun{void}{write0}{const char *s, GEN g}383384\fun{void}{write1}{const char *s, GEN g} use \kbd{fprintf}385386\fun{void}{writetex}{const char *s, GEN g} use \tet{GENtoTeXstr} and387\kbd{fprintf}.388389\fun{void}{printf0}{GEN fmt, GEN args} use \tet{pari_printf}.390391\fun{GEN}{strprintf}{GEN fmt, GEN args} use \tet{pari_sprintf}.392393\subsec{Input}394395\kbd{gp}'s input is read from the stream \tet{pari_infile}, which is changed396using397398\fun{FILE*}{switchin}{const char *name}399400Note that this function is quite complicated, maintaining stacks of files401to allow smooth error recovery and \kbd{gp} interaction. You will be better402off using \tet{gp_read_file}.403404\subsec{Control flow statements}405406\fun{GEN}{break0}{long n}. Use the C control statement \kbd{break}. Since407\kbd{break(2)} is invalid in C, either rework your code or use \kbd{goto}.408409\fun{GEN}{next0}{long n}. Use the C control statement \kbd{continue}. Since410\kbd{continue(2)} is invalid in C, either rework your code or use \kbd{goto}.411412\fun{GEN}{return0}{GEN x}. Use \kbd{return}!413414\fun{void}{error0}{GEN g}. Use \kbd{pari\_err(e\_USER,)}415416\fun{void}{warning0}{GEN g}. Use \kbd{pari\_warn(e\_USER,)}417418\subsec{Accessors}419420\fun{GEN}{vecslice0}{GEN A, long a, long b} implements $A[a..b]$.421422\fun{GEN}{matslice0}{GEN A, long a, long b, long c, long d}423implements $A[a..b, c..d]$.424425\subsec{Iterators}426427\fun{GEN}{apply0}{GEN f, GEN A} gp wrapper calling \tet{genapply}, where $f$428is a \typ{CLOSURE}, applied to $A$. Use \kbd{genapply} or a standard C loop.429430\fun{GEN}{select0}{GEN f, GEN A} gp wrapper calling \tet{genselect}, where $f$431is a \typ{CLOSURE} selecting from $A$. Use \kbd{genselect} or a standard C loop.432433\fun{GEN}{vecapply}{void *E, GEN (*f)(void* E, GEN x), GEN x} implements434\kbd{[a(x)|x<-b]}.435436\fun{GEN}{veccatapply}{void *E, GEN (*f)(void* E, GEN x), GEN x} implements437\kbd{concat([a(x)|x<-b])} which used to implement \kbd{[a0(x,y)|x<-b;y<-c(b)]}438which is equal to \kbd{concat([[a0(x,y)|y<-c(b)]|x<-b])}.439440\fun{GEN}{vecselect}{void *E, long (*f)(void* E, GEN x), GEN A}441implements \kbd{[x<-b,c(x)]}.442443\fun{GEN}{vecselapply}{void *Epred, long (*pred)(void* E, GEN x), void *Efun, GEN (*fun)(void* E, GEN x), GEN A}444implements \kbd{[a(x)|x<-b,c(x)]}.445446\subsec{Local precision}447448These functions allow to change \kbd{realprecision} locally when449calling the GP interpretor.450451\fun{void}{push_localprec}{long p} set the local precision to $p$.452453\fun{void}{push_localbitprec}{long b} set the local precision to $b$ bits.454455\fun{void}{pop_localprec}{void} reset the local precision to the previous456value.457458\fun{long}{get_localprec}{void} returns the current local precision.459460\fun{long}{get_localbitprec}{void} returns the current local precision in bits.461462\fun{void}{localprec}{long p} trivial wrapper around \kbd{push\_localprec}463(sanity checks and convert from decimal digits to a number of codewords).464Use \kbd{push\_localprec}.465466\fun{void}{localbitprec}{long p}467trivial wrapper around \kbd{push\_localbitprec}468(sanity checks). Use \kbd{push\_localbitprec}.469470These two function are used to implement \kbd{getlocalprec} and471\kbd{getlocalbitprec} for the GP interpreter and essentially return their472argument (the current dynamic precision, respectively in bits or as a473\kbd{prec} word count):474475\fun{long}{getlocalbitprec}{long bit}476477\fun{long}{getlocalprec}{long prec}478479480\subsec{Functions related to the GP evaluator}481482The prototype code \kbd{C} instructs the GP compiler to save the current483lexical context (pairs made of a lexical variable name and its value)484in a \kbd{GEN}, called \kbd{pack} in the sequel. This \kbd{pack} can be used485to evaluate expressions in the corresponding lexical context, providing it is486current.487488\fun{GEN}{localvars_read_str}{const char *s, GEN pack} evaluate the string $s$489in the lexical context given by \kbd{pack}. Used by \tet{geval_gp} in GP490to implement the behavior below:491\bprog492? my(z=3);eval("z=z^2");z493%1 = 9494@eprog495496\fun{long}{localvars_find}{GEN pack, entree *ep} does \kbd{pack} contain497a pair whose variable corresponds to \kbd{ep}? If so, where is the498corresponding value? (returns an offset on the value stack).499500\subsec{Miscellaneous}501502\fun{char*}{os_getenv}{const char *s} either calls \kbd{getenv}, or directly503return \kbd{NULL} if the \kbd{libc} does not provide it. Use \tet{getenv}.504505\fun{sighandler_t}{os_signal}{int sig, pari_sighandler_t fun} after a506\bprog507typedef void (*pari_sighandler_t)(int);508@eprog\noindent509(private type, not exported). Installs signal handler \kbd{fun} for510signal \kbd{sig}, using \tet{sigaction} with flag \tet{SA_NODEFER}. If511\kbd{sigaction} is not available use \tet{signal}. If even the latter is not512available, just return \tet{SIG_IGN}. Use \tet{sigaction}.513514\section{Embedded GP interpretor}515These function provide a simplified interface to embed a GP516interpretor in a program.517518\fun{void}{gp_embedded_init}{long rsize, long vsize}519Initialize the GP interpretor (like \kbd{pari\_init} does) with520\kbd{parisize=rsize} \kbd{rsize} and \kbd{parisizemax=vsize}.521522\fun{char *}{gp_embedded}{const char *s}523Evaluate the string \kbd{s} with GP and return the result as a string,524in a format similar to what GP displays (with the history index).525The resulting string is allocated on the PARI stack, so subsequent call526to \kbd{gp\_embedded} will destroy it.527528\section{Readline interface}529530Code which wants to use libpari readline (such as the Jupyter notebook)531needs to do the following:532\bprog533#include <readline.h>534#include <paripriv.h>535pari_rl_interface S;536...537pari_use_readline(S);538@eprog\noindent The variable $S$, as initialized above, encapsulates539the libpari readline interface. (And allow us to move gp's readline code540to libpari without introducing a mandatory dependency on readline in541libpari.) The following functions then become available:542543\fun{char**}{pari_completion_matches}{pari_rl_interface *pS, const char *s,544long pos, long *wordpos} given a command string $s$, where the cursor545is at index \kbd{pos}, return an array of completion matches.546547If \kbd{wordpos} is not \kbd{NULL}, set \kbd{*wordpos} to the index for the548start of the expression we complete.549550\fun{char**}{pari_completion}{pari_rl_interface *pS, char *text, int start,551int end} the low-level completer called by \tet{pari_completion_matches}.552The following wrapper553\bprog554char**555gp_completion(char *text, int START, int END)556{ return pari_completion(&S, text, START, END);)557@eprog\noindent is a valid value for558\tet{rl_attempted_completion_function}.559560\section{Constructors called by \kbd{pari\_init} functions}561562\fun{void}{pari_init_buffers}{}563564\fun{void}{pari_init_compiler}{}565566\fun{void}{pari_init_defaults}{}567568\fun{void}{pari_init_evaluator}{}569570\fun{void}{pari_init_files}{}571572\fun{void}{pari_init_floats}{}573574\fun{void}{pari_init_graphics}{}575576\fun{void}{pari_init_homedir}{}577578\fun{void}{pari_init_parser}{}579580\fun{void}{pari_init_paths}{}581582\fun{void}{pari_init_primetab}{}583584\fun{void}{pari_init_rand}{}585586\fun{void}{pari_init_seadata}{}587588\section{Destructors called by \kbd{pari\_close}}589590\fun{void}{pari_close_compiler}{}591592\fun{void}{pari_close_evaluator}{}593594\fun{void}{pari_close_files}{}595596\fun{void}{pari_close_floats}{}597598\fun{void}{pari_close_homedir}{}599600\fun{void}{pari_close_mf}{}601602\fun{void}{pari_close_parser}{}603604\fun{void}{pari_close_paths}{}605606\fun{void}{pari_close_primes}{}607608\section{Constructors and destructors used by the \kbd{pthreads} interface}609610\item Called by \tet{pari_thread_close}611612\fun{void}{pari_thread_close_files}{}613614\newpage615616\chapter{Regression tests, benches}617618This chapter documents how to write an automated test module, say \kbd{fun},619so that \kbd{make test-fun} executes the statements in the \kbd{fun} module620and times them, compares the output to a template, and prints an error621message if they do not match.622623\item Pick a \emph{new} name for your test, say \kbd{fun}, and write down a624GP script named \kbd{fun}. Make sure it produces some useful output and tests625adequately a set of routines.626627\item The script should not be too long: one minute runs should be enough.628Try to break your script into independent easily reproducible tests, this way629regressions are easier to debug; e.g. include \kbd{setrand(1)} statement before630a randomized computation. The expected output may be different on 32-bit and63164-bit machines but should otherwise be platform-independent. If possible, the632output shouldn't even depend on \kbd{sizeof(long)}; using a \kbd{realprecision}633that exists on both 32-bit and 64-bit architectures, e.g. \kbd{\bs p 38} is a634good first step. You can use \kbd{sizebyte(0)==16} to detect a 64-bit635architecture and \kbd{sizebyte(0)==8} for 32-bit.636637\item Dump your script into \kbd{src/test/in/} and run \kbd{Configure}.638639\item \kbd{make test-fun} now runs the new test, producing a \kbd{[BUG]} error640message and a \kbd{.dif} file in the relevant object directory \kbd{Oxxx}.641In fact, we compared the output to a nonexisting template, so this must fail.642643\item Now644\bprog645patch -p1 < Oxxx/fun.dif646@eprog\noindent647generates a template output in the right place \kbd{src/test/32/fun}, for648instance on a 32-bit machine.649650\item If different output is expected on 32-bit and 64-bit machines, run the651test on a 64-bit machine and patch again, thereby652producing \kbd{src/test/64/fun}. If, on the contrary, the output must be the653same (preferred behavior!), make sure the output template land in the654\kbd{src/test/32/} directory which provides a default template when the65564-bit output file is missing; in particular move the file from656\kbd{src/test/64/} to \kbd{src/test/32/} if the test was run on a 64-bit657machine.658659\item You can now re-run the test to check for regressions: no \kbd{[BUG]}660is expected this time! Of course you can at any time add some checks, and661iterate the test / patch phases. In particular, each time a bug in the662\kbd{fun} module is fixed, it is a good idea to add a minimal test case to663the test suite.664665\item By default, your new test is now included in \kbd{make test-all}. If666it is particularly annoying, e.g. opens tons of graphical windows as667\kbd{make test-ploth} or just much longer than the recommended minute, you668may edit \kbd{config/get\_tests} and add the \kbd{fun} test to the list of669excluded tests, in the \kbd{test\_extra\_out} variable.670671\item You can run a subset of existing tests by using the following idiom:672\bprog673cd Oxxx # call from relevant build directory674make TESTS="lfuntype lfun gamma" test-all675@eprog\noindent will run the \kbd{lfuntype}, \kbd{lfun} and \kbd{gamma} tests.676This produces a combined output whereas the alternative677\bprog678make test-lfuntype test-lfun test-gamma679@eprog\noindent would not.680681\item By default, the test is run on both the \kbd{gp-sta} and \kbd{gp-dyn}682binaries, making it twice as slow. If the test is somewhat long, it can683be annoying; you can restrict to one binary only using the \kbd{statest-all}684or \kbd{dyntest-all} targets. Both accept the \kbd{TESTS} argument seen above.685686\bprog687make test-lfuntype test-lfun gamma688@eprog\noindent would not.689690\item Finally, the \kbd{get\_tests} script also defines the recipe for691\kbd{make bench} timings, via the variable \kbd{test\_basic}. A test is692included as \kbd{fun} or \kbd{fun\_$n$}, where $n$ is an integer $\leq 1000$;693the latter means that the timing is weighted by a factor $n/1000$. (This was694introduced a long time ago, when the \kbd{nfields} bench was so much slower695than the others that it hid slowdowns elsewhere.)696697\section{Functions for GP2C}698699\subsec{Functions for safe access to components}700701These functions return the address of the requested component after checking702it is actually valid. This is used by GP2C -C.703704\fun{GEN*}{safegel}{GEN x, long l}, safe version of \kbd{gel(x,l)} for \typ{VEC},705\typ{COL} and \typ{MAT}.706707\fun{long*}{safeel}{GEN x, long l}, safe version of \kbd{x[l]} for \typ{VECSMALL}.708709\fun{GEN*}{safelistel}{GEN x, long l} safe access to \typ{LIST} component.710711\fun{GEN*}{safegcoeff}{GEN x, long a, long b} safe version of712\kbd{gcoeff(x,a, b)} for \typ{MAT}.713714\newpage715\chapter{Parallelism}716717PARI provides an abstraction, herafter called the MT engine, for doing718parallel computations. The exact same high level routines are used whether719the underlying communication protocol is POSIX threads or MPI and they behave720differently depending on how \kbd{libpari} was configured, specifically on721\kbd{Configure}'s \kbd{--mt} option. Sequential computation is also supported722(no \kbd{--mt} argument) which is helpful for debugging newly written723parallel code. The final section in this chapter comments a complete example.724725\section{The PARI multithread interface}726727\fun{void}{mt_queue_start}{struct pari_mt *pt, GEN worker} Let \kbd{worker}728be a \typ{CLOSURE} object of arity $1$. Initialize the opaque structure729\kbd{pt} to evaluate \kbd{worker} in parallel, using \kbd{nbthreads} threads.730This allocates data in731various ways, e.g., on the PARI stack or as malloc'ed objects: you may not732collect garbage on the PARI stack starting from an earlier \kbd{avma} point733until the parallel computation is over, it could destroy something in \kbd{pt}.734All ressources allocated outside the PARI stack are freed by735\kbd{mt\_queue\_end}.736737\fun{void}{mt_queue_start_lim}{struct pari_mt *pt, GEN worker, long lim}738as \kbd{mt\_queue\_start}, where \kbd{lim} is an upper bound on the number739of \kbd{tasks} to perform. Concretely the number of threads is the minimum740of \kbd{lim} and \kbd{nbthreads}. The values $0$ and $1$ of \kbd{lim} are741special:742743\item $0$: no limit, equivalent to \kbd{mt\_queue\_start} (use744\kbd{nbthreads} threads).745746\item $1$: no parallelism, evaluate the tasks sequentially.747748\fun{void}{mt_queue_submit}{struct pari_mt *pt, long taskid, GEN task} submit749\kbd{task} to be evaluated by \kbd{worker}; use \kbd{task = NULL} if no750further task needs to be submitted. The parameter \kbd{taskid} is attached to751the \kbd{task} but not used in any way by the \kbd{worker} or the MT engine,752it will be returned to you by \kbd{mt\_queue\_get} together with the result753for the task, allowing to match up results and submitted tasks if desired.754For instance, if the tasks $(t_1,\dots, t_m)$ are known in advance, stored in755a vector, and you want to recover the evaluation results in the same order as756in that vector, you may use consecutive integers $1, \dots, m$ as757\kbd{taskid}s. If you do not care about the ordering, on the other hand, you758can just use $\kbd{taskid} = 0$ for all tasks.759760The \kbd{taskid} parameter is ignored when \kbd{task} is \kbd{NULL}. It is761forbidden to call this function twice without an intervening762\kbd{mt\_queue\_get}.763764\fun{GEN}{mt_queue_get}{struct pari_mt *pt, long *taskid, long *pending}765return \kbd{NULL} until \kbd{mt\_queue\_submit} has submitted766tasks for the required number (\kbd{nbthreads}) of threads; then return the767result of the evaluation by \kbd{worker} of one of the previously submitted768tasks, in random order. Set \kbd{pending} to the number of remaining pending769tasks: if this is $0$ then no more tasks are pending and it is safe to call770\tet{mt_queue_end}. Set \kbd{*taskid} to the value attached to this task by771\kbd{mt\_queue\_submit}, unless the \kbd{taskid} pointer is \kbd{NULL}. It is772forbidden to call this function twice without an intervening773\kbd{mt\_queue\_submit}.774775\fun{void}{mt_queue_end}{struct pari_mt *pt} end the parallel execution776and free ressources attached to the opaque \kbd{pari\_mt} structure. For777instance malloc'ed data; in the \kbd{pthreads} interface, it would destroy778mutex locks, condition variables, etc. This must be called once there are no779longer pending tasks to avoid leaking ressources; but not before all tasks780have been processed else crashes will occur.781782\fun{long}{mt_nbthreads}{void} return the effective number of parallel threads783that would be started by \tet{mt_queue_start} if it has been called in place784of \tet{mt_nbthreads}.785786\section{Technical functions required by MPI}787788The functions in this section are needed when writing complex independent789programs in order to support the MPI MT engine, as more flexible790complement/variants of \kbd{pari\_init} and \kbd{pari\_close}.791792\fun{void}{mt_broadcast}{GEN code}: do nothing unless the MPI threading engine793is in use. In that case, evaluates the closure \kbd{code} on all secondary794nodes. This can be used to change the state of all MPI child nodes, e.g.,795in \tet{gpinstall} run in the main thread, which allows all nodes to use the796new function.797798\fun{void}{pari_mt_init}{void} \label{pari_mt_init}799when using MPI, it is often necessary to run initialization code on the child800nodes after PARI is initialized. This is done by calling successively:801802\item \tet{pari_init_opts} with the flag \tet{INIT_noIMTm}:803this initializes PARI, but not the MT engine;804805\item the required initialization code;806807\item \tet{pari_mt_init} to initialize the MT engine.808Note that under MPI, this function returns on the master node but enters809slave mode on the child nodes. Thus it is no longer possible to run810initialization code on the child nodes.811812\fun{void}{pari_mt_close}{void} \label{pari_mt_close}813when using MPI, calling \tet{pari_close} terminates the MPI execution814environment and it will not be possible to restart it. If this is815undesirable, call \tet{pari_close_opts} with the flag \tet{INIT_noIMTm}816instead of \kbd{pari\_close}: this closes PARI without terminating the MPI817execution environment. You may later call \kbd{pari\_mt\_close} to terminate818it. It is an error for a program to end without terminating the MPI execution819environment.820821\section{A complete example}822823We now proceed to an example exhibiting complex features of this824interface, in particular showing how to generate a valid \kbd{worker}.825Explanations and details follow.826827\bprogfile{../examples/pari-mt.c}828829We start from some arbitrary C function \kbd{Cworker} and create an830\kbd{entree} summarizing all that GP would need to know about it, in831particular832833\item a GP name \kbd{\_worker}; the leading \kbd{\_} is not necessary,834we use it as a namespace mechanism grouping private functions;835836\item the name of the C function;837838\item and its prototype, see \kbd{install} for an introduction to Prototype839Codes.840841\noindent The other three arguments ($0$, $20$ and \kbd{""}) are required in an842\kbd{entree} but not useful in our simple context: they are respectively a843valence ($0$ means ``nothing special''), a help section (20 is customary for844internal functions which need to be exported for technical reasons, see845\kbd{?20}), and a help text (no help).846847Then we initialize the MT engine; doing things in this order with a two part848initialization ensures that nodes have access to our \kbd{Cworker}. We849convert the \kbd{ep} data to a \typ{CLOSURE} using \kbd{strtofunction}, which850provides a valid \kbd{worker} to \kbd{mt\_queue\_start}. This creates a851parallel evaluation queue \kbd{mt}, and we proceed to submit all tasks,852recording all results. Results are stored in the right order853by making good use of the \kbd{taskid} label, although we have no control854over \emph{when} each result is returned. We finally free all ressources855attached to the \kbd{mt} structure. If needed, we could have collected all856garbage on the PARI stack using \kbd{gerepilecopy} on the \kbd{out} array and857gone on working instead of quitting.858859Note the argument passing convention for \kbd{Cworker}: the task consists of a860single vector containing all arguments as \kbd{GEN}s, which are interpreted861according to the function prototype, here \kbd{GL} so the first argument is862left as is and the second one is converted to a long integer. In more863complicated situations, this second (and possibly further) argument could864provide arbitrary evaluation contexts. In this example, we just used it as a865flag to indicate the kind of evaluation expected on the data: integer866factorization (0) or matrix determinant (1).867868Note also that869\bprog870gel(out, taskid) = mt_queue_get(&mt, &taskid, &pending);871@eprog \noindent instead of our use of a temporary \kbd{done} would have872undefined behaviour (\kbd{taskid} may be uninitialized in the left hand side).873\vfill\eject874\input index\end875876877