LOW LEVEL FUNCTIONS

OVERVIEW

The low-level functions, prefixed by ll, do not provide a coherent calling interface. They are directly accessible but one must be careful, for example, to allocate the right amount of memory for the result.

DESCRIPTION

The following holds for all low level functions except when explicitly stated otherwise: the memory to store the result must be allocated before the call to the low level function ; the allocated memory must have at least the size to hold the final result ; the memory for the result must not overlap with the memory occupied by the input ; the memory regions occupied by all or some elements of the input can be the same ; the input is not modified by the function.

Note that many details must be taken care of when using low level functions. For exemple if you have two polynomials and add them with ll_array_add then it is up to you to compute the degree of the resulting polynomial.

The prefix of functions indicates on what type of object it acts on.

VECTORS AND ARRAYS

Vectors and arrays are represented as simple C arrays of R, concretely they are simple pointers R*. There is no practical differences between vectors and arrays. A function prefixed by ll_vec manipulates arrays of the same size, like ll_vec_add or ll_vec_sub whereas a function prefixed by ll_array can manipulates vectors of different sizes, like ll_array_add.

MATRICES

Matrices are represented by a C array of R, namely they are simple pointers R*. Thus you must give to functions manipulating matrices the their sizes, the number of rows and the numbers of columns.

UNIVARIATE POLYNOMIALS

Univariate polynomials are represented as simple C arrays of R, concretely they are pointers R*. You must pass to functions manipulating polynomials their sizes. The size or length of a polynomial of degree d is defined to be d + 1.

DENSE BIVARIATE POLYNOMIALS

Bivariate polynomials are represented as simple C arrays of R, concretely they are pointers R*. You must pass their size to functions manipulating them. The size of a bivariate polynomial of degree dx in X and dy in Y is defined to be the pair (dx + 1,dy + 1).

LINKED LIST OF VECTORS

Linked lists of vectors are used for example by ll_dbpol_Xroots to return the list of roots.

struct ll_lili_vec {
  R *a;
  int na;
  struct ll_lili_vec *next,*prev;
};

You can manipulate directly the a and na members where na is thte size of a. The prev and next members shall not modified directly, use the appropriate functions instead (see below).

ALLOCATION FUNCTIONS

All allocation function, usually containing malloc in their name, return a pointer to allocated memory. In case the allocation fails, SIGABRT is raised. All function that use an *alloc(3) function from the standard library shall respect this behavior. You can catch this signal with sigaction(2). By default it will terminate silently your program.

LIST OF LOW LEVEL FUNCTIONS

VECTORS AND ARRAYS

void ll_vec_add(R *r,R *a,R *b,uma n);

Set r to a + b where a and b have size n. You can have r == a and/or r == b.

void ll_vec_neg(R *r,R *a,uma n);

Set r to -a where a has size n. You can have r == a.

void ll_vec_sub(R *r,R *a,R *b,uma n);

Set r to a - b where a and b have size n. You can have r == a and/or r == b.

void ll_vec_sub_sub(R *r,R *a,R *b,uma n);

Set r to r - a - b where r, a and b have size n. You can have r == a and/or r == b.

void ll_vec_add_mul(R *r,R *a,R *b,uma n);

Set r to r + (a .* b) where r, a and b have size n and ."* denote the component-wise product. You can have r == a and/or r == b.

void ll_vec_mul_by_cte(R *r,R *a,uma na,R b);

Set r to a * b where a has size na. You can have r == a.

void ll_vec_add_mul_by_cte(R *r,R *a,uma na,R b);

Set r to r + (a * b) where a has size na.

void ll_vec_sub_mul(R *r,R *a,R *b,uma n);

Set r to r - (a .* b) where r, a and b have size na and ."* denote the component-wise product.

void ll_vec_init(R *a,uma na);

Initialize and set to zero (with a call to R_init) each component of a of size na.

void ll_vec_init_random(R *a,uma na);

Initialize and randomize (with calls to R_init and R_random) each component of a of size na.

void ll_vec_clear(R *a,uma na);

Free the memory occupied by each component of a of size na.

