Pari/GP Reference Documentation  Contents - Index - Meta commands

Vectors, matrices, linear algebra and sets


algdep   charpoly   concat   lindep   listcreate   listinsert   listkill   listput   listsort   matadjoint   matcompanion   matdet   matdetint   matdiagonal   mateigen   matfrobenius   mathess   mathilbert   mathnf   mathnfmod   mathnfmodid   matid   matimage   matimagecompl   matindexrank   matintersect   matinverseimage   matisdiagonal   matker   matkerint   matmuldiagonal   matmultodiagonal   matpascal   matrank   matrix   matrixqz   matsize   matsnf   matsolve   matsolvemod   matsupplement   mattranspose   minpoly   qfgaussred   qfjacobi   qflll   qflllgram   qfminim   qfperfection   qfrep   qfsign   setintersect   setisset   setminus   setsearch   setunion   trace   vecextract   vecsort   vector   vectorsmall   vectorv  
 
algdep(x,k,{flag = 0})  

x being real/complex, or p-adic, finds a polynomial of degree at most k with integer coefficients having x as approximate root. Note that the polynomial which is obtained is not necessarily the "correct" one. In fact it is not even guaranteed to be irreducible. One can check the closeness either by a polynomial evaluation (use subst), or by computing the roots of the polynomial given by algdep (use polroots).

Internally, lindep([1,x,...,x^k], flag) is used. If lindep is not able to find a relation and returns a lower bound for the sup norm of the smallest relation, algdep returns that bound instead. A suitable non-zero value of flag may improve on the default behaviour:

\\\\\\\\\ LLL
 ? \p200
 ? algdep(2^(1/6)+3^(1/5), 30);      \\ wrong in 3.8s
 ? algdep(2^(1/6)+3^(1/5), 30, 100); \\ wrong in 1s
 ? algdep(2^(1/6)+3^(1/5), 30, 170); \\ right in 3.3s
 ? algdep(2^(1/6)+3^(1/5), 30, 200); \\ wrong in 2.9s
 ? \p250
 ? algdep(2^(1/6)+3^(1/5), 30);      \\ right in 2.8s
 ? algdep(2^(1/6)+3^(1/5), 30, 200); \\ right in 3.4s
 \\\\\\\\\ PSLQ
 ? \p200
 ? algdep(2^(1/6)+3^(1/5), 30, -3);  \\ failure in 14s.
 ? \p250
 ? algdep(2^(1/6)+3^(1/5), 30, -3);  \\ right in 18s

Proceeding by increments of 5 digits of accuracy, algdep with default flag produces its first correct result at 205 digits, and from then on a steady stream of correct results. Interestingly enough, our PSLQ also reliably succeeds from 205 digits on (and is 5 times slower at that accuracy).

The above example is the testcase studied in a 2000 paper by Borwein and Lisonek, Applications of integer relation algorithms, Discrete Math., 217, p.65--82. The paper conludes in the superiority of the PSLQ algorithm, which either shows that PARI's implementation of PSLQ is lacking, or that its LLL is extremely good. The version of PARI tested there was 1.39, which succeeded reliably from precision 265 on, in about 60 as much time as the current version.

The library syntax is algdep0(x,k,flag,prec), where k and flag are longs. Also available is algdep(x,k,prec) (flag = 0).

charpoly(A,{v = x},{flag = 0})  

characteristic polynomial of A with respect to the variable v, i.e.determinant of v*I-A if A is a square matrix. If A is not a square matrix, it returns the characteristic polynomial of the map "multiplication by A" if A is a scalar, in particular a polmod. E.g. charpoly(I) = x^2+1.

The value of flag is only significant for matrices.

If flag = 0, the method used is essentially the same as for computing the adjoint matrix, i.e.computing the traces of the powers of A.

If flag = 1, uses Lagrange interpolation which is almost always slower.

If flag = 2, uses the Hessenberg form. This is faster than the default when the coefficients are intmod a prime or real numbers, but is usually slower in other base rings.

