File: develop.tex

package info (click to toggle)
pari 2.17.3-1
  • links: PTS
  • area: main
  • in suites: forky, sid
  • size: 24,508 kB
  • sloc: ansic: 281,184; sh: 861; perl: 420; yacc: 214; makefile: 162; f90: 88
file content (1117 lines) | stat: -rw-r--r-- 42,523 bytes parent folder | download | duplicates (3)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
\def\TITLE{Developer's Guide to the PARI library}
\input parimacro.tex

% START TYPESET
\begintitle
\vskip 2.5truecm
\centerline{\mine Developer's Guide}
\vskip 1.truecm
\centerline{\mine to}
\vskip 1.truecm
\centerline{\mine the PARI library}
\vskip 1.truecm
\centerline{\sectiontitlebf (version \vers)}
\vskip 1.truecm
\authors
\endtitle

\copyrightpage
\tableofcontents
\openin\std=develop.aux
\ifeof\std
\else
  \input develop.aux
\fi
\chapno=0

\chapter{Work in progress}

This draft documents private internal functions and structures for hard-core
PARI developers. Anything in here is liable to change on short notice. Don't
use anything in the present document, unless you are implementing new
features for the PARI library. Try to fix the interfaces before using them,
or document them in a better way.
If you find an undocumented hack somewhere, add it here.

Hopefully, this will eventually document everything that we buried in
\kbd{paripriv.h} or even more private header files like \kbd{anal.h}.
Possibly, even implementation choices! Way to go.

\section{The type \typ{CLOSURE}}\kbdsidx{t_CLOSURE}\sidx{closure}
This type holds closures and functions in compiled form, so is deeply
linked to the internals of the GP compiler and evaluator.
The length of this type can be $6$, $7$ or $8$ depending whether the
object is an ``inline closure'', a ``function'' or a ``true closure''.

A function is a regular GP function. The GP input line is treated as a
function of arity $0$.

A true closure is a GP function defined in a nonempty lexical context.

An inline closure is a closure that appears in the code without
the preceding \kbd{->} token. They are generally attached to the prototype
code 'E' and 'I'. Inline closures can only exist as data of other closures,
see below.

In the following example,
\bprog
f(a=Euler)=x->sin(x+a);
g=f(Pi/2);
plot(x=0,2*Pi,g(x))
@eprog\noindent
\kbd{f} is a function, \kbd{g} is a true closure and both \kbd{Euler} and
\kbd{g(x)} are inline closures.

This type has a second codeword \kbd{z[1]}, which is the arity of the
function or closure. This is zero for inline closures. To access it, use

\fun{long}{closure_arity}{GEN C}

\item \kbd{z[2]} points to a \typ{STR} which holds the opcodes. To access it, use

\fun{GEN}{closure_get_code}{GEN C}.

\fun{const char *}{closure_codestr}{GEN C} returns as an array of \kbd{char}
starting at $1$.

\item \kbd{z[3]} points to a \typ{VECSMALL} which holds the operands of the opcodes.
To access it, use

\fun{GEN}{closure_get_oper}{GEN C}

\item \kbd{z[4]} points to a \typ{VEC} which hold the data referenced by the
\kbd{pushgen} opcodes, which can be \typ{CLOSURE}, and in particular
inline closures. To access it, use

\fun{GEN}{closure_get_data}{GEN C}

\item \kbd{z[5]} points to a \typ{VEC} which hold extra data needed for
error-reporting and debugging. See \secref{se:dbgclosure} for details.
To access it, use

\fun{GEN}{closure_get_dbg}{GEN C}

Additionally, for functions and true closures,

\item \kbd{z[6]} usually points to a \typ{VEC} with two components which are \typ{STR}.
The first one displays the list of arguments of the closure without the
enclosing parentheses, the second one the GP code of the function at the
right of the \kbd{->} token. They are used to display the closure, either in
implicit or explicit form. However for closures that were not generated from GP
code, \kbd{z[6]} can point to a \typ{STR} instead. To access it, use

\fun{GEN}{closure_get_text}{GEN C}

Additionally, for true closure,

\item \kbd{z[7]} points to a \typ{VEC} which holds the values of all lexical
variables defined in the scope the closure was defined. To access it, use

\fun{GEN}{closure_get_frame}{GEN C}

\subsec{Debugging information in closure}\label{se:dbgclosure}

Every \typ{CLOSURE} object \kbd{z} has a component \kbd{dbg=z[5]}
which hold extra data needed for error-reporting and debugging.
The object \kbd{dbg} is a \typ{VEC} with $3$ components:

\kbd{dbg[1]} is a \typ{VECSMALL} of the same length than \kbd{z[3]}. For each
opcode, it holds the position of the corresponding GP source code in the
strings stored in \kbd{z[6]} for function or true closures, positive indices
referring to the second strings, and negative indices referring to the first
strings, the last element being indexed as $-1$. For inline closures, the
string of the parent function or true closure is used instead.

\kbd{dbg[2]} is a \typ{VECSMALL} that lists opcodes index where new lexical
local variables are created. The value $0$ denotes the position before the
first offset and variables created by the prototype code 'V'.

\kbd{dbg[3]} is a \typ{VEC} of \typ{VECSMALL}s that give the list of
\kbd{entree*} of the lexical local variables created at a given index in
\kbd{dbg[2]}.

\section{The type \typ{LIST}}\kbdsidx{t_LIST}\sidx{list} This type needs to go
through various hoops to support GP's inconvenient memory model. Don't
use \typ{LIST}s in pure library mode, reimplement ordinary lists! This
dynamic type is implemented by a \kbd{GEN} of length 3: two codewords and a
vector containing the actual entries. In a normal setup (a finished list,
ready to be used),

\item the vector is malloc'ed, so that it can be realloc'ated without moving
the parent \kbd{GEN}.

\item all the entries are clones, possibly with cloned subcomponents; they
must be deleted with \tet{gunclone_deep}, not \tet{gunclone}.

The following macros are proper lvalues and access the components

\fun{long}{list_nmax}{GEN L}: current maximal number of elements. This grows
as needed.

\fun{GEN}{list_data}{GEN L}: the elements. If \kbd{v = list\_data(L)}, then
either \kbd{v} is \kbd{NULL} (empty list) or \kbd{l = lg(v)} is defined, and
the elements are \kbd{v[1]}, \dots, \kbd{v[l-1]}.

In most \kbd{gerepile} scenarios, the list components are not inspected
and a shallow copy of the malloc'ed vector is made. The functions
\kbd{gclone}, \kbd{copy\_bin\_canon} are exceptions, and make a full copy of
the list.

The main problem with lists is to avoid memory leaks; in the above setup,
a statement like \kbd{a = List(1)} would already leak memory, since
\kbd{List(1)} allocates memory, which is cloned (second allocation) when
assigned to \kbd{a}; and the original list is lost. The solution we
implemented is

\item to create anonymous lists (from \kbd{List}, \kbd{gtolist},
\kbd{concat} or \kbd{vecsort}) entirely on the stack, \emph{not} as described
above, and to set \kbd{list\_nmax} to $0$. Such a list is not yet proper and
trying to append elements to it fails:
\bprog
? listput(List(),1)
  ***   variable name expected: listput(List(),1)
  ***                                   ^----------------
@eprog\noindent
If we had been malloc'ing memory for the
\kbd{List([1,2,3])}, it would have leaked already.

\item as soon as a list is assigned to a variable (or a component thereof)
by the GP evaluator, the assigned list is converted to the proper format
(with \kbd{list\_nmax} set) previously described.

\fun{GEN}{listcopy}{GEN L} return a full copy of the \typ{LIST}~\kbd{L},
allocated on the \emph{stack} (hence \kbd{list\_nmax} is $0$). Shortcut for
\kbd{gcopy}.

\fun{GEN}{mklistcopy}{GEN x} returns a list with a single element $x$,
allocated on the stack. Used to implement most cases of \kbd{gtolist}
(except vectors and lists).

A typical low-level construct:
\bprog
  long l;
  /* assume L is a t_LIST */
  L = list_data(L); /* discard t_LIST wrapper */
  l = L? lg(L): 1;
  for (i = 1; i < l; i++) output( gel(L, i) );
  for (i = 1; i < l; i++) gel(L, i) = gclone( ... );
@eprog

\subsec{Maps as Lists}

GP's maps are implemented on top of \typ{LIST}s so as to benefit from
their peculiar memory models. Lists thus come in two subtypes: \typ{LIST\_RAW}
(actual lists) and \typ{LIST\_MAP} (a map).

\fun{GEN}{mklist_typ}{long t} create a list of subtype $t$.
\fun{GEN}{mklist}{void} is an alias for
\bprog
  mklist_typ(t_LIST_RAW);
@eprog
and
\fun{GEN}{mkmap}{void} is an alias for
\bprog
  mklist_typ(t_LIST_MAP);
@eprog

\fun{long}{list_typ}{GEN L} return the list subtype, either \typ{LIST\_RAW} or
\typ{LIST\_MAP}.

\fun{void}{listpop}{GEN L, long index} as \kbd{listpop0},
assuming that $L$ is a \typ{LIST\_RAW}.

\fun{GEN}{listput}{GEN list, GEN object, long index} as \kbd{listput0},
assuming that $L$ is a \typ{LIST\_RAW}. Return the element as inserted
in the list (a clone of \kbd{object}).

\fun{GEN}{listinsert}{GEN list, GEN object, long index} as \kbd{listinsert0},
assuming that $L$ is a \typ{LIST\_RAW}. Return the element as inserted
in the list (a clone of \kbd{object}).

\fun{GEN}{mapdomain}{GEN T} vector of keys of the map $T$.

\fun{GEN}{mapselect_shallow}{void *E, long (*f)(void* E, GEN x), GEN T}
vector of keys of the map $T$ whose value matches the predicate $f$.

\fun{GEN}{mapdomain_shallow}{GEN T} shallow version of \kbd{mapdomain}.

\fun{GEN}{maptomat}{GEN T} convert a map to a factorization matrix.

\fun{GEN}{maptomat_shallow}{GEN T} shallow version of \kbd{maptomat}.

\section{Protection of noninterruptible code}

GP allows the user to interrupt a computation by issuing SIGINT
(usually by entering control-C) or SIGALRM (usually using alarm()).
To avoid such interruption to occurs in section of code which are not
reentrant (in particular \kbd{malloc} and \kbd{free})
the following mechanism is provided:

\fun{}{BLOCK_SIGINT_START}{}
  Start a noninterruptible block code. Block both \kbd{SIGINT} and \kbd{SIGARLM}.

\fun{}{BLOCK_SIGALRM_START}{}
  Start a noninterruptible block code. Block only \kbd{SIGARLM}.
This is used in the \kbd{SIGINT} handler itself to delay an eventual pending
alarm.

\fun{}{BLOCK_SIGINT_END}{}
  End a noninterruptible block code

The above macros make use of the following global variables:

\tet{PARI_SIGINT_block}: set to $1$ (resp. $2$) by \kbd{BLOCK\_SIGINT\_START}
(resp. \kbd{BLOCK\_SIGALRM\_START}).

\tet{PARI_SIGINT_pending}: Either $0$ (no signal was blocked), \kbd{SIGINT}
(\kbd{SIGINT} was blocked) or \kbd{SIGALRM} (\kbd{SIGALRM} was blocked).
This need to be set by the signal handler.

Within a block, an automatic variable \kbd{int block} is defined which
records the value of \kbd{PARI\_SIGINT\_block} when entering the block.

\subsec{Multithread interruptions}

To support multithreaded programs, \kbd{BLOCK\_SIGINT\_START} and
\kbd{BLOCK\_SIGALRM\_START} call \kbd{MT\_SIGINT\_BLOCK(block)}, and
\kbd{BLOCK\_SIGINT\_END} calls \kbd{MT\_SIGINT\_UNBLOCK(block)}.

\tet{MT_SIGINT_BLOCK} and \tet{MT_SIGINT_UNBLOCK} are defined by the
multithread engine. They can calls the following public functions defined by
the multithread engine.

\fun{void}{mt_sigint_block}{void}

\fun{void}{mt_sigint_unblock}{void}

In practice this mechanism is used by the POSIX thread engine to protect against
asychronous cancellation.

\section{$\F_{l^2}$ field for small primes $l$}
Let $l>2$ be a prime \kbd{ulong}.  A \kbd{Fl2} is an element of the finite
field $\F_{l^2}$ represented (currently) by a \kbd{Flx} of degree at most $1$
modulo a polynomial of the form $x^2-D$ for some non square $0\leq D<p$.
Below \kbd{pi} denotes the pseudo inverse of \kbd{p}, see \kbd{Fl\_mul\_pre}

\fun{int}{Fl2_equal1}{GEN x} return $1$ if $x=1$, else return $0$.

\fun{GEN}{Fl2_mul_pre}{GEN x, GEN y, ulong D, ulong p, ulong pi} return $x\*y$.

\fun{GEN}{Fl2_sqr_pre}{GEN x, ulong D, ulong p, ulong pi} return $x^2$.

\fun{GEN}{Fl2_inv_pre}{GEN x, ulong D, ulong p, ulong pi} return $x^{-1}$.

\fun{GEN}{Fl2_pow_pre}{GEN x, GEN n, ulong D, ulong p, ulong pi} return
$x^n$.

\fun{GEN}{Fl2_sqrt_pre}{GEN a, ulong D, ulong p, ulong pi}, as \kbd{Flxq\_sqrt}.

\fun{GEN}{Fl2_sqrtn_pre}{GEN a, GEN n, ulong D, ulong p, ulong pi, GEN *zeta}
$n$-th root, as \kbd{Flxq\_sqrtn}.

\fun{GEN}{Fl2_norm_pre}{GEN x, GEN n, ulong D, ulong p, ulong pi} return the
norm of $x$.

\fun{GEN}{Flx_Fl2_eval_pre}{GEN P, GEN x, ulong D, ulong p, ulong pi}
return $P(x)$.

\section{Public functions useless outside of GP context}

These functions implement GP functionality for which the C language or
other libpari routines provide a better equivalent; or which are so tied
to the \kbd{gp} interpreter as to be virtually useless in \kbd{libpari}. Some
may be generated by \kbd{gp2c}. We document them here for completeness.

\subsec{Conversions}

\fun{GEN}{toser_i}{GEN x} internal shallow function, used to implement
automatic conversions to power series in GP (as in \kbd{cos(x)}).
Converts a \typ{POL} or a \typ{RFRAC} to a \typ{SER} in the same variable and
precision \kbd{precdl} (the global variable corresponding to
\kbd{seriesprecision}). Returns $x$ itself for a \typ{SER}, and \kbd{NULL}
for other argument types. The fact that it uses a global variable makes it
awkward whenever you're not implementing a new transcendental function in GP.
Use \tet{RgX_to_ser} or \tet{rfrac_to_ser} for a fast clean alternative to
\kbd{gtoser}.

\fun{GEN}{listinit}{GEN x} a \typ{LIST} (from \kbd{List} or \kbd{Map}) may
exist in two different forms due to GP memory model:

\item an ordinary \emph{read-only} copy on the PARI stack (as produced by
\kbd{gtolist} or \kbd{gtomap}) to which one may not assign elements
(\kbd{listput} will fail) unless the list is empty.

\item a feature-complete GP list using (malloc'ed) \kbd{block}s to
allow dynamic insertions. An empty list is automaticaly promoted to this
status on first insertion.

The \kbd{listinit} function creates a copy of existing \typ{SER} $x$ and
makes sure it is of the second kind. Variants of this are automatically
called by \kbd{gp} when assigning a \typ{LIST} to a GP variable; the
mecanism avoid memory leaks when creating a constant list, e.g.
\kbd{List([1,2,3])} (read-only), without assigning it to a variable. Whereas
after \kbd{L = List([1,2,3])} (GP list), we keep a pointer to the object and
may delete it when $L$ goes out of scope.

This \kbd{libpari} function allows \kbd{gp2c} to simulate this process by
generating \kbd{listinit} calls at appropriate places.

\subsec{Output}

\fun{void}{out_print1}{PariOUT *out, const char *sep, GEN g, long flag}
internal function underlying the
\kbd{print} GP function. Prints the entries of the \typ{VEC} $g$, one by one,
without any separator; entries of type \typ{STR} are printed without enclosing
quotes. \fl is one of \tet{f_RAW}, \tet{f_PRETTYMAT} or \tet{f_TEX}, using the
output context \kbd{out} and separator \kbd{sep} between
successive entries (no separator if \kbd{NULL}).

\fun{void}{out_print0}{PariOUT *out, const char *sep, GEN g, long flag} as
\tet{out_print1}, but also output a terminating newline.

\fun{char*}{pari_sprint0}{const char *s, GEN g, long flag} displays $s$,
then \kbd{print0(g, flag)}.

\subsec{Input}

\kbd{gp}'s input is read from the stream \tet{pari_infile}, which is changed
using

\fun{FILE*}{switchin}{const char *name}

Note that this function is quite complicated, maintaining stacks of files
to allow smooth error recovery and \kbd{gp} interaction. You will be better
off using \tet{gp_read_file}.

\subsec{Control flow statements}

\fun{GEN}{break0}{long n}. Use the C control statement \kbd{break}. Since
\kbd{break(2)} is invalid in C, either rework your code or use \kbd{goto}.

\fun{GEN}{next0}{long n}. Use the C control statement \kbd{continue}. Since
\kbd{continue(2)} is invalid in C, either rework your code or use \kbd{goto}.

\fun{GEN}{return0}{GEN x}. Use \kbd{return}!

\subsec{Accessors}

\fun{GEN}{vecslice0}{GEN A, long a, long b} implements $A[a..b]$.

\fun{GEN}{matslice0}{GEN A, long a, long b, long c, long d}
implements $A[a..b, c..d]$.

\subsec{Iterators}

\fun{GEN}{apply0}{GEN f, GEN A} gp wrapper calling \tet{genapply}, where $f$
is a \typ{CLOSURE}, applied to $A$. Use \kbd{genapply} or a standard C loop.

\fun{GEN}{select0}{GEN f, GEN A} gp wrapper calling \tet{genselect}, where $f$
is a \typ{CLOSURE} selecting from $A$. Use \kbd{genselect} or a standard C loop.

\fun{GEN}{vecapply}{void *E, GEN (*f)(void* E, GEN x), GEN x} implements
\kbd{[a(x)|x<-b]}.

\fun{GEN}{veccatapply}{void *E, GEN (*f)(void* E, GEN x), GEN x} implements
\kbd{concat([a(x)|x<-b])} which used to implement \kbd{[a0(x,y)|x<-b;y<-c(b)]}
which is equal to \kbd{concat([[a0(x,y)|y<-c(b)]|x<-b])}.

\fun{GEN}{vecselect}{void *E, long (*f)(void* E, GEN x), GEN A}
implements \kbd{[x<-b,c(x)]}.

\fun{GEN}{vecselapply}{void *Epred, long (*pred)(void* E, GEN x), void *Efun, GEN (*fun)(void* E, GEN x), GEN A}
implements \kbd{[a(x)|x<-b,c(x)]}.

\subsec{Local precision}

These functions allow to change \kbd{realprecision} locally when
calling the GP interpretor.

\fun{void}{push_localprec}{long p} set the local precision to $p$.

\fun{void}{push_localbitprec}{long b} set the local precision to $b$ bits.

\fun{void}{pop_localprec}{void} reset the local precision to the previous
value.

\fun{long}{get_localprec}{void} returns the current local precision.

\fun{long}{get_localbitprec}{void} returns the current local precision in bits.

\fun{void}{localprec}{long p} trivial wrapper around \kbd{push\_localprec}
(sanity checks and convert from decimal digits to a number of codewords).
Use \kbd{push\_localprec}.

\fun{void}{localbitprec}{long p}
 trivial wrapper around \kbd{push\_localbitprec}
(sanity checks). Use \kbd{push\_localbitprec}.

These two function are used to implement \kbd{getlocalprec} and
\kbd{getlocalbitprec} for the GP interpreter and essentially return their
argument (the current dynamic precision, respectively in bits or as a
\kbd{prec} word count):

\fun{long}{getlocalbitprec}{long bit}

\fun{long}{getlocalprec}{long prec}


\subsec{Functions related to the GP evaluator}

The prototype code \kbd{C} instructs the GP compiler to save the current
lexical context (pairs made of a lexical variable name and its value)
in a \kbd{GEN}, called \kbd{pack} in the sequel. This \kbd{pack} can be used
to evaluate expressions in the corresponding lexical context, providing it is
current.

\fun{GEN}{localvars_read_str}{const char *s, GEN pack} evaluate the string $s$
in the lexical context given by \kbd{pack}.  Used by \tet{geval_gp} in GP
to implement the behavior below:
\bprog
? my(z=3);eval("z=z^2");z
%1 = 9
@eprog

\fun{long}{localvars_find}{GEN pack, entree *ep} does \kbd{pack} contain
a pair whose variable corresponds to \kbd{ep}? If so, where is the
corresponding value? (returns an offset on the value stack).

\subsec{Miscellaneous}

\fun{char*}{os_getenv}{const char *s} either calls \kbd{getenv}, or directly
return \kbd{NULL} if the \kbd{libc} does not provide it. Use \tet{getenv}.

\fun{sighandler_t}{os_signal}{int sig, pari_sighandler_t fun} after a
\bprog
  typedef void (*pari_sighandler_t)(int);
@eprog\noindent
(private type, not exported). Installs signal handler \kbd{fun} for
signal \kbd{sig}, using \tet{sigaction} with flag \tet{SA_NODEFER}. If
\kbd{sigaction} is not available use \tet{signal}. If even the latter is not
available, just return \tet{SIG_IGN}. Use \tet{sigaction}.

\section{Embedded GP interpretor}
These function provide a simplified interface to embed a GP
interpretor in a program.

\fun{void}{gp_embedded_init}{long rsize, long vsize}
Initialize the GP interpretor (like \kbd{pari\_init} does) with
\kbd{parisize=rsize} and \kbd{parisizemax=vsize}.

\fun{long}{gp_embedded}{const char *s}
Evaluate the string \kbd{s} as if it was input to the GP interpretor,
displaying results, timings and error, filling the history etc, as
needed. return $1$ if an error was triggered, $0$ otherwise.

\section{Readline interface}

Code which wants to use libpari readline (such as the Jupyter notebook)
needs to do the following:
\bprog
#include <readline.h>
#include <paripriv.h>
pari_rl_interface S;
...
pari_use_readline(S);
@eprog\noindent The variable $S$, as initialized above, encapsulates
the libpari readline interface. (And allow us to move gp's readline code
to libpari without introducing a mandatory dependency on readline in
libpari.) The following functions then become available:

\fun{char**}{pari_completion_matches}{pari_rl_interface *pS, const char *s,
long pos, long *wordpos} given a command string $s$, where the cursor
is at index \kbd{pos}, return an array of completion matches.

If \kbd{wordpos} is not \kbd{NULL}, set \kbd{*wordpos} to the index for the
start of the expression we complete.

\fun{char**}{pari_completion}{pari_rl_interface *pS, char *text, int start,
int end} the low-level completer called by \tet{pari_completion_matches}.
The following wrapper
\bprog
char**
gp_completion(char *text, int START, int END)
{ return pari_completion(&S, text, START, END);)
@eprog\noindent is a valid value for
\tet{rl_attempted_completion_function}.

\section{Constructors called by \kbd{pari\_init} functions}

\fun{void}{pari_init_buffers}{}

\fun{void}{pari_init_compiler}{}

\fun{void}{pari_init_defaults}{}

\fun{void}{pari_init_evaluator}{}

\fun{void}{pari_init_files}{}

\fun{void}{pari_init_floats}{}

\fun{void}{pari_init_graphics}{}

\fun{void}{pari_init_homedir}{}

\fun{void}{pari_init_parser}{}

\fun{void}{pari_init_paths}{}

\fun{void}{pari_init_primetab}{}

\fun{void}{pari_init_rand}{}

\section{Destructors called by \kbd{pari\_close}}

\fun{void}{pari_close_compiler}{}

\fun{void}{pari_close_evaluator}{}

\fun{void}{pari_close_files}{}

\fun{void}{pari_close_floats}{}

\fun{void}{pari_close_homedir}{}

\fun{void}{pari_close_mf}{}

\fun{void}{pari_close_parser}{}

\fun{void}{pari_close_paths}{}

\fun{void}{pari_close_primes}{}

\section{Constructors and destructors used by the \kbd{pthreads} interface}

\item Called by \tet{pari_thread_close}

\fun{void}{pari_thread_close_files}{}

\section{Functions for GP2C}

\subsec{Functions for safe access to components}

These functions return the address of the requested component after checking
it is actually valid. This is used by GP2C -C.

\fun{GEN*}{safegel}{GEN x, long l}, safe version of \kbd{gel(x,l)} for \typ{VEC},
\typ{COL} and \typ{MAT}.

\fun{long*}{safeel}{GEN x, long l}, safe version of \kbd{x[l]} for \typ{VECSMALL}.

\fun{GEN*}{safelistel}{GEN x, long l} safe access to \typ{LIST} component.

\fun{GEN*}{safegcoeff}{GEN x, long a, long b} safe version of
\kbd{gcoeff(x,a, b)} for \typ{MAT}.

\newpage

\chapter{Regression tests, benches}

This chapter documents how to write an automated test module, say \kbd{fun},
so that \kbd{make test-fun} executes the statements in the \kbd{fun} module
and times them, compares the output to a template, and prints an error
message if they do not match.

\item Pick a \emph{new} name for your test, say \kbd{fun}, and write down a
GP script named \kbd{fun}. Make sure it produces some useful output and tests
adequately a set of routines.

\item The script should not be too long: one minute runs should be enough.
Try to break your script into independent easily reproducible tests, this way
regressions are easier to debug; e.g. include \kbd{setrand(1)} statement before
a randomized computation. The expected output may be different on 32-bit and
64-bit machines but should otherwise be platform-independent. If possible, the
output shouldn't even depend on \kbd{sizeof(long)}; using a \kbd{realprecision}
that exists on both 32-bit and 64-bit architectures, e.g. the
default~\kbd{\bs p 38} is a good first step. You can use
\kbd{sizebyte(0)==16} to detect a 64-bit architecture and
\kbd{sizebyte(0)==8} for 32-bit.

\item Dump your script into \kbd{src/test/in/} and run \kbd{Configure}.

\item \kbd{make test-fun} now runs the new test, producing a \kbd{[BUG]} error
message and a \kbd{.dif} file in the relevant object directory \kbd{Oxxx}.
In fact, we compared the output to a nonexisting template, so this must fail.

\item Now
\bprog
  patch -p1 < Oxxx/fun.dif
@eprog\noindent
generates a template output in the right place \kbd{src/test/32/fun}, for
instance on a 32-bit machine.

\item If different output is expected on 32-bit and 64-bit machines, run the
test on a 64-bit machine and patch again, thereby
producing \kbd{src/test/64/fun}. If, on the contrary, the output must be the
same (preferred behavior!), make sure the output template land in the
\kbd{src/test/32/} directory which provides a default template when the
64-bit output file is missing; in particular move the file from
\kbd{src/test/64/} to \kbd{src/test/32/} if the test was run on a 64-bit
machine.

\item You can now re-run the test to check for regressions: no \kbd{[BUG]}
is expected this time! Of course you can at any time add some checks, and
iterate the test / patch phases. In particular, each time a bug in the
\kbd{fun} module is fixed, it is a good idea to add a minimal test case to
the test suite.

\item By default, your new test is now included in \kbd{make test-all}. If
it is particularly annoying, e.g. opens tons of graphical windows as
\kbd{make test-ploth} or just much longer than the recommended minute, you
may edit \kbd{config/get\_tests} and add the \kbd{fun} test to the list of
excluded tests, in the \kbd{test\_extra\_out} variable.

\item You can run a subset of existing tests by using the following idiom:
\bprog
  cd Oxxx     # call from relevant build directory
  make TESTS="lfuntype lfun gamma" test-all
@eprog\noindent will run the \kbd{lfuntype}, \kbd{lfun} and \kbd{gamma} tests.
This produces a combined output whereas the alternative
\bprog
  make test-lfuntype test-lfun test-gamma
@eprog\noindent would not.

\item By default, the test is run on both the \kbd{gp-sta} and \kbd{gp-dyn}
binaries, making it twice as slow. If the test is somewhat long, it can
be annoying; you can restrict to one binary only using the \kbd{statest-all}
or \kbd{dyntest-all} targets. Both accept the \kbd{TESTS} argument seen
above; again, the alternative

\bprog
  make test-lfuntype test-lfun test-gamma
@eprog\noindent would not.

\item Finally, the \kbd{get\_tests} script also defines the recipe for
\kbd{make bench} timings, via the variable \kbd{test\_basic}. A test is
included as \kbd{fun} or \kbd{fun\_$n$}, where $n$ is an integer $\leq 1000$;
the latter means that the timing is weighted by a factor $n/1000$. (This was
introduced a long time ago, when the \kbd{nfields} bench was so much slower
than the others that it hid slowdowns elsewhere.)

\newpage
\chapter{Tips and tricks to edit and debug the PARI sources}

First, fetch the pari sources from our \kbd{git} repository:
\bprog
  git clone https://pari.math.u-bordeaux.fr/git/pari.git
@eprog\noindent
The following sections assume you defined an environment variable
\kbd{PARIDIR}, pointing to the toplevel where you install the pari sources,
e.g. in \kbd{.bashrc},
\bprog
  export PARIDIR=$HOME/pari
@eprog

\section{Tags}
To allow seamless navigation in PARI sources, we will use a
\kbd{ctags} program. For instance, on Debian/Ubuntu systems,
\bprog
  sudo apt-get install exuberant-ctags
@eprog

\subsec{For \kbd{vim}}
Type \kbd{make ctags} in PARI's toplevel, then add to your \kbd{.vimrc}
\bprog
  set tags=./tags,$PARIDIR/src/tags
@eprog\noindent Now you can navigate the PARI sources in the usual \kbd{vim}
way. Try \kbd{vi -t gadd}, then move the cursor to an identifier in the
neighborhood, \kbd{typ} say, then type \kbd{Ctrl-]}  (and \kbd{Ctrl-t} to come
back).

\subsec{For \kbd{emacs}}
Type \kbd{make etags} in PARI's toplevel, then add to your \kbd{.emacs}
\bprog
  (setq tags-table-list '("$PARIDIR/src"))
@eprog

\section{Editor customization}

You can define special-purpose editor options to help you edit PARI code. We
restrict to the \kbd{vim} editor; \kbd{emacs} and other programmer's editors
give you analogous possibilities.

\subsec{Create and customize \kbd{\$HOME/.vim/ftplugin/c.vim}}

There you change \kbd{vim} settings whenever you're editing a C file.
(Cleaner would be to call this \kbd{pari.vim} and load it conditionally from
\kbd{c.vim} but this goes beyond the scope of these notes.) For instance

\bprog
  setlocal path+=$PARIDIR/src/headers
  setlocal path+=$PARIDIR/Olinux-x86_64
@eprog\noindent This allows the editor to find and process our include files;
for instance typing \kbd{gf} on \kbd{\#include "paripriv.h"} will open that
file!

By the way, you can define more complex macros there, such as
\bprog
map <buffer> <Esc>r i  return NULL; /*LCOV_EXCL_LINE*/<Esc>
@eprog\noindent When typing \kbd{<M-i>}, this inserts at the cursor the string
\bprog
  return NULL; /*LCOV_EXCL_LINE*/
@eprog\noindent which is used in \kbd{libpari} code to mark unreachable
lines, inserted to avoid warnings from the compiler, for instance after an
error which we know will not return, when the function has a non \kbd{void}
return type. The \kbd{LCOV\_EXCL\_LINE} comment instructs our coverage tool
not to report that line as not reached by our testing suite. See
\kbd{https://pari.math.u-bordeaux.fr/lcov-report/}.

\bprog
map <buffer> <Esc>a :if expand("%:t") == 'paridecl.h'\
  <Bar> edit #\
  <Bar> else\
  <Bar> edit $PARIDIR/src/headers/paridecl.h\
  <Bar> endif<C-M>
@eprog\noindent This opens \kbd{paridecl.h} when hitting \kbd{<M-a>}
from a PARI source file. Typing \kbd{<M-a>} once again moves back to where
you where. That allows to edit a function, declare it or fix its declaration,
and come back working. The same can be done for \kbd{paripriv.h}, obviously.

\subsec{Create and customize \kbd{\$HOME/.vim/after/syntax/c.vim}}

This allows to colorize in an appropriate way PARI-specific types and
constants. For instance
\bprog
  syntax keyword cType GEN pari_sp ulong uchar hashentry hashtable
  syntax keyword cNumber avma
  syntax keyword cNumber gen_0 gen_1 gen_m1 gen_2 gen_m2 ghalf
@eprog\noindent
We also like to flag spaces at the end of a line as an error:
\bprog
  syntax match cError " \+$"
@eprog

\section{GDB configuration}
Edit \kbd{\$HOME/.gdbinit} to teach \kbd{gdb} a number of useful macros.
Here are some of the ones we use. Elaborate according to your needs.

\subsec{Printing GENs}

\bprog
define o
  call output((GEN)$arg0)
end
document o
  Print the GEN value the same way as gp's print().
  For instance o x prints the GEN x, o gadd(gen_1, gen_2) prints '3', etc.
end
@eprog

\bprog
define om
  call outmat((GEN)$arg0)
end
document om
  Pretty-print the GEN value, very suitable for matrices. For instance
  'o matid(2)' prints [1,0;0,1] but 'om matid(2)' prints
    [1 0]

    [0 1]
end
@eprog

\bprog
define v
  if $argc > 1
    call dbgGEN((GEN)$arg0, $arg1)
  else
    call dbgGEN((GEN)$arg0, 2)
  end
end
document v
  Emit the GEN value the same way as gp's \x would. By default,
  truncates the leaves at 2 words. This is changed by the 2nd optional
  argument; special case -1 = no truncation.
end
@eprog

The following GDB macros apply some preprocessing before printing so write and
leave intermediate results (e.g., \kbd{liftall} result) on the PARI stack.
\bprog
define olm
  call outmat(liftall((GEN)$arg0))
end
document olm
  Pretty-print the GEN value lifting whenever possible.
end
@eprog

\bprog
define osm # prec_w: shallow copy with precision decreased
  call outmat(gprec_w((GEN)$arg0,3))
end
document osm
  Pretty-print the GEN value, printing floats with minimal precision,
  suitable for large matrices or polynomials with floating point
  entries of huge accuracy.
end
@eprog

\subsec{Advanced use}

\bprog
define be
  break pari_err
end
document be
  Set a breakpoint in the PARI error handler.
end

define db
  call setalldebug($argv0)
end
document db
  Set 'DEBUGLEVEL' variables in all modules to the value given as argument
  (between 0 and 10). Emulates gp's \g n.
end

define fs
  call dbg_fill_stack()
end
document fs
  Overwrite the unused portion of PARI's stack by junk. Allows to diagnose
  garbage collection mistakes.
end
@eprog

The final one is particularly useful.
\bprog
define w1
  shell rm -f /var/tmp/gp.tmp1
  call gpwritebin("/var/tmp/gp.tmp1",$arg0)
end
document w1
  Write GEN argument to temp file in binary PARI format
end
@eprog\noindent The following GP function can now work its magic:
\bprog
  rtmp({n = 1}) = read( Str("/var/tmp/gp.tmp", n) );
@eprog\noindent In particular, in any GP shell, \kbd{rtmp()} will read in
the argument that was saved under \kbd{gdb} by \kbd{w1}, allowing to test and
inspect it with GP functions and start again on error. This is in sharp
contrast to the same inspection under GDB, which is error prone and
unforgiving: an error when manipulation C functions is likely to crash the
process, hitting a garbage collection point might destroy the object, etc.
This is particularly useful when the saved value takes a non-negligible time
to compute.

Another possibility is to compare the behavior of a newly modified \kbd{gp}
with a reference version. One can find and save internal values from the two
different \kbd{gp} binaries from \kbd{gdb} and compare or test them in either
one.

One can obviously duplicate the definition to have \kbd{w2} write to
\kbd{gp.tmp2}, etc. Here is a more generic solution:
\bprog
define W
  set $tmp=stack_sprintf("/var/tmp/gp.tmp%d", $arg1)
  eval "shell rm -f %s", $tmp
  call gpwritebin($tmp, $arg0)
end
document W
  'W x 7' writes the GEN x to /var/tmp/gp.tmp7. Note that 'rtmp(7)'
  in a GP session would then recover the value of 'x'
end
@eprog

\newpage
\chapter{Parallelism}

PARI provides an abstraction, herafter called the MT engine, for doing
parallel computations. The exact same high level routines are used whether
the underlying communication protocol is POSIX threads or MPI and they behave
differently depending on how \kbd{libpari} was configured, specifically on
\kbd{Configure}'s \kbd{--mt} option. Sequential computation is also supported
(no \kbd{--mt} argument) which is helpful for debugging newly written
parallel code. The final section in this chapter comments a complete example.

\section{The PARI multithread interface}

\fun{void}{mt_queue_start}{struct pari_mt *pt, GEN worker} Let \kbd{worker}
be a \typ{CLOSURE} object of arity $1$.  Initialize the opaque structure
\kbd{pt} to evaluate \kbd{worker} in parallel, using \kbd{nbthreads} threads.
This allocates data in
various ways, e.g., on the PARI stack or as malloc'ed objects: you may not
collect garbage on the PARI stack starting from an earlier \kbd{avma} point
until the parallel computation is over, it could destroy something in \kbd{pt}.
All ressources allocated outside the PARI stack are freed by
\kbd{mt\_queue\_end}.

\fun{void}{mt_queue_start_lim}{struct pari_mt *pt, GEN worker, long lim}
as \kbd{mt\_queue\_start}, where \kbd{lim} is an upper bound on the number
of \kbd{tasks} to perform. Concretely the number of threads is the minimum
of \kbd{lim} and \kbd{nbthreads}. The values $0$ and $1$ of \kbd{lim} are
special:

\item $0$: no limit, equivalent to \kbd{mt\_queue\_start} (use
\kbd{nbthreads} threads).

\item $1$: no parallelism, evaluate the tasks sequentially.

\fun{void}{mt_queue_submit}{struct pari_mt *pt, long taskid, GEN task} submit
\kbd{task} to be evaluated by \kbd{worker}; use \kbd{task = NULL} if no
further task needs to be submitted. The parameter \kbd{taskid} is attached to
the \kbd{task} but not used in any way by the \kbd{worker} or the MT engine,
it will be returned to you by \kbd{mt\_queue\_get} together with the result
for the task, allowing to match up results and submitted tasks if desired.
For instance, if the tasks $(t_1,\dots, t_m)$ are known in advance, stored in
a vector, and you want to recover the evaluation results in the same order as
in that vector, you may use consecutive integers $1, \dots, m$ as
\kbd{taskid}s. If you do not care about the ordering, on the other hand, you
can just use $\kbd{taskid} = 0$ for all tasks.

The \kbd{taskid} parameter is ignored when \kbd{task} is \kbd{NULL}. It is
forbidden to call this function twice without an intervening
\kbd{mt\_queue\_get}.

\fun{GEN}{mt_queue_get}{struct pari_mt *pt, long *taskid, long *pending}
return \kbd{NULL} until \kbd{mt\_queue\_submit} has submitted
tasks for the required number (\kbd{nbthreads}) of threads; then return the
result of the evaluation by \kbd{worker} of one of the previously submitted
tasks, in random order. Set \kbd{pending} to the number of remaining pending
tasks: if this is $0$ then no more tasks are pending and it is safe to call
\tet{mt_queue_end}. Set \kbd{*taskid} to the value attached to this task by
\kbd{mt\_queue\_submit}, unless the \kbd{taskid} pointer is \kbd{NULL}. It is
forbidden to call this function twice without an intervening
\kbd{mt\_queue\_submit}.

\fun{void}{mt_queue_end}{struct pari_mt *pt} end the parallel execution
and free ressources attached to the opaque \kbd{pari\_mt} structure. For
instance malloc'ed data; in the \kbd{pthreads} interface, it would destroy
mutex locks, condition variables, etc. This must be called once there are no
longer pending tasks to avoid leaking ressources; but not before all tasks
have been processed else crashes will occur.

\fun{long}{mt_nbthreads}{void} return the effective number of parallel threads
that would be started by \tet{mt_queue_start} if it has been called in place
of \tet{mt_nbthreads}.

\section{Technical functions required by MPI}

The functions in this section are needed when writing complex independent
programs in order to support the MPI MT engine, as more flexible
complement/variants of \kbd{pari\_init} and \kbd{pari\_close}.

\fun{void}{mt_broadcast}{GEN code}: do nothing unless the MPI threading engine
is in use. In that case, evaluates the closure  \kbd{code} on all secondary
nodes. This can be used to change the state of all MPI child nodes, e.g.,
in \tet{gpinstall} run in the main thread, which allows all nodes to use the
new function.

\fun{void}{pari_mt_init}{void}
when using MPI, it is often necessary to run initialization code on the child
nodes after PARI is initialized. This is done by calling successively:

\item \tet{pari_init_opts} with the flag \tet{INIT_noIMTm}:
this initializes PARI, but not the MT engine;

\item the required initialization code;

\item \tet{pari_mt_init} to initialize the MT engine.
Note that under MPI, this function returns on the master node but enters
slave mode on the child nodes. Thus it is no longer possible to run
initialization code on the child nodes.

\fun{void}{pari_mt_close}{void}
when using MPI, calling \tet{pari_close} terminates the MPI execution
environment and it will not be possible to restart it. If this is
undesirable, call \tet{pari_close_opts} with the flag \tet{INIT_noIMTm}
instead of \kbd{pari\_close}: this closes PARI without terminating the MPI
execution environment. You may later call \kbd{pari\_mt\_close} to terminate
it. It is an error for a program to end without terminating the MPI execution
environment.

\section{A complete example}

We now proceed to an example exhibiting complex features of this
interface, in particular showing how to generate a valid \kbd{worker}.
Explanations and details follow.

\bprogfile{../examples/pari-mt.c}

We start from some arbitrary C function \kbd{Cworker} and create an
\kbd{entree} summarizing all that GP would need to know about it, in
particular

\item a GP name \kbd{\_worker}; the leading \kbd{\_} is not necessary,
we use it as a namespace mechanism grouping private functions;

\item the name of the C function;

\item and its prototype, see \kbd{install} for an introduction to Prototype
Codes.

\noindent The other three arguments ($0$, $20$ and \kbd{""}) are required in an
\kbd{entree} but not useful in our simple context: they are respectively a
valence ($0$ means ``nothing special''), a help section (20 is customary for
internal functions which need to be exported for technical reasons, see
\kbd{?20}), and a help text (no help).

Then we initialize the MT engine; doing things in this order with a two part
initialization ensures that nodes have access to our \kbd{Cworker}. We
convert the \kbd{ep} data to a \typ{CLOSURE} using \kbd{strtofunction}, which
provides a valid \kbd{worker} to \kbd{mt\_queue\_start}. This creates a
parallel evaluation queue \kbd{mt}, and we proceed to submit all tasks,
recording all results. Results are stored in the right order
by making good use of the \kbd{taskid} label, although we have no control
over \emph{when} each result is returned. We finally free all ressources
attached to the \kbd{mt} structure. If needed, we could have collected all
garbage on the PARI stack using \kbd{gerepilecopy} on the \kbd{out} array and
gone on working instead of quitting.

Note the argument passing convention for \kbd{Cworker}: the task consists of a
single vector containing all arguments as \kbd{GEN}s, which are interpreted
according to the function prototype, here \kbd{GL} so the first argument is
left as is and the second one is converted to a long integer. In more
complicated situations, this second (and possibly further) argument could
provide arbitrary evaluation contexts. In this example, we just used it as a
flag to indicate the kind of evaluation expected on the data: integer
factorization (0) or matrix determinant (1).

Note also that
\bprog
  gel(out, taskid) = mt_queue_get(&mt, &taskid, &pending);
@eprog \noindent instead of our use of a temporary \kbd{done} would have
undefined behaviour (\kbd{taskid} may be uninitialized in the left hand side).

\chapter{Cross-compiling PARI for Windows from Linux}

We use \kbd{mingw}. Please use the kit at

\kbd{https://pari.math.u-bordeaux.fr/pub/pari/windows/paricrossmingwkit.tgz}

\noindent This kit provides helper scripts and support binary to
cross-compile PARI/GP for \kbd{mingw} and \kbd{mingw64} with \kbd{readline},
\kbd{gmp} and \kbd{perl} support. To start:

\item Install the cross-compiling environment. On Debian and Ubuntu this can
be achieved by using the scripts \kbd{install32} and \kbd{install64}.

\item In this directory, type
\bprog
. ./setup
@eprog
to set the environment variables, especially \kbd{PARIKIT}. This assumes a
\kbd{bash} (or \kbd{zsh}) shell.

\item The following commands can then be executed \emph{from the toplevel of
a pari source tree}:

\bprog
mkwine32    : build the installer package (32-bit)
mkwine64    : build the installer package (64-bit)
mkwinebin32 : built the stand-alone GP binary (32-bit)
mkwinebin64 : built the stand-alone GP binary (64-bit)
@eprog

\noindent For convenience, the following variants are available for the
standalone binaries (both \kbd{mkwinebin32} and \kbd{mkwinebin64}), the
default is \kbd{rl}:
\bprog
mkwinebin64 norl : disable readline
mkwinebin64 rl   : only build readline
mkwinebin64 all  : build both readline and noreadline version
@eprog

\vfill\eject
\input index\end