uma ll_vec_hamming_weight(R *a,uma na);

Return the Hamming weight (the number of nonzero components) of a of size na.

uma ll_vec_hamming_distance(R *a,R *b,uma n);

Return the Hamming distance between a and b, both of size n. You can have a == b, in which case it will return 0.

void ll_vec_zero(R *a,uma na);

Set each component of a to zero where a has size na.

uma ll_vec_iszero(R *a,uma na);

Return 0 if and only if a has a nonzero component where a has size na.

void ll_vec_random(R *a,uma na);

Randomize (with R_random) each component of a of size na.

void ll_vec_random_hamming_weight(R *r,uma nr,uma w);

Set r to a random vector of size nr of Hamming weight w.

void ll_vec_copy(R *r,R *a,uma na);

Perform a deep copy of a of size na into r.

uma ll_vec_true_size(R *a,uma na);

Return the greatest integer i <= na such that a[i - 1] is nonzero. If all the components of a are zero, 1 is returned.

R *ll_vec_malloc_init(uma nr);

Return a pointer to a vector of size nr where each component is set to zero.

R *ll_vec_realloc(R *a,uma na,uma nna);

Changes the size nna of the memory block occupied by the vector a of size na and return a pointer to a vector of size nna. The content of a will be unchanged. If nna > na the new vector will be completed with zeros.

R *ll_vec_malloc_copy(R *a,uma na);

Perform a deep copy of a of size na and return a pointer

R *ll_vec_malloc_random(uma nr);

Return a pointer to a random vector of size nr.

void ll_vec_free(R *r,uma nr);

Free the memory occupied by the vector r of size nr.

uma ll_vec_cmp(R *a,R *b,uma n);

Return 0 if and only if a and b, both of size n, are equal.

void ll_vec_print(FILE *out,R *a,uma na);

Pruma to out the vector a of size na.

void ll_array_add(R *r,R *a,uma na,R *b,uma nb);

Set r to a + b where na and nb are the sizes of a and b. If a and b have not the same size the shortest is considered to be extended with zeros. You can have r == a and/or r == b.

void ll_array_sub(R *r,R *a,uma na,R *b,uma nb);

Set r to a - b where na and nb are the sizes of a and b. If a and b have not the same size the shortest is considered to be extended with zeros. You can have r == a and/or r == b.

void ll_array_add_mul_by_cte(R *r,uma nr,R *a,uma na,R b);

Set r to r + (a * b) where nr and na are the sizes of r and a. If r and a have not the same size the shortest is considered extended with zeros. You can have r == a, b must not be in r or a.

MATRICES

void ll_mat_copy(R *r,R *a,uma ra,uma ca);

Do a deep copy of a of size (ar,ac) to r.

R *ll_mat_malloc_init(uma rm,uma cm);

Return a pointer to a region of (rm x cm) R's.

R *ll_mat_malloc_random(uma mr,uma mc);

Return a pointer to a region of (rm x cm) random R's.

void ll_mat_free(R *m,uma rm,uma cm);

Free m of size rm x cm.

void ll_mat_transpose(R *r,R *a,uma ra,uma ca);

Set r to be the transpose of a which has size (ra,ca).

void ll_mat_print(FILE *out,R *a,uma ra,uma ca);

Pruma to out the matrix a of size (ra,ca).

void ll_mat_put(R *r,uma cr,uma i0,uma j0,R *a,uma ra,uma ca);

Put a of size (ra,ca) into r with cr columns at position (i0,j0).

UNIVARIATE POLYNOMIALS

int ll_upol_iszero(R *a,uma na);

Return nonzero if and only if na = 1 and a[0] = 0. Note that if the degree of a is na - 1 this function returns nonzero if and only if a with size na is the zero polynomial.

void ll_upol_binpow_mod(R *r,uma *nr,R *a,uma na,mpz_t e,R *m,uma nm);

Set r to a(X) ^ e mod m where a has size na and m has size nm and set nr to the size of r. Note that e must be a GMP integer. The repeated squaring algorithm is used.

void ll_upol_derivative_naive(R *r,R *a,uma na);