The library syntax is charpoly0(A,v,flag), where v is the variable number. Also available are the functions caract(A,v) (flag = 1), carhess(A,v) (flag = 2), and caradj(A,v,pt) where, in this last case, pt is a GEN* which, if not equal to NULL, will receive the address of the adjoint matrix of A (see matadjoint), so both can be obtained at once.

concat(x,{y})  

concatenation of x and y. If x or y is not a vector or matrix, it is considered as a one-dimensional vector. All types are allowed for x and y, but the sizes must be compatible. Note that matrices are concatenated horizontally, i.e.the number of rows stays the same. Using transpositions, it is easy to concatenate them vertically.

To concatenate vectors sideways (i.e.to obtain a two-row or two-column matrix), use Mat instead (see the example there). Concatenating a row vector to a matrix having the same number of columns will add the row to the matrix (top row if the vector is x, i.e.comes first, and bottom row otherwise).

The empty matrix [;] is considered to have a number of rows compatible with any operation, in particular concatenation. (Note that this is definitely not the case for empty vectors [] or []~.)

If y is omitted, x has to be a row vector or a list, in which case its elements are concatenated, from left to right, using the above rules.

? concat([1,2], [3,4])
 %1 = [1, 2, 3, 4]
 ? a = [[1,2]~, [3,4]~]; concat(a)
 %2 =
 [1 3]
 
 [2 4]
 
 ? concat([1,2; 3,4], [5,6]~)
 %3 =
 [1 2 5]
 
 [3 4 6]
 ? concat([%, [7,8]~, [1,2,3,4]])
 %5 =
 [1 2 5 7]
 
 [3 4 6 8]
 
 [1 2 3 4]

The library syntax is concat(x,y).

lindep(x,{flag = 0})  

x being a vector with p-adic or real/complex coefficients, finds a small integral linear combination among these coefficients.

If x is p-adic, flag is meaningless and the algorithm LLL-reduces a suitable (dual) lattice.

Otherwise, the value of flag determines the algorithm used; in the current version of PARI, we suggest to use non-negative values, since it is by far the fastest and most robust implementation. See the detailed example in Section [Label: se:algdep] ( algdep).

If flag >= 0, uses a floating point (variable precision) LLL algorithm. This is in general much faster than the other variants. If flag = 0 the accuracy is chosen internally using a crude heuristic. If flag > 0 the computation is done with an accuracy of flag decimal digits. In that case, the parameter flag should be between 0.6 and 0.9 times the number of correct decimal digits in the input.

If flag = -1, uses a variant of the LLL algorithm due to Hastad, Lagarias and Schnorr (STACS 1986). If the precision is too low, the routine may enter an infinite loop.

If flag = -2, x is allowed to be (and in any case interpreted as) a matrix. Returns a non trivial element of the kernel of x, or 0 if x has trivial kernel. The element is defined over the field of coefficients of x, and is in general not integral.

If flag = -3, uses the PSLQ algorithm. This may return a real number B, indicating that the input accuracy was exhausted and that no relation exist whose sup norm is less than B.

If flag = -4, uses an experimental 2-level PSLQ, which does not work at all. (Should be rewritten.)

The library syntax is lindep0(x,flag,prec). Also available is lindep(x,prec) (flag = 0).

listcreate(n)  

creates an empty list of maximal length n.

This function is useless in library mode.

listinsert(list,x,n)  

inserts the object x at position n in list (which must be of type t_LIST). All the remaining elements of list (from position n+1 onwards) are shifted to the right. This and listput are the only commands which enable you to increase a list's effective length (as long as it remains under the maximal length specified at the time of the listcreate).

This function is useless in library mode.

listkill(list)  

kill list. This deletes all elements from list and sets its effective length to 0. The maximal length is not affected.

This function is useless in library mode.

listput(list,x,{n})  

sets the n-th element of the list list (which must be of type t_LIST) equal to x. If n is omitted, or greater than the list current effective length, just appends x. This and listinsert are the only commands which enable you to increase a list's effective length (as long as it remains under the maximal length specified at the time of the listcreate).

If you want to put an element into an occupied cell, i.e.if you don't want to change the effective length, you can consider the list as a vector and use the usual list[n] = x construct.

This function is useless in library mode.

listsort(list,{flag = 0})  

sorts list (which must be of type t_LIST) in place. If flag is non-zero, suppresses all repeated coefficients. This is much faster than the vecsort command since no copy has to be made.

This function is useless in library mode.

matadjoint(x)  

adjoint matrix of x, i.e.the matrix y of cofactors of x, satisfying x*y = det(x)*\Id. x must be a (non-necessarily invertible) square matrix.

The library syntax is adj(x).

matcompanion(x)  

the left companion matrix to the polynomial x.

The library syntax is assmat(x).

matdet(x,{flag = 0})  

determinant of x. x must be a square matrix.

If flag = 0, uses Gauss-Bareiss.

If flag = 1, uses classical Gaussian elimination, which is better when the entries of the matrix are reals or integers for example, but usually much worse for more complicated entries like multivariate polynomials.

The library syntax is det(x) (flag = 0) and det2(x) (flag = 1).

matdetint(x)  

x being an m x n matrix with integer coefficients, this function computes a multiple of the determinant of the lattice generated by the columns of x if it is of rank m, and returns zero otherwise. This function can be useful in conjunction with the function mathnfmod which needs to know such a multiple. To obtain the exact determinant (assuming the rank is maximal), you can compute matdet(mathnfmod(x, matdetint(x))).

Note that as soon as one of the dimensions gets large (m or n is larger than 20, say), it will often be much faster to use mathnf(x, 1) or mathnf(x, 4) directly.

The library syntax is detint(x).

matdiagonal(x)  

x being a vector, creates the diagonal matrix whose diagonal entries are those of x.

The library syntax is diagonal(x).

mateigen(x)  

gives the eigenvectors of x as columns of a matrix.

The library syntax is eigen(x).

matfrobenius(M,{flag = 0},{v = x})  

returns the Frobenius form of the square matrix M. If flag = 1, returns only the elementary divisors as a vectr of polynomials in the variable v. If flag = 2, returns a two-components vector [F,B] where F is the Frobenius form and B is the basis change so that M = B^{-1}FB.

The library syntax is matfrobenius(M,flag,v), where v is the variable number.

mathess(x)  

Hessenberg form of the square matrix x.

The library syntax is hess(x).

mathilbert(x)  

x being a long, creates the Hilbert matrixof order x, i.e.the matrix whose coefficient (i,j) is 1/ (i+j-1).

The library syntax is mathilbert(x).

mathnf(x,{flag = 0})  

if x is a (not necessarily square) matrix with integer entries, finds the upper triangular Hermite normal form of x. If the rank of x is equal to its number of rows, the result is a square matrix. In general, the columns of the result form a basis of the lattice spanned by the columns of x.

If flag = 0, uses the naive algorithm. This should never be used if the dimension is at all large (larger than 10, say). It is recommanded to use either mathnfmod(x, matdetint(x)) (when x has maximal rank) or mathnf(x, 1). Note that the latter is in general faster than mathnfmod, and also provides a base change matrix.

If flag = 1, uses Batut's algorithm, which is much faster than the default. Outputs a two-component row vector [H,U], where H is the upper triangular Hermite normal form of x defined as above, and U is the unimodular transformation matrix such that xU = [0|H]. U has in general huge coefficients, in particular when the kernel is large.

If flag = 3, uses Batut's algorithm, but outputs [H,U,P], such that H and U are as before and P is a permutation of the rows such that P applied to xU gives H. The matrix U is smaller than with flag = 1, but may still be large.

If flag = 4, as in case 1 above, but uses a heuristic variant of LLL reduction along the way. The matrix U is in general close to optimal (in terms of smallest L_2 norm), but the reduction is slower than in case 1.

The library syntax is mathnf0(x,flag). Also available are hnf(x) (flag = 0) and hnfall(x) (flag = 1). To reduce huge (say 400 x 400 and more) relation matrices (sparse with small entries), you can use the pair hnfspec / hnfadd. Since this is rather technical and the calling interface may change, they are not documented yet. Look at the code in basemath/alglin1.c.

mathnfmod(x,d)  