Set r to be the formal derivative of a of size na. The naive algorithm is used.

void ll_upol_derivative(R *r,R *a,uma na);

Set r to be the formal derivative of a of size na.

void ll_upol_discrepancy(R r,R *a,uma na,uma m,R x);

Set r to be the evaluation at x of the m-th Hasse derivative of a with size na.

void ll_upol_div_by_affine(R *r,R *a,uma na,R b);

Set r to the quotient of the Euclidean division of a by X + b where a has size na.

void ll_upol_eea(R *gcd,uma *ngcd,R *u,uma *nu,R *v,uma *nv,R *a,uma na,R *b,uma nb);

Set gcd of size ngcd, u of size nu and v of size nv to be the GCD and the Bezout coefficients of a of size na and b of size nb. Namely it returns gcd, u and v such that u(X) a(X) + v(X) b(X) = gcd(X). Note that u must point to a region of at least nb allocated R's and v must point to a region of at least na allocated R's.

struct ll_vec_lili *ll_upol_eqdegfact(R *a,uma na,uma d);

Return all the irredusible factors of a of size na of degree d. Note that all irreducible factors of a must have degree d and that a must have degree na - 1.

void ll_upol_trace_poly(R *r,uma *nr,uma d,R *a,uma na,R *m,uma nm)

Set r of size nr to be the trace polynomial of degree d of a of size na in the sense of exercise 14.16, page 413 of Modern Computer Algebra, 2nd edition by Von zur Gathen and Gerhard.

void ll_upol_cantor_zassenhaus(R *r,R *a,uma na,uma d);

Set r to be a proper factor of a. All irreducible factors of a must have degree d. Note that r may not be irreducible and that a must be monic and has degree na - 1. The algorithm of Cantor and Zassenhaus is used.

void ll_upol_eqdegsplit(R *r,R *a,uma na);

Set r to be a proper factor of a. Note that r may not be irreducible. Note that a must be monic and have degree na - 1.

void ll_upol_quorem_naive(R *r,uma *nr,R *q,R *a,uma na,R *b,uma nb);

Set r of size nr and q to be the rest and the quotient of the Euclidean division of a of size na by b of size nb. Note that the degree of b must be nb - 1. The schoolbook algorithm is used. You can have r == a.

void ll_upol_quo_naive(R *q,R *a,uma na,R *b,uma nb);

Set q to be the quotient of the Euclidean division of a of size na by b of size nb. Note that the degree of b must be nb - 1. The schoolbook algorithm is used.

void ll_upol_rem_naive(R *r,uma *nr,R *a,uma na,R *b,uma nb);

Set r of size nr to be the rest of the Euclidean division of a of size na by b of size nb. Note that the degree of b must be nb - 1. The schoolbook algorithm is used. You can have r == a.

void ll_upol_quorem(R *r,uma *nr,R *q,R *a,uma na,R *b,uma nb);

Set r of size nr and q to be the rest and the quotient of the Euclidean division of a of size na by b of size nb. Note that the degree of b must be nb - 1 and that r must point to a region of at least na allocated R's.

void ll_upol_quo(R *q,R *a,uma na,R *b,uma nb);

Set q to be the quotient of the Euclidean division of a of size na by b of size nb. Note that the degree of b must be nb - 1 and that r must point to a region of at least na allocated R's.

void ll_upol_rem(R *r,uma *nr,R *a,uma na,R *b,uma nb);

Set r of size nr to be the rest of the Euclidean division of a of size na by b of size nb. Note that the degree of b must be nb - 1 and that r must point to a region of at least na allocated R's.

void ll_upol_eval(R r,R *a,uma na,R x);

Set r to be a(x) where a is a univariate polynomial of size na.

void ll_upol_multi_eval(R *r,R *a,uma na,R *b,uma nb);

Set r such that r[i] is a(b[i]) where a is a univariate polynomial of size na and b is a vector of size nb.

void ll_upol_eval_horner(R r,R *a,uma na,R x);

Set r to be a(x) where a is a univariate polynomial of size na. The Horner algorithm is used.

void ll_upol_multi_eval_horner(R *r,R *a,uma na,R *b,uma nb);