if x is a (not necessarily square) matrix of maximal rank with integer entries, and d is a multiple of the (non-zero) determinant of the lattice spanned by the columns of x, finds the upper triangular Hermite normal form of x.

If the rank of x is equal to its number of rows, the result is a square matrix. In general, the columns of the result form a basis of the lattice spanned by the columns of x. This is much faster than mathnf when d is known.

The library syntax is hnfmod(x,d).

mathnfmodid(x,d)  

outputs the (upper triangular) Hermite normal form of x concatenated with d times the identity matrix. Assumes that x has integer entries.

The library syntax is hnfmodid(x,d).

matid(n)  

creates the n x n identity matrix.

The library syntax is matid(n) where n is a long.

Related functions are gscalmat(x,n), which creates x times the identity matrix (x being a GEN and n a long), and gscalsmat(x,n) which is the same when x is a long.

matimage(x,{flag = 0})  

gives a basis for the image of the matrix x as columns of a matrix. A priori the matrix can have entries of any type. If flag = 0, use standard Gauss pivot. If flag = 1, use matsupplement.

The library syntax is matimage0(x,flag). Also available is image(x) (flag = 0).

matimagecompl(x)  

gives the vector of the column indices which are not extracted by the function matimage. Hence the number of components of matimagecompl(x) plus the number of columns of matimage(x) is equal to the number of columns of the matrix x.

The library syntax is imagecompl(x).

matindexrank(x)  

x being a matrix of rank r, gives two vectors y and z of length r giving a list of rows and columns respectively (starting from 1) such that the extracted matrix obtained from these two vectors using vecextract(x,y,z) is invertible.

The library syntax is indexrank(x).

matintersect(x,y)  

x and y being two matrices with the same number of rows each of whose columns are independent, finds a basis of the Q-vector space equal to the intersection of the spaces spanned by the columns of x and y respectively. See also the function idealintersect, which does the same for free Z-modules.

The library syntax is intersect(x,y).

matinverseimage(M,y)  

gives a column vector belonging to the inverse image z of the column vector or matrix y by the matrix M if one exists (i.e such that Mz = y), the empty vector otherwise. To get the complete inverse image, it suffices to add to the result any element of the kernel of x obtained for example by matker.

The library syntax is inverseimage(x,y).

matisdiagonal(x)  

returns true (1) if x is a diagonal matrix, false (0) if not.

The library syntax is isdiagonal(x), and this returns a long integer.

matker(x,{flag = 0})  

gives a basis for the kernel of the matrix x as columns of a matrix. A priori the matrix can have entries of any type.

If x is known to have integral entries, set flag = 1.

Note: The library function FpM_ker(x, p), where x has integer entries reduced mod p and p is prime, is equivalent to, but orders of magnitude faster than, matker(x*Mod(1,p)) and needs much less stack space. To use it under gp, type install(FpM_ker, GG) first.

The library syntax is matker0(x,flag). Also available are ker(x) (flag = 0), keri(x) (flag = 1).

matkerint(x,{flag = 0})  

gives an LLL-reduced Z-basis for the lattice equal to the kernel of the matrix x as columns of the matrix x with integer entries (rational entries are not permitted).

If flag = 0, uses a modified integer LLL algorithm.

If flag = 1, uses matrixqz(x,-2). If LLL reduction of the final result is not desired, you can save time using matrixqz(matker(x),-2) instead.

The library syntax is matkerint0(x,flag). Also available is kerint(x) (flag = 0).

matmuldiagonal(x,d)  

product of the matrix x by the diagonal matrix whose diagonal entries are those of the vector d. Equivalent to, but much faster than x* matdiagonal(d).

The library syntax is matmuldiagonal(x,d).

matmultodiagonal(x,y)  

product of the matrices x and y assuming that the result is a diagonal matrix. Much faster than x*y in that case. The result is undefined if x*y is not diagonal.

The library syntax is matmultodiagonal(x,y).

matpascal(x,{q})  

creates as a matrix the lower triangular Pascal triangle of order x+1 (i.e.with binomial coefficients up to x). If q is given, compute the q-Pascal triangle (i.e.using q-binomial coefficients).

The library syntax is matqpascal(x,q), where x is a long and q = NULL is used to omit q. Also available is matpascal(x).

matrank(x)  

rank of the matrix x.

The library syntax is rank(x), and the result is a long.

matrix(m,n,{X},{Y},{expr = 0})  

creation of the m x n matrix whose coefficients are given by the expression expr. There are two formal parameters in expr, the first one (X) corresponding to the rows, the second (Y) to the columns, and X goes from 1 to m, Y goes from 1 to n. If one of the last 3 parameters is omitted, fill the matrix with zeroes.

The library syntax is matrice(GEN nlig,GEN ncol,entree *e1,entree *e2,char *expr).

matrixqz(x,p)  

x being an m x n matrix with m >= n with rational or integer entries, this function has varying behaviour depending on the sign of p:

If p >= 0, x is assumed to be of maximal rank. This function returns a matrix having only integral entries, having the same image as x, such that the GCD of all its n x n subdeterminants is equal to 1 when p is equal to 0, or not divisible by p otherwise. Here p must be a prime number (when it is non-zero). However, if the function is used when p has no small prime factors, it will either work or give the message "impossible inverse modulo" and a non-trivial divisor of p.

If p = -1, this function returns a matrix whose columns form a basis of the lattice equal to Z^n intersected with the lattice generated by the columns of x.

If p = -2, returns a matrix whose columns form a basis of the lattice equal to Z^n intersected with the Q-vector space generated by the columns of x.

The library syntax is matrixqz0(x,p).

matsize(x)  

x being a vector or matrix, returns a row vector with two components, the first being the number of rows (1 for a row vector), the second the number of columns (1 for a column vector).

The library syntax is matsize(x).

matsnf(X,{flag = 0})  

if X is a (singular or non-singular) matrix outputs the vector of elementary divisors of X (i.e.the diagonal of the Smith normal form of X).

The binary digits of flag mean:

1 (complete output): if set, outputs [U,V,D], where U and V are two unimodular matrices such that UXV is the diagonal matrix D. Otherwise output only the diagonal of D.

2 (generic input): if set, allows polynomial entries, in which case the input matrix must be square. Otherwise, assume that X has integer coefficients with arbitrary shape.

4 (cleanup): if set, cleans up the output. This means that elementary divisors equal to 1 will be deleted, i.e.outputs a shortened vector D' instead of D. If complete output was required, returns [U',V',D'] so that U'XV' = D' holds. If this flag is set, X is allowed to be of the form D or [U,V,D] as would normally be output with the cleanup flag unset.

The library syntax is matsnf0(X,flag). Also available is smith(X) (flag = 0).

matsolve(x,y)  

x being an invertible matrix and y a column vector, finds the solution u of x*u = y, using Gaussian elimination. This has the same effect as, but is a bit faster, than x^{-1}*y.

The library syntax is gauss(x,y).

matsolvemod(m,d,y,{flag = 0})  

m being any integral matrix, d a vector of positive integer moduli, and y an integral column vector, gives a small integer solution to the system of congruences sum_i m_{i,j}x_j = y_i (mod d_i) if one exists, otherwise returns zero. Shorthand notation: y (resp.d) can be given as a single integer, in which case all the y_i (resp.d_i) above are taken to be equal to y (resp.d).

  ? m = [1,2;3,4];
   ? matsolvemod(m, [3,4], [1,2]~)
   %2 = [-2, 0]~
   ? matsolvemod(m, 3, 1) \\ m X = [1,1]~ over F_3
   %3 = [-1, 1]~

If flag = 1, all solutions are returned in the form of a two-component row vector [x,u], where x is a small integer solution to the system of congruences and u is a matrix whose columns give a basis of the homogeneous system (so that all solutions can be obtained by adding x to any linear combination of columns of u). If no solution exists, returns zero.

The library syntax is matsolvemod0(m,d,y,flag). Also available are gaussmodulo(m,d,y) (flag = 0) and gaussmodulo2(m,d,y) (flag = 1).

matsupplement(x)  

assuming that the columns of the matrix x are linearly independent (if they are not, an error message is issued), finds a square invertible matrix whose first columns are the columns of x, i.e.supplement the columns of x to a basis of the whole space.