Set r such that r[i] is a(b[i]) where a is a univariate polynomial of size na and b is a vector of size nb.

void ll_upol_gcd_naive(R *r,uma *nr,R *a,uma na,R *b,uma nb);

Set r of size nr to be GCD of a of size na and b of size nb. Note that the degree of a must be na - 1 and the degree of b must be nb - 1. The schoolbook algorithm is used.

void ll_upol_gcd(R *r,uma *nr,R *a,uma na,R *b,uma nb);

Set r of size nr to be GCD of a of size na and b of size nb. Note that the degree of a must be na - 1 and the degree of b must be nb - 1.

void ll_upol_interpolation_naive(R *r,uma *nr,R *x,R *y,uma n);

Set r of size nr to be the polynomial of degree at most n - 1 which interpolates the n points (x[i],y[i]), in other words such that r(x[i]) = y[i]. The naive algorithm with Lagrange polynomials is used.

void ll_upol_interpolation(R *r,uma *nr,R *x,R *y,uma n);

Set r of size nr to be the polynomial of degree at most n - 1 which interpolates the n points (x[i],y[i]), in other words such that r(x[i]) = y[i].

void ll_upol_mul_naive(R *r,R *a,uma na,R *b,uma nb);

Set r to a * b where a has size na and b has size nb. The shoolbook algorithm is used.

void ll_upol_mul(R *r,R *a,uma na,R *b,uma nb);

Set r to a * b where a has size na and b has size nb. By default the shcoolbook algorithm is used, ll_upol_naive is called. You can provide your own multiplication function with MY_UPOL_MUL, see decoding-defines(3) for more details.

void ll_upol_mul_by_X(R *r,R *a,uma na);

Set r to be a(X) * X where a has size na. You can have r == a.

void ll_upol_mul_by_Xn(R *r,R *a,uma na,uma n);

Set r to be a(X) * X^n where a has size na. You can have r == a.

void ll_upol_mul_by_affine(R *r,R *a,uma na,R b)

Set r to be a * (X + b) where a is a univariate polynomial of size na. You can have r == a.

void ll_upol_normalize(R *r,uma nr);

Normalize r of size nr. note that r must have degree nr - 1.

void ll_upol_partial_eea_naive(R *r,uma *nr,R *u,uma *nu,R *v,uma *nv,R *a,uma na,R *b,uma nb,uma k);

Set r of size nr, u of size nu and v of size nv to be the polynomials involved in the Bezout's identity produced by the extended Euclidean algorithm when the remainder has degree at most k with input a of size na and b of size nb. Namely it returns r, u and v such that u(X) a(X) + v(X) b(X) = r(X) such that $r$ is the first remainder which has degree at most k. Note that u must point to a region of at least nb allocated R's and v must point to a region of at least na allocated R's. The naive partial extended Euclidean algorithm is used.

void ll_upol_partial_eea(R *r,uma *nr,R *u,uma *nu,R *v,uma *nv,R *a,uma na,R *b,uma nb);

Set r of size nr, u of size nu and v of size nv to be the polynomials involved in the Bezout's identity produced by the extended Euclidean algorithm when the remainder has degree at most k with input a of size na and b of size nb. Namely it returns r, u and v such that u(X) a(X) + v(X) b(X) = r(X) such that $r$ is the first remainder which has degree at most k. Note that u must point to a region of at least nb allocated R's and v must point to a region of at least na allocated R's.

void ll_upol_pow_mod(R *r,uma *nr,R *a,uma na,mpz_t e,R *m,uma nm);

Set r to a(X) ^ e mod m where a has size na and m has size nm and set nr to the size of r. Note that e must be a GMP integer. If e is of type uma, see the ll_upol_pow_mod_uma function instead.

void ll_upol_pow_mod_uma(R *r,uma *nr,R *a,uma na,uma e,R *m,uma nm);

Set r to a(X) ^ e mod m where a has size na and m has size nm and set nr to the size of r. Note that e must be a uma integer. If e is a multiprecision integer, see the ll_upol_pow_mod function instead.