The library syntax is suppl(x).

mattranspose(x)  

or x~: transpose of x. This has an effect only on vectors and matrices.

The library syntax is gtrans(x).

minpoly(A,{v = x},{flag = 0})  

minimal polynomial of A with respect to the variable v., i.e. the monic polynomial P of minimal degree (in the variable v) such that P(A) = 0.

The library syntax is minpoly(A,v), where v is the variable number.

qfgaussred(q)  

decomposition into squares of the quadratic form represented by the symmetric matrix q. The result is a matrix whose diagonal entries are the coefficients of the squares, and the non-diagonal entries represent the bilinear forms. More precisely, if (a_{ij}) denotes the output, one has q(x) = sum_i a_{ii} (x_i + sum_{j > i} a_{ij} x_j)^2

The library syntax is sqred(x).

qfjacobi(x)  

x being a real symmetric matrix, this gives a vector having two components: the first one is the vector of eigenvalues of x, the second is the corresponding orthogonal matrix of eigenvectors of x. The method used is Jacobi's method for symmetric matrices.

The library syntax is jacobi(x).

qflll(x,{flag = 0})  

LLL algorithm applied to the columns of the matrix x. The columns of x must be linearly independent, unless specified otherwise below. The result is a unimodular transformation matrix T such that x.T is an LLL-reduced basis of the lattice generated by the column vectors of x.