void ll_upol_print(FILE *out,R *a,uma na);

Pruma into out the univariate polynomial a of size na.

void ll_upol_prod_affine(R *r,R *a,uma na);

Set r to be the product of the X + a[i] where a is a vector of size na. The dichotomy algorithm is used.

void ll_upol_prod_affine_neg(R *r,R *a,uma na);

Set r to be the product of the X - a[i] where a is a vector of size na. The dichotomy algorithm is used.

void ll_upol_random(R *r,uma nr);

Set r to be a random polynomial of size nr. The resulting polynomial will have degree nr - 1. If you do not require this property prefer ll_vec_random instead.

R *ll_upol_malloc_random(uma nr);

Return a pointer to a random polynomial of size nr. The resulting polynomial will have degree nr - 1. If you do not require this property prefer ll_vec_malloc_random instead.

void ll_upol_roots(R *r,uma *nr,R *a,uma na);

Set r of size nr to be the vector constituted by the roots of a of size na. By default no algorithm is provided and the function will terminate your program. Usually the ring provides its own root finding function, see decoding-alphabets(3) for details. You can provide your own root finding function with MY_UPOL_ROOTS, see decoding-defines(3) for details.

void ll_upol_roots_fq(R *r,uma *nr,R *a,uma na);

Assuming that R is a finite field, set r to the nr roots in R of a of size na. Note that a must have degree na - 1.

void ll_upol_shift_naive(R *r,R *a,uma na,R b);

Set r to a(X + b) where a has size na. The schoolbook algorithm is used.

void ll_upol_shift(R *r,R *a,uma na,R b);

Set r to a(X + b) where a has size na. By default, if the ring has caracteristic 2 (RING_CHAR2) ll_upol_shift_car2 is called, and ll_upol_shift_fast is called otherwise. In both cases the dichotomy algorithm is used. You can also provide your own shifting function, see decoding-defines(3) for details.

void ll_upol_shift_fast2(R *r,R *a,uma na,R b,struct ll_upol_shift_fast_pc *pc);

Set r to a(X + b) where a has size na and pc is a pointer to a C structure where precomputations for the algorithm have been made by ll_upol_shift_fast_prepare. The memory occupied by pc is freed with ll_upol_shift_fast_finish. For example:

struct ll_upol_shift_fast_pc *pc;;
pc = ll_upol_shift_fast_prepare(na,x);
[ ... ]
ll_upol_shift_fast2(r,a,na,x,pc);
[ ... ]
ll_upol_shift_fast_finish(pc);
struct ll_upol_shift_fast_pc *ll_upol_shift_fast_prepare(uma na,R x);

Return a pointer to a structure containing precomputations for the shift of univariate polynomials of size na at x.

void ll_upol_shift_fast_finish(struct ll_upol_shift_fast_pc *pc);

Free the memory occupied by pc where pc was returned by a call to ll_upol_shift_fast_prepare.

void ll_upol_shift_fast_finish(struct ll_upol_shift_fast_pc *pc);

Set r to a(X + x) where a has size nax. The dichotomy algorithm is used. Note that if your ring has caracteristic 2 prefer ll_upol_shift_car22.

void ll_upol_shift_car22(R *r,R *a,uma na,R b,struct ll_upol_shift_car2_pc *pc);

Set r to a(X + b) where a has size na and pc is a pointer to a C structure where precomputations for the algorithm have been made by ll_upol_shift_car2_prepare. The memory occupied by pc is freed with ll_upol_shift_car2_finish. Note that this algorithm is specific for rings of caracterstic 2. A call to ll_upol_shift_car22 with a ring of caracteristic not 2 will return a mathematically false result. The dichotomy algorithm is used. You can have r == a. For example:

struct ll_upol_shift_car2_pc *pc;
pc = ll_upol_shift_car2_prepare(na,x);
[ ... ]
ll_upol_shift_car22(r,a,na,x,pc);
[ ... ]
ll_upol_shift_car2_finish(pc);
struct ll_upol_shift_car2_pc *ll_upol_shift_car2_prepare(uma na,R x);