If flag = 0 (default), the computations are done with floating point numbers, using Householder matrices for orthogonalization. If x has integral entries, then computations are nonetheless approximate, with precision varying as needed (Lehmer's trick, as generalized by Schnorr).

If flag = 1, it is assumed that x is integral. The computation is done entirely with integers. In this case, x needs not be of maximal rank, but if it is not, T will not be square. This is slower and no more accurate than flag = 0 above if x has small dimension (say 100 or less).

If flag = 2, x should be an integer matrix whose columns are linearly independent. Returns a partially reduced basis for x, using an unpublished algorithm by Peter Montgomery: a basis is said to be partially reduced if |v_i ± v_j| >= |v_i| for any two distinct basis vectors v_i, v_j.

This is significantly faster than flag = 1, esp. when one row is huge compared to the other rows. Note that the resulting basis is not LLL-reduced in general.

If flag = 4, x is assumed to have integral entries, but needs not be of maximal rank. The result is a two-component vector of matrices: the columns of the first matrix represent a basis of the integer kernel of x (not necessarily LLL-reduced) and the second matrix is the transformation matrix T such that x.T is an LLL-reduced Z-basis of the image of the matrix x.

If flag = 5, case as case 4, but x may have polynomial coefficients.

If flag = 8, same as case 0, but x may have polynomial coefficients.

The library syntax is qflll0(x,flag,prec). Also available are lll(x,prec) (flag = 0), lllint(x) (flag = 1), and lllkerim(x) (flag = 4).

qflllgram(G,{flag = 0})  

same as qflll, except that the matrix G = x~ * x is the Gram matrix of some lattice vectors x, and not the coordinates of the vectors themselves. In particular, G must now be a square symmetric real matrix, corresponding to a positive definite quadratic form. The result is a unimodular transformation matrix T such that x.T is an LLL-reduced basis of the lattice generated by the column vectors of x.

If flag = 0 (default): the computations are done with floating point numbers, using Householder matrices for orthogonalization. If G has integral entries, then computations are nonetheless approximate, with precision varying as needed (Lehmer's trick, as generalized by Schnorr).

If flag = 1: G has integer entries, still positive but not necessarily definite (i.ex needs not have maximal rank). The computations are all done in integers and should be slower than the default, unless the latter triggers accuracy problems.

flag = 4: G has integer entries, gives the kernel and reduced image of x.

flag = 5: same as case 4, but G may have polynomial coefficients.

The library syntax is qflllgram0(G,flag,prec). Also available are lllgram(G,prec) (flag = 0), lllgramint(G) (flag = 1), and lllgramkerim(G) (flag = 4).

qfminim(x,{b},{m},{flag = 0})  

x being a square and symmetric matrix representing a positive definite quadratic form, this function deals with the vectors of x whose norm is less than or equal to b, enumerated using the Fincke-Pohst algorithm. The function searches for the minimal non-zero vectors if b is omitted. The precise behaviour depends on flag.

If flag = 0 (default), seeks at most 2m vectors. The result is a three-component vector, the first component being the number of vectors found, the second being the maximum norm found, and the last vector is a matrix whose columns are the vectors found, only one being given for each pair ± v (at most m such pairs). The vectors are returned in no particular order. In this variant, an explicit m must be provided.

If flag = 1, ignores m and returns the first vector whose norm is less than b. In this variant, an explicit b must be provided.

In both these cases, x is assumed to have integral entries. The implementation uses low precision floating point computations for maximal speed, which gives incorrect result when x has large entries. (The condition is checked in the code and the routine will raise an error if large rounding errors occur.) A more robust, but much slower, implementation is chosen if the following flag is used:

If flag = 2, x can have non integral real entries. In this case, if b is omitted, the "minimal" vectors only have approximately the same norm. If b is omitted, m is an upper bound for the number of vectors that will be stored and returned, but all minimal vectors are nevertheless enumerated. If m is omitted, all vectors found are stored and returned; note that this may be a huge vector!

The library syntax is qfminim0(x,b,m,flag,prec), also available are minim(x,b,m) (flag = 0), minim2(x,b,m) (flag = 1). In all cases, an omitted b or m is coded as NULL.

qfperfection(x)  

x being a square and symmetric matrix with integer entries representing a positive definite quadratic form, outputs the perfection rank of the form. That is, gives the rank of the family of the s symmetric matrices v_iv_i^t, where s is half the number of minimal vectors and the v_i (1 <= i <= s) are the minimal vectors.

As a side note to old-timers, this used to fail bluntly when x had more than 5000 minimal vectors. Beware that the computations can now be very lengthy when x has many minimal vectors.

The library syntax is perf(x).

qfrep(q, B, {flag = 0})  

q being a square and symmetric matrix with integer entries representing a positive definite quadratic form, outputs the vector whose i-th entry, 1 <= i <= B is half the number of vectors v such that q(v) = i. This routine uses a naive algorithm based on qfminim, and will fail if any entry becomes larger than 2^{31}.

The binary digits of flag mean:

* 1: count vectors of even norm from 1 to 2B.

* 2: return a t_VECSMALL instead of a t_GEN

The library syntax is qfrep0(q, B, flag).

qfsign(x)  

signature of the quadratic form represented by the symmetric matrix x. The result is a two-component vector.

The library syntax is signat(x).

setintersect(x,y)  

intersection of the two sets x and y.

The library syntax is setintersect(x,y).

setisset(x)  

returns true (1) if x is a set, false (0) if not. In PARI, a set is simply a row vector whose entries are strictly increasing. To convert any vector (and other objects) into a set, use the function Set.

The library syntax is setisset(x), and this returns a long.

setminus(x,y)  

difference of the two sets x and y, i.e.set of elements of x which do not belong to y.

The library syntax is setminus(x,y).

setsearch(x,y,{flag = 0})  

searches if y belongs to the set x. If it does and flag is zero or omitted, returns the index j such that x[j] = y, otherwise returns 0. If flag is non-zero returns the index j where y should be inserted, and 0 if it already belongs to x (this is meant to be used in conjunction with listinsert).

This function works also if x is a sorted list (see listsort).

The library syntax is setsearch(x,y,flag) which returns a long integer.

setunion(x,y)  

union of the two sets x and y.

The library syntax is setunion(x,y).

trace(x)  

this applies to quite general x. If x is not a matrix, it is equal to the sum of x and its conjugate, except for polmods where it is the trace as an algebraic number.

For x a square matrix, it is the ordinary trace. If x is a non-square matrix (but not a vector), an error occurs.

The library syntax is gtrace(x).

vecextract(x,y,{z})  

extraction of components of the vector or matrix x according to y. In case x is a matrix, its components are as usual the columns of x. The parameter y is a component specifier, which is either an integer, a string describing a range, or a vector.

If y is an integer, it is considered as a mask: the binary bits of y are read from right to left, but correspond to taking the components from left to right. For example, if y = 13 = (1101)_2 then the components 1,3 and 4 are extracted.

If y is a vector, which must have integer entries, these entries correspond to the component numbers to be extracted, in the order specified.

If y is a string, it can be

* a single (non-zero) index giving a component number (a negative index means we start counting from the end).

* a range of the form "a..b", where a and b are indexes as above. Any of a and b can be omitted; in this case, we take as default values a = 1 and b = -1, i.e.the first and last components respectively. We then extract all components in the interval [a,b], in reverse order if b < a.

In addition, if the first character in the string is ^, the complement of the given set of indices is taken.

If z is not omitted, x must be a matrix. y is then the line specifier, and z the column specifier, where the component specifier is as explained above.

? v = [a, b, c, d, e];
 ? vecextract(v, 5)          \\ mask

 %1 = [a, c]  ? vecextract(v, [4, 2, 1]) \\ component list

 %2 = [d, b, a]  ? vecextract(v, "2..4") \\ interval

 %3 = [b, c, d]  ? vecextract(v, "-1..-3") \\ interval + reverse order

 %4 = [e, d, c]  ? vecextract(v, "^2") \\ complement

 %5 = [a, c, d, e]  ? vecextract(matid(3), "2..", "..")  %6 =  [0 1 0]    [0 0 1]

The library syntax is extract(x,y) or matextract(x,y,z).

vecsort(x,{k},{flag = 0})  

sorts the vector x in ascending order, using a mergesort method. x must be a vector, and its components integers, reals, or fractions.

If k is present and is an integer, sorts according to the value of the k-th subcomponents of the components ofx. Note that mergesort is stable, hence is the initial ordering of "equal" entries (with respect to the sorting criterion) is not changed.

k can also be a vector, in which case the sorting is done lexicographically according to the components listed in the vector k. For example, if k = [2,1,3], sorting will be done with respect to the second component, and when these are equal, with respect to the first, and when these are equal, with respect to the third.

The binary digits of flag mean:

* 1: indirect sorting of the vector x, i.e.if x is an n-component vector, returns a permutation of [1,2,...,n] which applied to the components of x sorts x in increasing order. For example, vecextract(x, vecsort(x,,1)) is equivalent to vecsort(x).

* 2: sorts x by ascending lexicographic order (as per the lex comparison function).

* 4: use descending instead of ascending order.

The library syntax is vecsort0(x,k,flag). To omit k, use NULL instead. You can also use the simpler functions

sort(x) ( = vecsort0(x, NULL,0)).

indexsort(x) ( = vecsort0(x, NULL,1)).

lexsort(x) ( = vecsort0(x, NULL,2)).

Also available are sindexsort(x) and sindexlexsort(x) which return a t_VECSMALL v, where v[1]...v[n] contain the indices.

vector(n,{X},{expr = 0})  

creates a row vector (type t_VEC) with n components whose components are the expression expr evaluated at the integer points between 1 and n. If one of the last two arguments is omitted, fill the vector with zeroes.

Avoid modifying X within expr; if you do, the formal variable still runs from 1 to n. In particular, vector(n,i,expr) is not equivalent to

    v = vector(n)
     for (i = 1, n, v[i] = expr)

as the following example shows:

    n = 3
     v = vector(n); vector(n, i, i++)            ----> [2, 3, 4]
     v = vector(n); for (i = 1, n, v[i] = i++)   ----> [2, 0, 4]

The library syntax is vecteur(GEN nmax, entree *ep, char *expr).

vectorsmall(n,{X},{expr = 0})  

creates a row vector of small integers (type t_VECSMALL) with n components whose components are the expression expr evaluated at the integer points between 1 and n. If one of the last two arguments is omitted, fill the vector with zeroes.

The library syntax is vecteursmall(GEN nmax, entree *ep, char *expr).

vectorv(n,X,expr)  

as vector, but returns a column vector (type t_COL).

The library syntax is vvecteur(GEN nmax, entree *ep, char *expr).