Return a structure containing precomputations for the shift of univariate polynomials of size nax at x. Note that this function is specific for caracteristic 2.

void ll_upol_shift_car2_finish(struct ll_upol_shift_car2_pc *pc);

Free the memory occupied by pc where pc was returned by a call to ll_upol_shift_fast_prepare. Note that this function is specific for caracteristic 2.

void ll_upol_shift_car2(R *r,R *a,uma na,R x);

Set r to a(X + x) where a has size na. Note that this function is specific for caracteristic 2. A call to ll_upol_shift_car2 with a ring of caracteristic not 2 will return a mathematically false result. The dichotomy algorithm is used.

void ll_upol_square_naive(R *r,R *a,uma na);

Set r to a(X)^2 where a has size na. The schoolbook algorithm is used.

void ll_upol_square(R *r,R *a,uma na);

Set r to a(X)^2 where a has size na. The schoolbook algorithm is used, ll_upol_square_naive is called. You can provide your own squaring function with MY_UPOL_SQUARE, see decoding-defines(3) for more details.

DENSE BIVARIATE POLYNOMIALS

struct ll_vec_lili *ll_dbpol_Xroots(R *a,uma nax,uma nay,uma n);

Return a linked list containing the polynomials f such that a(f(X),X) = 0 mod X^n where a has size (nax,nay). The algorithm is a mix between Rock and Ruckenstein's and Berthomieu, Lecerf and Quintin's.

void ll_dbpol_discrepancy(R r,R *a,uma nax,uma nay,uma mx,uma my,R x,R y);

Set r to be the evaluation at (x,y) of the (mx,my)-Hasse derivative of a of size (nax,nay).

void ll_dbpol_eval(R z,R *a,uma nax,uma nay,R x,R y);

Evaluate a of size (nax,nay) at (x,y).

void ll_dbpol_horner_upol(R *r,uma *nr,R *a,uma nax,uma nay,R *b,uma nb);

Evaluate a of size (nax,nay) at b of size nb.

void ll_dbpol_eval_horner(R z,R *a,uma nax,uma nay,R x,R y);

Evaluate a of size (nax,nay) at (x,y). The Horner algorihtm is used.

void ll_dbpol_eval_horner_upol(R *r,uma *nr,R *a,uma nax,uma nay,R *b,uma nb);

Evaluate a of size (nax,nay) at b of size nb. The Horner algorithm is used.

void ll_dbpol_mul_naive(R *r,R *a,uma nax,uma nay,R *b,uma nbx,uma nby);

Set r to a * b where a and b have size (nax,nay) and (nbx,nby). The schoolbook algorithm is used.

void ll_dbpol_mul(R *r,R *a,uma nax,uma nay,R *b,uma nbx,uma nby);

Set r to a * b where a and b have size (nax,nay) and (nbx,nby). By default ll_dbpol_mul_naive is called. You can provide a better algorithm with MY_DBPOL_MUL. See decoding-defines(3) for more details. Some rings can provide their own function. See decoding-alphabets(3) for more details.

void ll_dbpol_mul_by_Y(R *r,R *a,uma nax,uma nay);

Set r to a * Y where a has size (nax,nay). You can have r == a.

void ll_dbpol_mul_by_Yn(R *r,R *a,uma nax,uma nay,uma n);

Set r to a * Y^n where a has size (nax,nay). You can have r == a.

void ll_dbpol_mul_by_affine(R *r,R *a,uma nax,uma nay,R b);

Set r to a * (Y + b) where a has size (nax,nay). You can have r == a.

void ll_dbpol_mul_by_upolX(R *r,R *a,uma nax,uma nay,R *b,uma nbx);

Set r to a * b where a has size (nax,nay).

void ll_dbpol_mul_by_upolY(R *r,R *a,uma nax,uma nay,R *b,uma nby)

Set r to a(X,Y) * b(Y) where a has size (nax,nay) and b must be represented as univariate polynomial of size nby.

void ll_dbpol_print_factorized(FILE *out,R *a,uma nax,uma nay);

Pruma to out a of size (nax,nay). It prints the bivariate polynomial according to the powers of Y: f0(X) + f1(X).Y + ... + fn(X).Y^n.

void ll_dbpol_replace_X_by_XYn(R *r,uma *nry,R *a,uma nax,uma nay,uma n);

Set r to a(X,XY^n) where a has size (nax,nay). You can have r == a.

void ll_dbpol_shift_naive(R *r,R *a,uma nax,uma nay,R x,R y);

Set r to a(X + x,Y + y) where a has size (nax,nay). The horner algorithm is used.

void ll_dbpol_shift(R *r,R *a,uma nax,uma nay,R x,R y);

Set r to a(X + x,Y + y) where a has size (nax,nay). You can provide your own function with MY_DBPOL_SHIFT. If your ring has caracteristic 2 (RING_COMM_CHAR2) then ll_dbpol_shift_car2 is called and ll_dbpol_shift_fast otherwise.

void ll_dbpol_shift_fast2(R *r,R *a,uma nax,uma nay,R x,R y,struct ll_dbpol_shift_fast_pc *pc);

Set r to a(X + x,Y + y) where a has size (nax,nay) and pc is a pointer to a C structure where precomputations for the algorithm have been made by ll_dbpol_shift_fast_prepare. The memory occupied by pc is freed with ll_dbpol_shift_fast_finish. For example:

struct ll_dbpol_shift_fast_pc *pc;
pc = ll_dbpol_shift_fast_prepare(nax,nay,x,y);
[ ... ]
ll_dbpol_shift_fast2(r,a,nax,nay,x,y,pc);
[ ... ]
ll_dbpol_shift_fast_finish(pc);
struct ll_dbpol_shift_fast_pc *ll_dbpol_shift_fast_prepare(uma nax,uma nay,R x,R y);

Return a pointer to a structure containing precomputations for the shift of dense bivariate polynomials of size (nax,nay) at (x,y).

void ll_dbpol_shift_fast_finish(struct ll_dbpol_shift_fast_pc *pc);

Free the memory occupied by pc where pc was returned by a call to ll_dbpol_shift_fast_prepare.

void ll_dbpol_shift_fast(R *r,R *a,uma nax,uma nay,R x,R y);

Set r to a(X + x,Y + y) where a has size (nax,nay). The dichotomy algorithm is used. Note that if your ring has caracteristic 2 prefer ll_dbpol_shift_car22.

void ll_dbpol_shift_car22(R *r,R *a,uma nax,uma nay,R x,R y,struct ll_dbpol_shift_car2_pc *pc);

Set r to a(X + x,Y + y) where a has size (nax,nay) and pc is a pointer to a C structure where precomputations for the algorithm have been made by ll_dbpol_shift_car2_prepare. The memory occupied by pc is freed with ll_dbpol_shift_car2_finish. Note that this algorithm is specific for rings of caracterstic 2. A call to ll_dbpol_shift_car22 with a ring of caracteristic not 2 will return a mathematically false result. The dichotomy algorithm is used. For example:

struct ll_dbpol_shift_car2_pc *pc;
pc = ll_dbpol_shift_car2_prepare(nax,nay,x,y);
[ ... ]
ll_dbpol_shift_car22(r,a,nax,nay,x,y,pc);
[ ... ]
ll_dbpol_shift_car2_finish(pc);
struct ll_dbpol_shift_car2_pc *ll_dbpol_shift_car2_prepare(uma nax,uma nay,R x,R y);

Return a structure containing precomputations for the shift of dense bivariate polynomials of size (nax,nay) at (x,y). Note that this function is specific for caracteristic 2.

void ll_dbpol_shift_car2(R *r,R *a,uma nax,uma nay,R x,R y);

Set r to a(X + x,Y + y) where a has size (nax,nay). Note that this function is specific for caracteristic 2. A call to ll_dbpol_shift_car2 with a ring of caracteristic not 2 will return a mathematically false result. The dichotomy algorithm is used.

void ll_dbpol_shift_car2_finish(struct ll_dbpol_shift_car2_pc *pc)

Free the memory occupied by pc where pc was returned by a call to ll_dbpol_shift_fast_prepare. Note that this function is specific for caracteristic 2.

void ll_dbpol_shift_in_X(R *r,R *a,uma nax,uma nay,R x);

Set r to a(X + x,Y) where a has size (nax,nay). If your ring has caracteristic 2 (RING_COMM_CHAR2) then ll_upol_shift_car22 is called, otherwise ll_upol_shift_fast2 is called. The dichotomy algorithm is used.

void ll_dbpol_upol_map_XtoY(R *r,uma nrx,R *a,uma na);

Set r to a(Y) where a is a univariate polynomial of size na. The size or r will be (nrx,na).

INTERPOLATION ALGORITHMS

uma ll_koetter_single(R *q,uma *nqx,uma *nqy,R *supp,uma n,uma k,uma L,uma m,R *y,uma tau);

Set q of size (nqx,nqy) to be the interpolation polynomial of the Guruswami-Sudan list decoding algorithm for a Reed-Solomon code of support supp, length n, dimension k, list-size L and for a received word y with at most tau errors. The multiplicity m is the same for all points (supp[i],y[i]). If q is found 0 is returned otherwise -1 is returned. The Koetter algorithm is used.

uma ll_koetter_multi(R *q,uma *nqx,uma *nqy,R *supp,uma n,uma k,uma L,uma *m,R *y,uma tau);

Set q of size (nqx,nqy) to be the interpolation polynomial of the Guruswami-Sudan list decoding algorithm for a Reed-Solomon code of support supp, length n, dimension k, list-size L and for a received word y with at most tau errors. The multiplicities for each pouma (supp[i],y[i]) is given by m[i]. If q is found 0 is returned otherwise -1 is returned. The Koetter algorithm is used.

LINKED LIST OF VECTORS

struct ll_vec_lili *ll_vec_lili_node_malloc(uma na);

Return a ll_vec_lili struct containing an allocated vector of size na.

struct ll_vec_lili *ll_vec_lili_node_malloc_copy(R *a,uma na);

Return a ll_vec_lili struct containing a deep copy of the vector a of size na.

void ll_vec_lili_node_free(struct ll_vec_lili *node);

Free the ll_vec_lili struct node.

void ll_vec_lili_node_attach_after(struct ll_vec_lili *after,struct ll_vec_lili *node);

Attach the node after after.

void ll_vec_lili_node_attach_before(struct ll_vec_lili *before,struct ll_vec_lili *node);

Attach node before before.

void ll_vec_lili_node_detach(struct ll_vec_lili *node);

Detach node.

void ll_vec_lili_list_free(struct ll_vec_lili *start);

Free the linked list whose first node is start.

int ll_vec_lili_list_in(struct ll_vec_lili *start,R *a,uma na);

Return nonzero if and only if a of size na is within the linked list pointed to by start.

MISC FUNCTIONS

void ll_die(const char *msg);

Pruma msg to stderr and terminate the program with exit code -1.

void ll_print(FILE *out,const char *s,...);

Produce output into out according to the format s. The format is (for the moment) much simpler than printf's.

%e

Print an element of the ring.

%v

Print a vector, its size must be passed as an argument.

%m

Print a matrix, its size must be passed as an argument.

%p

Print a univariate polynomial, its size must be passed as an argument.

%q

Print a dense bivariate polynomial, its size must passed as an argument.

%d, %c, %s, %f

Print an int, a char, a string (char*) and a double.

For example:

R *a;
int na,nax,nay,ra,ca;
int n;
char b,*str;
double d;

ll_print(stdout,"A vector %v\n",a,na);
ll_print(stdout,"A matrix %m\n",a,ra,ca);
ll_print(stdout,"A univariate polynomial %p\n",a,na);
ll_print(stdout,"A dense bivariate polynomial %q\n",a,nax,nay);
ll_print(stdout,"And the usual integer %d,"
                "char %c, string %s, double %f\n",n,b,str,d);

AUTHOR

Written by Guillaume Quintin (quintin@lix.polytechnique.fr).

SEE ALSO

decoding-defines(3), decoding-alphabets(3), malloc(3), abort(3), signal(7), sigaction(2), raise(3), printf(3)