/** * Computes the overall period of the variables I for (MI) mod |d|, where M is * a matrix and |d| a vector. Produce a diagonal matrix S = (s_k) where s_k is * the overall period of i_k * @param M the set of affine functions of I (row-vectors) * @param d the column-vector representing the modulos */ Matrix * affine_periods(Matrix * M, Matrix * d) { Matrix * S; unsigned int i,j; Value tmp; Value * periods = (Value *)malloc(sizeof(Value) * M->NbColumns); value_init(tmp); for(i=0; i< M->NbColumns; i++) { value_init(periods[i]); value_set_si(periods[i], 1); } for (i=0; i<M->NbRows; i++) { for (j=0; j< M->NbColumns; j++) { value_gcd(tmp, d->p[i][0], M->p[i][j]); value_divexact(tmp, d->p[i][0], tmp); value_lcm(periods[j], periods[j], tmp); } } value_clear(tmp); /* 2- build S */ S = Matrix_Alloc(M->NbColumns, M->NbColumns); for (i=0; i< M->NbColumns; i++) for (j=0; j< M->NbColumns; j++) if (i==j) value_assign(S->p[i][j],periods[j]); else value_set_si(S->p[i][j], 0); /* 3- clean up */ for(i=0; i< M->NbColumns; i++) value_clear(periods[i]); free(periods); return S; } /* affine_periods */
static Polyhedron *partition2polyhedron(Matrix *A, struct barvinok_options *options) { int i; unsigned nvar, nparam; Matrix *M; Polyhedron *P; nvar = A->NbColumns; nparam = A->NbRows; M = Matrix_Alloc(nvar + nparam, 1 + nvar + nparam + 1); assert(M); for (i = 0; i < nparam; ++i) { Vector_Copy(A->p[i], M->p[i] + 1, nvar); value_set_si(M->p[i][1 + nvar + i], -1); } for (i = 0; i < nvar; ++i) { value_set_si(M->p[nparam + i][0], 1); value_set_si(M->p[nparam + i][1 + i], 1); } P = Constraints2Polyhedron(M, options->MaxRays); Matrix_Free(M); return P; }
static Matrix *Polyhedron2standard_form(Polyhedron *P, Matrix **T) { int i, j; int rows; unsigned dim = P->Dimension; Matrix *M2; Matrix *H, *U; Matrix M; assert(P->NbEq == 0); Polyhedron_Remove_Positivity_Constraint(P); for (i = 0; i < P->NbConstraints; ++i) assert(value_zero_p(P->Constraint[i][1+dim])); Polyhedron_Matrix_View(P, &M, P->NbConstraints); H = standard_constraints(&M, 0, &rows, &U); *T = homogenize(U); Matrix_Free(U); M2 = Matrix_Alloc(rows, 2+dim+rows); for (i = dim; i < H->NbRows; ++i) { Vector_Copy(H->p[i], M2->p[i-dim]+1, dim); value_set_si(M2->p[i-dim][1+i], -1); } for (i = 0, j = H->NbRows-dim; i < dim; ++i) { if (First_Non_Zero(H->p[i], i) == -1) continue; Vector_Oppose(H->p[i], M2->p[j]+1, dim); value_set_si(M2->p[j][1+j+dim], 1); ++j; } Matrix_Free(H); return M2; }
int Vector_IsZero(Value * v, unsigned length) { unsigned i; if (value_notzero_p(v[0])) return 0; else { value_set_si(v[0], 1); for (i=length-1; value_zero_p(v[i]); i--); value_set_si(v[0], 0); return (i==0); } }
/** * Computes the validity lattice of a set of equalities. I.e., the lattice * induced on the last <tt>b</tt> variables by the equalities involving the * first <tt>a</tt> integer existential variables. The submatrix of Eqs that * concerns only the existential variables (so the first a columns) is assumed * to be full-row rank. * @param Eqs the equalities * @param a the number of existential integer variables, placed as first * variables * @param vl the (returned) validity lattice, in homogeneous form. It is * allocated if initially set to null, or reused if already allocated. */ void Equalities_validityLattice(Matrix * Eqs, int a, Matrix** vl) { unsigned int b = Eqs->NbColumns-2-a; unsigned int r = Eqs->NbRows; Matrix * A=NULL, * B=NULL, *I = NULL, *Lb=NULL, *sol=NULL; Matrix *H, *U, *Q; unsigned int i; if (dbgCompParm) { printf("Computing validity lattice induced by the %d first variables of:" ,a); show_matrix(Eqs); } if (b==0) { ensureMatrix((*vl), 1, 1); value_set_si((*vl)->p[0][0], 1); return; } /* 1- check that there is an integer solution to the equalities */ /* OPT: could change integerSolution's profile to allocate or not*/ Equalities_integerSolution(Eqs, &sol); /* if there is no integer solution, there is no validity lattice */ if (sol==NULL) { if ((*vl)!=NULL) Matrix_Free(*vl); return; } Matrix_subMatrix(Eqs, 0, 1, r, 1+a, &A); Matrix_subMatrix(Eqs, 0, 1+a, r, 1+a+b, &B); linearInter(A, B, &I, &Lb); Matrix_Free(A); Matrix_Free(B); Matrix_Free(I); if (dbgCompParm) { show_matrix(Lb); } /* 2- The linear part of the validity lattice is the left HNF of Lb */ left_hermite(Lb, &H, &Q, &U); Matrix_Free(Lb); Matrix_Free(Q); Matrix_Free(U); /* 3- build the validity lattice */ ensureMatrix((*vl), b+1, b+1); Matrix_copySubMatrix(H, 0, 0, b, b, (*vl), 0,0); Matrix_Free(H); for (i=0; i< b; i++) { value_assign((*vl)->p[i][b], sol->p[0][a+i]); } Matrix_Free(sol); Vector_Set((*vl)->p[b],0, b); value_set_si((*vl)->p[b][b], 1); } /* validityLattice */
static Matrix *VectorArray2Matrix(VectorArray array, unsigned cols) { int i, j; Matrix *M = Matrix_Alloc(array->Size, cols+1); for (i = 0; i < array->Size; ++i) { for (j = 0; j < cols; ++j) value_set_si(M->p[i][j], array->Data[i][j]); value_set_si(M->p[i][cols], 1); } return M; }
/** * Computes the intersection of two linear lattices, whose base vectors are * respectively represented in A and B. * If I and/or Lb is set to NULL, then the matrix is allocated. * Else, the matrix is assumed to be allocated already. * I and Lb are rk x rk, where rk is the rank of A (or B). * @param A the full-row rank matrix whose column-vectors are the basis for the * first linear lattice. * @param B the matrix whose column-vectors are the basis for the second linear * lattice. * @param Lb the matrix such that B.Lb = I, where I is the intersection. * @return their intersection. */ static void linearInter(Matrix * A, Matrix * B, Matrix ** I, Matrix **Lb) { Matrix * AB=NULL; int rk = A->NbRows; int a = A->NbColumns; int b = B->NbColumns; int i,j, z=0; Matrix * H, *U, *Q; /* ensure that the spanning vectors are in the same space */ assert(B->NbRows==rk); /* 1- build the matrix * (A 0 1) * (0 B 1) */ AB = Matrix_Alloc(2*rk, a+b+rk); Matrix_copySubMatrix(A, 0, 0, rk, a, AB, 0, 0); Matrix_copySubMatrix(B, 0, 0, rk, b, AB, rk, a); for (i=0; i< rk; i++) { value_set_si(AB->p[i][a+b+i], 1); value_set_si(AB->p[i+rk][a+b+i], 1); } if (dbgCompParm) { show_matrix(AB); } /* 2- Compute its left Hermite normal form. AB.U = [H 0] */ left_hermite(AB, &H, &Q, &U); Matrix_Free(AB); Matrix_Free(Q); /* count the number of non-zero colums in H */ for (z=H->NbColumns-1; value_zero_p(H->p[H->NbRows-1][z]); z--); z++; if (dbgCompParm) { show_matrix(H); printf("z=%d\n", z); } Matrix_Free(H); /* if you split U in 9 submatrices, you have: * A.U_13 = -U_33 * B.U_23 = -U_33, * where the nb of cols of U_{*3} equals the nb of zero-cols of H * U_33 is a (the smallest) combination of col-vectors of A and B at the same * time: their intersection. */ Matrix_subMatrix(U, a+b, z, U->NbColumns, U->NbColumns, I); Matrix_subMatrix(U, a, z, a+b, U->NbColumns, Lb); if (dbgCompParm) { show_matrix(U); } Matrix_Free(U); } /* linearInter */
/*---------------------------------------------------------------------*/ static int exist_points(int pos,Polyhedron *Pol,Value *context) { Value LB, UB, k,tmp; value_init(LB); value_init(UB); value_init(k); value_init(tmp); value_set_si(LB,0); value_set_si(UB,0); /* Problem if UB or LB is INFINITY */ if (lower_upper_bounds(pos,Pol,context,&LB,&UB) !=0) { errormsg1("exist_points", "infdom", "infinite domain"); value_clear(LB); value_clear(UB); value_clear(k); value_clear(tmp); return -1; } value_set_si(context[pos],0); if(value_lt(UB,LB)) { value_clear(LB); value_clear(UB); value_clear(k); value_clear(tmp); return 0; } if (!Pol->next) { value_subtract(tmp,UB,LB); value_increment(tmp,tmp); value_clear(UB); value_clear(LB); value_clear(k); return (value_pos_p(tmp)); } for (value_assign(k,LB);value_le(k,UB);value_increment(k,k)) { /* insert k in context */ value_assign(context[pos],k); if (exist_points(pos+1,Pol->next,context) > 0 ) { value_clear(LB); value_clear(UB); value_clear(k); value_clear(tmp); return 1; } } /* Reset context */ value_set_si(context[pos],0); value_clear(UB); value_clear(LB); value_clear(k); value_clear(tmp); return 0; }
/* * Compute n! */ void Factorial(int n, Value *fact) { int i; Value tmp; value_init(tmp); value_set_si(*fact,1); for (i=1;i<=n;i++) { value_set_si(tmp,i); value_multiply(*fact,*fact,tmp); } value_clear(tmp); } /* Factorial */
void insert_constraint_into_matrix (CloogMatrix *m, int row, ppl_const_Constraint_t cstr) { ppl_Coefficient_t c; ppl_dimension_type i, dim, nb_cols = m->NbColumns; ppl_Constraint_space_dimension (cstr, &dim); ppl_new_Coefficient (&c); for (i = 0; i < dim; i++) { ppl_Constraint_coefficient (cstr, i, c); ppl_Coefficient_to_mpz_t (c, m->p[row][i + 1]); } for (i = dim; i < nb_cols - 1; i++) value_set_si (m->p[row][i + 1], 0); ppl_Constraint_inhomogeneous_term (cstr, c); ppl_Coefficient_to_mpz_t (c, m->p[row][nb_cols - 1]); value_set_si (m->p[row][0], 1); switch (ppl_Constraint_type (cstr)) { case PPL_CONSTRAINT_TYPE_LESS_THAN: oppose_constraint (m, row); case PPL_CONSTRAINT_TYPE_GREATER_THAN: value_sub_int (m->p[row][nb_cols - 1], m->p[row][nb_cols - 1], 1); break; case PPL_CONSTRAINT_TYPE_LESS_OR_EQUAL: oppose_constraint (m, row); case PPL_CONSTRAINT_TYPE_GREATER_OR_EQUAL: break; case PPL_CONSTRAINT_TYPE_EQUAL: value_set_si (m->p[row][0], 0); break; default: /* Not yet implemented. */ gcc_unreachable(); } ppl_delete_Coefficient (c); }
/* * Given matrices 'Mat1' and 'Mat2', compute the matrix product and store in * matrix 'Mat3' */ void Matrix_Product(Matrix *Mat1,Matrix *Mat2,Matrix *Mat3) { int Size, i, j, k; unsigned NbRows, NbColumns; Value **q1, **q2, *p1, *p3,sum; NbRows = Mat1->NbRows; NbColumns = Mat2->NbColumns; Size = Mat1->NbColumns; if(Mat2->NbRows!=Size||Mat3->NbRows!=NbRows||Mat3->NbColumns!=NbColumns) { fprintf(stderr, "? Matrix_Product : incompatable matrix dimension\n"); return; } value_init(sum); p3 = Mat3->p_Init; q1 = Mat1->p; q2 = Mat2->p; /* Mat3[i][j] = Sum(Mat1[i][k]*Mat2[k][j] where sum is over k = 1..nbrows */ for (i=0;i<NbRows;i++) { for (j=0;j<NbColumns;j++) { p1 = *(q1+i); value_set_si(sum,0); for (k=0;k<Size;k++) { value_addmul(sum, *p1, *(*(q2+k)+j)); p1++; } value_assign(*p3,sum); p3++; } } value_clear(sum); return; } /* Matrix_Product */
static double compute_enode(enode *p, Value *list_args) { int i; Value m, param; double res=0.0; if (!p) return(0.); value_init(m); value_init(param); if (p->type == polynomial) { if (p->size > 1) value_assign(param,list_args[p->pos-1]); /* Compute the polynomial using Horner's rule */ for (i=p->size-1;i>0;i--) { res +=compute_evalue(&p->arr[i],list_args); res *=VALUE_TO_DOUBLE(param); } res +=compute_evalue(&p->arr[0],list_args); } else if (p->type == periodic) { value_assign(m,list_args[p->pos-1]); /* Choose the right element of the periodic */ value_set_si(param,p->size); value_pmodulus(m,m,param); res = compute_evalue(&p->arr[VALUE_TO_INT(m)],list_args); } value_clear(m); value_clear(param); return res; } /* compute_enode */
/* * Return the component of 'p' with minimum non-zero absolute value. 'index' * points to the component index that has the minimum value. If no such value * and index is found, Value 1 is returned. */ void Vector_Min_Not_Zero(Value *p,unsigned length,int *index,Value *min) { Value aux; int i; i = First_Non_Zero(p, length); if (i == -1) { value_set_si(*min,1); return; } *index = i; value_absolute(*min, p[i]); value_init(aux); for (i = i+1; i < length; i++) { if (value_zero_p(p[i])) continue; value_absolute(aux, p[i]); if (value_lt(aux,*min)) { value_assign(*min,aux); *index = i; } } value_clear(aux); } /* Vector_Min_Not_Zero */
/*--------------------------------------------------------------*/ int Polyhedron_Not_Empty(Polyhedron *P,Polyhedron *C,int MAXRAYS) { int res,i; Value *context; Polyhedron *L; POL_ENSURE_FACETS(P); POL_ENSURE_VERTICES(P); POL_ENSURE_FACETS(C); POL_ENSURE_VERTICES(C); /* Create a context vector size dim+2 and set it to all zeros */ context = (Value *) malloc((P->Dimension+2)*sizeof(Value)); /* Initialize array 'context' */ for (i=0;i<(P->Dimension+2);i++) value_init(context[i]); Vector_Set(context,0,(P->Dimension+2)); /* Set context[P->Dimension+1] = 1 (the constant) */ value_set_si(context[P->Dimension+1],1); L = Polyhedron_Scan(P,C,MAXRAYS); res = exist_points(1,L,context); Domain_Free(L); /* Clear array 'context' */ for (i=0;i<(P->Dimension+2);i++) value_clear(context[i]); free(context); return res; }
static ppl_Pointset_Powerset_C_Polyhedron_t dr_equality_constraints (graphite_dim_t dim, graphite_dim_t pos, graphite_dim_t nb_subscripts) { ppl_Polyhedron_t subscript_equalities; ppl_Pointset_Powerset_C_Polyhedron_t res; Value v, v_op; graphite_dim_t i; value_init (v); value_init (v_op); value_set_si (v, 1); value_set_si (v_op, -1); ppl_new_C_Polyhedron_from_space_dimension (&subscript_equalities, dim, 0); for (i = 0; i < nb_subscripts; i++) { ppl_Linear_Expression_t expr; ppl_Constraint_t cstr; ppl_Coefficient_t coef; ppl_new_Coefficient (&coef); ppl_new_Linear_Expression_with_dimension (&expr, dim); ppl_assign_Coefficient_from_mpz_t (coef, v); ppl_Linear_Expression_add_to_coefficient (expr, pos + i, coef); ppl_assign_Coefficient_from_mpz_t (coef, v_op); ppl_Linear_Expression_add_to_coefficient (expr, pos + i + nb_subscripts, coef); ppl_new_Constraint (&cstr, expr, PPL_CONSTRAINT_TYPE_EQUAL); ppl_Polyhedron_add_constraint (subscript_equalities, cstr); ppl_delete_Linear_Expression (expr); ppl_delete_Constraint (cstr); ppl_delete_Coefficient (coef); } ppl_new_Pointset_Powerset_C_Polyhedron_from_C_Polyhedron (&res, subscript_equalities); value_clear (v); value_clear (v_op); ppl_delete_Polyhedron (subscript_equalities); return res; }
/* Return * T 0 * 0 1 */ static Matrix *homogenize(Matrix *T) { int i; Matrix *H = Matrix_Alloc(T->NbRows+1, T->NbColumns+1); for (i = 0; i < T->NbRows; ++i) Vector_Copy(T->p[i], H->p[i], T->NbColumns); value_set_si(H->p[T->NbRows][T->NbColumns], 1); return H; }
/* * Return the inner product of the two Vectors 'p1' and 'p2' */ void Inner_Product(Value *p1, Value *p2, unsigned length, Value *ip) { int i; if (length != 0) value_multiply(*ip, p1[0], p2[0]); else value_set_si(*ip, 0); for(i = 1; i < length; i++) value_addmul(*ip, p1[i], p2[i]); } /* Inner_Product */
/* * Assign 'n' to each component of Vector 'p' */ void Vector_Set(Value *p,int n,unsigned length) { Value *cp; int i; cp = p; for (i=0;i<length;i++) { value_set_si(*cp,n); cp++; } return; } /* Vector_Set */
/** * Given a matrix that defines a full-dimensional affine lattice, returns the * affine sub-lattice spanned in the k first dimensions. * Useful for instance when you only look for the parameters' validity lattice. * @param lat the original full-dimensional lattice * @param subLat the sublattice */ void Lattice_extractSubLattice(Matrix * lat, unsigned int k, Matrix ** subLat) { Matrix * H, *Q, *U, *linLat = NULL; unsigned int i; dbgStart(Lattice_extractSubLattice); /* if the dimension is already good, just copy the initial lattice */ if (k==lat->NbRows-1) { if (*subLat==NULL) { (*subLat) = Matrix_Copy(lat); } else { Matrix_copySubMatrix(lat, 0, 0, lat->NbRows, lat->NbColumns, (*subLat), 0, 0); } return; } assert(k<lat->NbRows-1); /* 1- Make the linear part of the lattice triangular to eliminate terms from other dimensions */ Matrix_subMatrix(lat, 0, 0, lat->NbRows, lat->NbColumns-1, &linLat); /* OPT: any integer column-vector elimination is ok indeed. */ /* OPT: could test if the lattice is already in triangular form. */ left_hermite(linLat, &H, &Q, &U); if (dbgCompParmMore) { show_matrix(H); } Matrix_Free(Q); Matrix_Free(U); Matrix_Free(linLat); /* if not allocated yet, allocate it */ if (*subLat==NULL) { (*subLat) = Matrix_Alloc(k+1, k+1); } Matrix_copySubMatrix(H, 0, 0, k, k, (*subLat), 0, 0); Matrix_Free(H); Matrix_copySubMatrix(lat, 0, lat->NbColumns-1, k, 1, (*subLat), 0, k); for (i=0; i<k; i++) { value_set_si((*subLat)->p[k][i], 0); } value_set_si((*subLat)->p[k][k], 1); dbgEnd(Lattice_extractSubLattice); } /* Lattice_extractSubLattice */
Value *compute_poly(Enumeration *en,Value *list_args) { Value *tmp; /* double d; int i; */ tmp = (Value *) malloc (sizeof(Value)); assert(tmp != NULL); value_init(*tmp); value_set_si(*tmp,0); if(!en) return(tmp); /* no ehrhart polynomial */ if(en->ValidityDomain) { if(!en->ValidityDomain->Dimension) { /* no parameters */ value_set_double(*tmp,compute_evalue(&en->EP,list_args)+.25); return(tmp); } } else return(tmp); /* no Validity Domain */ while(en) { if(in_domain(en->ValidityDomain,list_args)) { #ifdef EVAL_EHRHART_DEBUG Print_Domain(stdout,en->ValidityDomain,NULL); print_evalue(stdout,&en->EP,NULL); #endif /* d = compute_evalue(&en->EP,list_args); i = d; printf("(double)%lf = %d\n", d, i ); */ value_set_double(*tmp,compute_evalue(&en->EP,list_args)+.25); return(tmp); } else en=en->next; } value_set_si(*tmp,0); return(tmp); /* no compatible domain with the arguments */ } /* compute_poly */
static ppl_Constraint_t build_pairwise_constraint (graphite_dim_t dim, graphite_dim_t pos1, graphite_dim_t pos2, int c, enum ppl_enum_Constraint_Type cstr_type) { ppl_Linear_Expression_t expr; ppl_Constraint_t cstr; ppl_Coefficient_t coef; Value v, v_op, v_c; value_init (v); value_init (v_op); value_init (v_c); value_set_si (v, 1); value_set_si (v_op, -1); value_set_si (v_c, c); ppl_new_Coefficient (&coef); ppl_new_Linear_Expression_with_dimension (&expr, dim); ppl_assign_Coefficient_from_mpz_t (coef, v); ppl_Linear_Expression_add_to_coefficient (expr, pos1, coef); ppl_assign_Coefficient_from_mpz_t (coef, v_op); ppl_Linear_Expression_add_to_coefficient (expr, pos2, coef); ppl_assign_Coefficient_from_mpz_t (coef, v_c); ppl_Linear_Expression_add_to_inhomogeneous (expr, coef); ppl_new_Constraint (&cstr, expr, cstr_type); ppl_delete_Linear_Expression (expr); ppl_delete_Coefficient (coef); value_clear (v); value_clear (v_op); value_clear (v_c); return cstr; }
/* * Return the number of ways to choose 'b' items from 'a' items, that is, * return a!/(b!(a-b)!) */ void CNP(int a,int b, Value *result) { int i; Value tmp; value_init(tmp); value_set_si(*result,1); /* If number of items is less than the number to be choosen, return 1 */ if(a <= b) { value_clear(tmp); return; } for(i=a;i>b;--i) { value_set_si(tmp,i); value_multiply(*result,*result,tmp); } for(i=1;i<=(a-b);++i) { value_set_si(tmp,i); value_division(*result,*result,tmp); } value_clear(tmp); } /* CNP */
/* * Compute n!/(p!(n-p)!) */ void Binomial(int n, int p, Value *result) { int i; Value tmp; Value f; value_init(tmp);value_init(f); if (n-p<p) p=n-p; if (p!=0) { value_set_si(*result,(n-p+1)); for (i=n-p+2;i<=n;i++) { value_set_si(tmp,i); value_multiply(*result,*result,tmp); } Factorial(p,&f); value_division(*result,*result,f); } else value_set_si(*result,1); value_clear(f);value_clear(tmp); } /* Binomial */
static void extend_scattering (poly_bb_p pbb, int max_scattering) { ppl_dimension_type nb_old_dims, nb_new_dims; int nb_added_dims, i; ppl_Coefficient_t coef; Value one; nb_added_dims = max_scattering - pbb_nb_scattering_transform (pbb); value_init (one); value_set_si (one, 1); ppl_new_Coefficient (&coef); ppl_assign_Coefficient_from_mpz_t (coef, one); gcc_assert (nb_added_dims >= 0); nb_old_dims = pbb_nb_scattering_transform (pbb) + pbb_dim_iter_domain (pbb) + scop_nb_params (PBB_SCOP (pbb)); nb_new_dims = nb_old_dims + nb_added_dims; ppl_insert_dimensions (PBB_TRANSFORMED_SCATTERING (pbb), pbb_nb_scattering_transform (pbb), nb_added_dims); PBB_NB_SCATTERING_TRANSFORM (pbb) += nb_added_dims; /* Add identity matrix for the added dimensions. */ for (i = max_scattering - nb_added_dims; i < max_scattering; i++) { ppl_Constraint_t cstr; ppl_Linear_Expression_t expr; ppl_new_Linear_Expression_with_dimension (&expr, nb_new_dims); ppl_Linear_Expression_add_to_coefficient (expr, i, coef); ppl_new_Constraint (&cstr, expr, PPL_CONSTRAINT_TYPE_EQUAL); ppl_Polyhedron_add_constraint (PBB_TRANSFORMED_SCATTERING (pbb), cstr); ppl_delete_Constraint (cstr); ppl_delete_Linear_Expression (expr); } ppl_delete_Coefficient (coef); value_clear (one); }
/** * Given a matrix with m parameterized equations, compress the nb_parms * parameters and n-m variables so that m variables are integer, and transform * the variable space into a n-m space by eliminating the m variables (using * the equalities) the variables to be eliminated are chosen automatically by * the function. * <b>Deprecated.</b> Try to use Constraints_fullDimensionize instead. * @param M the constraints * @param the number of parameters * @param validityLattice the the integer lattice underlying the integer * solutions. */ Matrix * full_dimensionize(Matrix const * M, int nbParms, Matrix ** validityLattice) { Matrix * Eqs, * Ineqs; Matrix * permutedEqs, * permutedIneqs; Matrix * Full_Dim; Matrix * WVL; /* The Whole Validity Lattice (vars+parms) */ unsigned int i,j; int nbElimVars; unsigned int * permutation, * permutationInv; /* 0- Split the equalities and inequalities from each other */ split_constraints(M, &Eqs, &Ineqs); /* 1- if the polyhedron is already full-dimensional, return it */ if (Eqs->NbRows==0) { Matrix_Free(Eqs); (*validityLattice) = Identity_Matrix(nbParms+1); return Ineqs; } nbElimVars = Eqs->NbRows; /* 2- put the vars to be eliminated at the first positions, and compress the other vars/parms -> [ variables to eliminate / parameters / variables to keep ] */ permutation = find_a_permutation(Eqs, nbParms); if (dbgCompParm) { printf("Permuting the vars/parms this way: [ "); for (i=0; i< Eqs->NbColumns; i++) { printf("%d ", permutation[i]); } printf("]\n"); } permutedEqs = mpolyhedron_permute(Eqs, permutation); WVL = compress_parms(permutedEqs, Eqs->NbColumns-2-Eqs->NbRows); if (dbgCompParm) { printf("Whole validity lattice: "); show_matrix(WVL); } mpolyhedron_compress_last_vars(permutedEqs, WVL); permutedIneqs = mpolyhedron_permute(Ineqs, permutation); if (dbgCompParm) { show_matrix(permutedEqs); } Matrix_Free(Eqs); Matrix_Free(Ineqs); mpolyhedron_compress_last_vars(permutedIneqs, WVL); if (dbgCompParm) { printf("After compression: "); show_matrix(permutedIneqs); } /* 3- eliminate the first variables */ if (!mpolyhedron_eliminate_first_variables(permutedEqs, permutedIneqs)) { fprintf(stderr,"full-dimensionize > variable elimination failed. \n"); return NULL; } if (dbgCompParm) { printf("After elimination of the variables: "); show_matrix(permutedIneqs); } /* 4- get rid of the first (zero) columns, which are now useless, and put the parameters back at the end */ Full_Dim = Matrix_Alloc(permutedIneqs->NbRows, permutedIneqs->NbColumns-nbElimVars); for (i=0; i< permutedIneqs->NbRows; i++) { value_set_si(Full_Dim->p[i][0], 1); for (j=0; j< nbParms; j++) value_assign(Full_Dim->p[i][j+Full_Dim->NbColumns-nbParms-1], permutedIneqs->p[i][j+nbElimVars+1]); for (j=0; j< permutedIneqs->NbColumns-nbParms-2-nbElimVars; j++) value_assign(Full_Dim->p[i][j+1], permutedIneqs->p[i][nbElimVars+nbParms+j+1]); value_assign(Full_Dim->p[i][Full_Dim->NbColumns-1], permutedIneqs->p[i][permutedIneqs->NbColumns-1]); } Matrix_Free(permutedIneqs); /* 5- Keep only the the validity lattice restricted to the parameters */ *validityLattice = Matrix_Alloc(nbParms+1, nbParms+1); for (i=0; i< nbParms; i++) { for (j=0; j< nbParms; j++) value_assign((*validityLattice)->p[i][j], WVL->p[i][j]); value_assign((*validityLattice)->p[i][nbParms], WVL->p[i][WVL->NbColumns-1]); } for (j=0; j< nbParms; j++) value_set_si((*validityLattice)->p[nbParms][j], 0); value_assign((*validityLattice)->p[nbParms][nbParms], WVL->p[WVL->NbColumns-1][WVL->NbColumns-1]); /* 6- Clean up */ Matrix_Free(WVL); return Full_Dim; } /* full_dimensionize */
/** Removes the equalities that involve only parameters, by eliminating some * parameters in the polyhedron's constraints and in the context.<p> * <b>Updates M and Ctxt.</b> * @param M1 the polyhedron's constraints * @param Ctxt1 the constraints of the polyhedron's context * @param renderSpace tells if the returned equalities must be expressed in the * parameters space (renderSpace=0) or in the combined var/parms space * (renderSpace = 1) * @param elimParms the list of parameters that have been removed: an array * whose 1st element is the number of elements in the list. (returned) * @return the system of equalities that involve only parameters. */ Matrix * Constraints_Remove_parm_eqs(Matrix ** M1, Matrix ** Ctxt1, int renderSpace, unsigned int ** elimParms) { int i, j, k, nbEqsParms =0; int nbEqsM, nbEqsCtxt, allZeros, nbTautoM = 0, nbTautoCtxt = 0; Matrix * M = (*M1); Matrix * Ctxt = (*Ctxt1); int nbVars = M->NbColumns-Ctxt->NbColumns; Matrix * Eqs; Matrix * EqsMTmp; /* 1- build the equality matrix(ces) */ nbEqsM = 0; for (i=0; i< M->NbRows; i++) { k = First_Non_Zero(M->p[i], M->NbColumns); /* if it is a tautology, count it as such */ if (k==-1) { nbTautoM++; } else { /* if it only involves parameters, count it */ if (k>= nbVars+1) nbEqsM++; } } nbEqsCtxt = 0; for (i=0; i< Ctxt->NbRows; i++) { if (value_zero_p(Ctxt->p[i][0])) { if (First_Non_Zero(Ctxt->p[i], Ctxt->NbColumns)==-1) { nbTautoCtxt++; } else { nbEqsCtxt ++; } } } nbEqsParms = nbEqsM + nbEqsCtxt; /* nothing to do in this case */ if (nbEqsParms+nbTautoM+nbTautoCtxt==0) { (*elimParms) = (unsigned int*) malloc(sizeof(int)); (*elimParms)[0] = 0; if (renderSpace==0) { return Matrix_Alloc(0,Ctxt->NbColumns); } else { return Matrix_Alloc(0,M->NbColumns); } } Eqs= Matrix_Alloc(nbEqsParms, Ctxt->NbColumns); EqsMTmp= Matrix_Alloc(nbEqsParms, M->NbColumns); /* copy equalities from the context */ k = 0; for (i=0; i< Ctxt->NbRows; i++) { if (value_zero_p(Ctxt->p[i][0]) && First_Non_Zero(Ctxt->p[i], Ctxt->NbColumns)!=-1) { Vector_Copy(Ctxt->p[i], Eqs->p[k], Ctxt->NbColumns); Vector_Copy(Ctxt->p[i]+1, EqsMTmp->p[k]+nbVars+1, Ctxt->NbColumns-1); k++; } } for (i=0; i< M->NbRows; i++) { j=First_Non_Zero(M->p[i], M->NbColumns); /* copy equalities that involve only parameters from M */ if (j>=nbVars+1) { Vector_Copy(M->p[i]+nbVars+1, Eqs->p[k]+1, Ctxt->NbColumns-1); Vector_Copy(M->p[i]+nbVars+1, EqsMTmp->p[k]+nbVars+1, Ctxt->NbColumns-1); /* mark these equalities for removal */ value_set_si(M->p[i][0], 2); k++; } /* mark the all-zero equalities for removal */ if (j==-1) { value_set_si(M->p[i][0], 2); } } /* 2- eliminate parameters until all equalities are used or until we find a contradiction (overconstrained system) */ (*elimParms) = (unsigned int *) malloc((Eqs->NbRows+1) * sizeof(int)); (*elimParms)[0] = 0; allZeros = 0; for (i=0; i< Eqs->NbRows; i++) { /* find a variable that can be eliminated */ k = First_Non_Zero(Eqs->p[i], Eqs->NbColumns); if (k!=-1) { /* nothing special to do for tautologies */ /* if there is a contradiction, return empty matrices */ if (k==Eqs->NbColumns-1) { printf("Contradiction in %dth row of Eqs: ",k); show_matrix(Eqs); Matrix_Free(Eqs); Matrix_Free(EqsMTmp); (*M1) = Matrix_Alloc(0, M->NbColumns); Matrix_Free(M); (*Ctxt1) = Matrix_Alloc(0,Ctxt->NbColumns); Matrix_Free(Ctxt); free(*elimParms); (*elimParms) = (unsigned int *) malloc(sizeof(int)); (*elimParms)[0] = 0; if (renderSpace==1) { return Matrix_Alloc(0,(*M1)->NbColumns); } else { return Matrix_Alloc(0,(*Ctxt1)->NbColumns); } } /* if we have something we can eliminate, do it in 3 places: Eqs, Ctxt, and M */ else { k--; /* k is the rank of the variable, now */ (*elimParms)[0]++; (*elimParms)[(*elimParms[0])]=k; for (j=0; j< Eqs->NbRows; j++) { if (i!=j) { eliminate_var_with_constr(Eqs, i, Eqs, j, k); eliminate_var_with_constr(EqsMTmp, i, EqsMTmp, j, k+nbVars); } } for (j=0; j< Ctxt->NbRows; j++) { if (value_notzero_p(Ctxt->p[i][0])) { eliminate_var_with_constr(Eqs, i, Ctxt, j, k); } } for (j=0; j< M->NbRows; j++) { if (value_cmp_si(M->p[i][0], 2)) { eliminate_var_with_constr(EqsMTmp, i, M, j, k+nbVars); } } } } /* if (k==-1): count the tautologies in Eqs to remove them later */ else { allZeros++; } } /* elimParms may have been overallocated. Now we know how many parms have been eliminated so we can reallocate the right amount of memory. */ if (!realloc((*elimParms), ((*elimParms)[0]+1)*sizeof(int))) { fprintf(stderr, "Constraints_Remove_parm_eqs > cannot realloc()"); } Matrix_Free(EqsMTmp); /* 3- remove the "bad" equalities from the input matrices and copy the equalities involving only parameters */ EqsMTmp = Matrix_Alloc(M->NbRows-nbEqsM-nbTautoM, M->NbColumns); k=0; for (i=0; i< M->NbRows; i++) { if (value_cmp_si(M->p[i][0], 2)) { Vector_Copy(M->p[i], EqsMTmp->p[k], M->NbColumns); k++; } } Matrix_Free(M); (*M1) = EqsMTmp; EqsMTmp = Matrix_Alloc(Ctxt->NbRows-nbEqsCtxt-nbTautoCtxt, Ctxt->NbColumns); k=0; for (i=0; i< Ctxt->NbRows; i++) { if (value_notzero_p(Ctxt->p[i][0])) { Vector_Copy(Ctxt->p[i], EqsMTmp->p[k], Ctxt->NbColumns); k++; } } Matrix_Free(Ctxt); (*Ctxt1) = EqsMTmp; if (renderSpace==0) {/* renderSpace=0: equalities in the parameter space */ EqsMTmp = Matrix_Alloc(Eqs->NbRows-allZeros, Eqs->NbColumns); k=0; for (i=0; i<Eqs->NbRows; i++) { if (First_Non_Zero(Eqs->p[i], Eqs->NbColumns)!=-1) { Vector_Copy(Eqs->p[i], EqsMTmp->p[k], Eqs->NbColumns); k++; } } } else {/* renderSpace=1: equalities rendered in the combined space */ EqsMTmp = Matrix_Alloc(Eqs->NbRows-allZeros, (*M1)->NbColumns); k=0; for (i=0; i<Eqs->NbRows; i++) { if (First_Non_Zero(Eqs->p[i], Eqs->NbColumns)!=-1) { Vector_Copy(Eqs->p[i], &(EqsMTmp->p[k][nbVars]), Eqs->NbColumns); k++; } } } Matrix_Free(Eqs); Eqs = EqsMTmp; return Eqs; } /* Constraints_Remove_parm_eqs */
/** * Eliminates all the equalities in a set of constraints and returns the set of * constraints defining a full-dimensional polyhedron, such that there is a * bijection between integer points of the original polyhedron and these of the * resulting (projected) polyhedron). * If VL is set to NULL, this funciton allocates it. Else, it assumes that * (*VL) points to a matrix of the right size. * <p> The following things are done: * <ol> * <li> remove equalities involving only parameters, and remove as many * parameters as there are such equalities. From that, the list of * eliminated parameters <i>elimParms</i> is built. * <li> remove equalities that involve variables. This requires a compression * of the parameters and of the other variables that are not eliminated. * The affine compresson is represented by matrix VL (for <i>validity * lattice</i>) and is such that (N I 1)^T = VL.(N' I' 1), where N', I' * are integer (they are the parameters and variables after compression). *</ol> *</p> */ void Constraints_fullDimensionize(Matrix ** M, Matrix ** C, Matrix ** VL, Matrix ** Eqs, Matrix ** ParmEqs, unsigned int ** elimVars, unsigned int ** elimParms, int maxRays) { unsigned int i, j; Matrix * A=NULL, *B=NULL; Matrix * Ineqs=NULL; unsigned int nbVars = (*M)->NbColumns - (*C)->NbColumns; unsigned int nbParms; int nbElimVars; Matrix * fullDim = NULL; /* variables for permutations */ unsigned int * permutation; Matrix * permutedEqs=NULL, * permutedIneqs=NULL; /* 1- Eliminate the equalities involving only parameters. */ (*ParmEqs) = Constraints_removeParmEqs(M, C, 0, elimParms); /* if the polyehdron is empty, return now. */ if ((*M)->NbColumns==0) return; /* eliminate the columns corresponding to the eliminated parameters */ if (elimParms[0]!=0) { Constraints_removeElimCols(*M, nbVars, (*elimParms), &A); Matrix_Free(*M); (*M) = A; Constraints_removeElimCols(*C, 0, (*elimParms), &B); Matrix_Free(*C); (*C) = B; if (dbgCompParm) { printf("After false parameter elimination: \n"); show_matrix(*M); show_matrix(*C); } } nbParms = (*C)->NbColumns-2; /* 2- Eliminate the equalities involving variables */ /* a- extract the (remaining) equalities from the poyhedron */ split_constraints((*M), Eqs, &Ineqs); nbElimVars = (*Eqs)->NbRows; /* if the polyhedron is already full-dimensional, return */ if ((*Eqs)->NbRows==0) { Matrix_identity(nbParms+1, VL); return; } /* b- choose variables to be eliminated */ permutation = find_a_permutation((*Eqs), nbParms); if (dbgCompParm) { printf("Permuting the vars/parms this way: [ "); for (i=0; i< (*Eqs)->NbColumns-2; i++) { printf("%d ", permutation[i]); } printf("]\n"); } Constraints_permute((*Eqs), permutation, &permutedEqs); Equalities_validityLattice(permutedEqs, (*Eqs)->NbRows, VL); if (dbgCompParm) { printf("Validity lattice: "); show_matrix(*VL); } Constraints_compressLastVars(permutedEqs, (*VL)); Constraints_permute(Ineqs, permutation, &permutedIneqs); if (dbgCompParmMore) { show_matrix(permutedIneqs); show_matrix(permutedEqs); } Matrix_Free(*Eqs); Matrix_Free(Ineqs); Constraints_compressLastVars(permutedIneqs, (*VL)); if (dbgCompParm) { printf("After compression: "); show_matrix(permutedIneqs); } /* c- eliminate the first variables */ assert(Constraints_eliminateFirstVars(permutedEqs, permutedIneqs)); if (dbgCompParmMore) { printf("After elimination of the variables: "); show_matrix(permutedIneqs); } /* d- get rid of the first (zero) columns, which are now useless, and put the parameters back at the end */ fullDim = Matrix_Alloc(permutedIneqs->NbRows, permutedIneqs->NbColumns-nbElimVars); for (i=0; i< permutedIneqs->NbRows; i++) { value_set_si(fullDim->p[i][0], 1); for (j=0; j< nbParms; j++) { value_assign(fullDim->p[i][j+fullDim->NbColumns-nbParms-1], permutedIneqs->p[i][j+nbElimVars+1]); } for (j=0; j< permutedIneqs->NbColumns-nbParms-2-nbElimVars; j++) { value_assign(fullDim->p[i][j+1], permutedIneqs->p[i][nbElimVars+nbParms+j+1]); } value_assign(fullDim->p[i][fullDim->NbColumns-1], permutedIneqs->p[i][permutedIneqs->NbColumns-1]); } Matrix_Free(permutedIneqs); } /* Constraints_fullDimensionize */
/** * Given a system of equalities, looks if it has an integer solution in the * combined space, and if yes, returns one solution. * <p>pre-condition: the equalities are full-row rank (without the constant * part)</p> * @param Eqs the system of equations (as constraints) * @param I a feasible integer solution if it exists, else NULL. Allocated if * initially set to NULL, else reused. */ void Equalities_integerSolution(Matrix * Eqs, Matrix **I) { Matrix * Hm, *H=NULL, *U, *Q, *M=NULL, *C=NULL, *Hi; Matrix *Ip; int i; Value mod; unsigned int rk; if (Eqs==NULL){ if ((*I)!=NULL) Matrix_Free(*I); I = NULL; return; } /* we use: AI = C = (Ha 0).Q.I = (Ha 0)(I' 0)^T */ /* with I = Qinv.I' = U.I'*/ /* 1- compute I' = Hainv.(-C) */ /* HYP: the equalities are full-row rank */ rk = Eqs->NbRows; Matrix_subMatrix(Eqs, 0, 1, rk, Eqs->NbColumns-1, &M); left_hermite(M, &Hm, &Q, &U); Matrix_Free(M); Matrix_subMatrix(Hm, 0, 0, rk, rk, &H); if (dbgCompParmMore) { show_matrix(Hm); show_matrix(H); show_matrix(U); } Matrix_Free(Q); Matrix_Free(Hm); Matrix_subMatrix(Eqs, 0, Eqs->NbColumns-1, rk, Eqs->NbColumns, &C); Matrix_oppose(C); Hi = Matrix_Alloc(rk, rk+1); MatInverse(H, Hi); if (dbgCompParmMore) { show_matrix(C); show_matrix(Hi); } /* put the numerator of Hinv back into H */ Matrix_subMatrix(Hi, 0, 0, rk, rk, &H); Ip = Matrix_Alloc(Eqs->NbColumns-2, 1); /* fool Matrix_Product on the size of Ip */ Ip->NbRows = rk; Matrix_Product(H, C, Ip); Ip->NbRows = Eqs->NbColumns-2; Matrix_Free(H); Matrix_Free(C); value_init(mod); for (i=0; i< rk; i++) { /* if Hinv.C is not integer, return NULL (no solution) */ value_pmodulus(mod, Ip->p[i][0], Hi->p[i][rk]); if (value_notzero_p(mod)) { if ((*I)!=NULL) Matrix_Free(*I); value_clear(mod); Matrix_Free(U); Matrix_Free(Ip); Matrix_Free(Hi); I = NULL; return; } else { value_pdivision(Ip->p[i][0], Ip->p[i][0], Hi->p[i][rk]); } } /* fill the rest of I' with zeros */ for (i=rk; i< Eqs->NbColumns-2; i++) { value_set_si(Ip->p[i][0], 0); } value_clear(mod); Matrix_Free(Hi); /* 2 - Compute the particular solution I = U.(I' 0) */ ensureMatrix((*I), Eqs->NbColumns-2, 1); Matrix_Product(U, Ip, (*I)); Matrix_Free(U); Matrix_Free(Ip); if (dbgCompParm) { show_matrix(*I); } }
/* * Given a rational matrix 'Mat'(k x k), compute its inverse rational matrix * 'MatInv' k x k. * The output is 1, * if 'Mat' is non-singular (invertible), otherwise the output is 0. Note:: * (1) Matrix 'Mat' is modified during the inverse operation. * (2) Matrix 'MatInv' must be preallocated before passing into this function. */ int Matrix_Inverse(Matrix *Mat,Matrix *MatInv ) { int i, k, j, c; Value x, gcd, piv; Value m1,m2; Value *den; if(Mat->NbRows != Mat->NbColumns) { fprintf(stderr,"Trying to invert a non-square matrix !\n"); return 0; } /* Initialize all the 'Value' variables */ value_init(x); value_init(gcd); value_init(piv); value_init(m1); value_init(m2); k = Mat->NbRows; /* Initialise MatInv */ Vector_Set(MatInv->p[0],0,k*k); /* Initialize 'MatInv' to Identity matrix form. Each diagonal entry is set*/ /* to 1. Last column of each row (denominator of each entry in a row) is */ /* also set to 1. */ for(i=0;i<k;++i) { value_set_si(MatInv->p[i][i],1); /* value_set_si(MatInv->p[i][k],1); /* denum */ } /* Apply Gauss-Jordan elimination method on the two matrices 'Mat' and */ /* 'MatInv' in parallel. */ for(i=0;i<k;++i) { /* Check if the diagonal entry (new pivot) is non-zero or not */ if(value_zero_p(Mat->p[i][i])) { /* Search for a non-zero pivot down the column(i) */ for(j=i;j<k;++j) if(value_notzero_p(Mat->p[j][i])) break; /* If no non-zero pivot is found, the matrix 'Mat' is non-invertible */ /* Return 0. */ if(j==k) { /* Clear all the 'Value' variables */ value_clear(x); value_clear(gcd); value_clear(piv); value_clear(m1); value_clear(m2); return 0; } /* Exchange the rows, row(i) and row(j) so that the diagonal element */ /* Mat->p[i][i] (pivot) is non-zero. Repeat the same operations on */ /* matrix 'MatInv'. */ for(c=0;c<k;++c) { /* Interchange rows, row(i) and row(j) of matrix 'Mat' */ value_assign(x,Mat->p[j][c]); value_assign(Mat->p[j][c],Mat->p[i][c]); value_assign(Mat->p[i][c],x); /* Interchange rows, row(i) and row(j) of matrix 'MatInv' */ value_assign(x,MatInv->p[j][c]); value_assign(MatInv->p[j][c],MatInv->p[i][c]); value_assign(MatInv->p[i][c],x); } } /* Make all the entries in column(i) of matrix 'Mat' zero except the */ /* diagonal entry. Repeat the same sequence of operations on matrix */ /* 'MatInv'. */ for(j=0;j<k;++j) { if (j==i) continue; /* Skip the pivot */ value_assign(x,Mat->p[j][i]); if(value_notzero_p(x)) { value_assign(piv,Mat->p[i][i]); value_gcd(gcd, x, piv); if (value_notone_p(gcd) ) { value_divexact(x, x, gcd); value_divexact(piv, piv, gcd); } for(c=((j>i)?i:0);c<k;++c) { value_multiply(m1,piv,Mat->p[j][c]); value_multiply(m2,x,Mat->p[i][c]); value_subtract(Mat->p[j][c],m1,m2); } for(c=0;c<k;++c) { value_multiply(m1,piv,MatInv->p[j][c]); value_multiply(m2,x,MatInv->p[i][c]); value_subtract(MatInv->p[j][c],m1,m2); } /* Simplify row(j) of the two matrices 'Mat' and 'MatInv' by */ /* dividing the rows with the common GCD. */ Vector_Gcd(&MatInv->p[j][0],k,&m1); Vector_Gcd(&Mat->p[j][0],k,&m2); value_gcd(gcd, m1, m2); if(value_notone_p(gcd)) { for(c=0;c<k;++c) { value_divexact(Mat->p[j][c], Mat->p[j][c], gcd); value_divexact(MatInv->p[j][c], MatInv->p[j][c], gcd); } } } } } /* Find common denom for each row */ den = (Value *)malloc(k*sizeof(Value)); value_set_si(x,1); for(j=0 ; j<k ; ++j) { value_init(den[j]); value_assign(den[j],Mat->p[j][j]); /* gcd is always positive */ Vector_Gcd(&MatInv->p[j][0],k,&gcd); value_gcd(gcd, gcd, den[j]); if (value_neg_p(den[j])) value_oppose(gcd,gcd); /* make denominator positive */ if (value_notone_p(gcd)) { for (c=0; c<k; c++) value_divexact(MatInv->p[j][c], MatInv->p[j][c], gcd); /* normalize */ value_divexact(den[j], den[j], gcd); } value_gcd(gcd, x, den[j]); value_divexact(m1, den[j], gcd); value_multiply(x,x,m1); } if (value_notone_p(x)) for(j=0 ; j<k ; ++j) { for (c=0; c<k; c++) { value_division(m1,x,den[j]); value_multiply(MatInv->p[j][c],MatInv->p[j][c],m1); /* normalize */ } } /* Clear all the 'Value' variables */ for(j=0 ; j<k ; ++j) { value_clear(den[j]); } value_clear(x); value_clear(gcd); value_clear(piv); value_clear(m1); value_clear(m2); free(den); return 1; } /* Matrix_Inverse */
/* Compute integer hull of truncated linear cone C, i.e., of C with * the origin removed. * Here, we do this by first computing the Hilbert basis of C * and then discarding elements from this basis that are rational * overconvex combinations of other elements in the basis. */ Matrix *Cone_Hilbert_Integer_Hull(Polyhedron *C, struct barvinok_options *options) { int i, j, k; Matrix *hilbert = Cone_Hilbert_Basis(C, options->MaxRays); Matrix *rays, *hull; unsigned dim = C->Dimension; Value tmp; unsigned MaxRays = options->MaxRays; /* When checking for redundant points below, we want to * check if there are any _rational_ solutions. */ POL_UNSET(options->MaxRays, POL_INTEGER); POL_ENSURE_VERTICES(C); rays = Matrix_Alloc(C->NbRays-1, C->Dimension); for (i = 0, j = 0; i < C->NbRays; ++i) { if (value_notzero_p(C->Ray[i][1+C->Dimension])) continue; Vector_Copy(C->Ray[i]+1, rays->p[j++], C->Dimension); } /* We only sort the pointers into the big Value array */ qsort(rays->p, rays->NbRows, sizeof(Value *), lex_cmp); qsort(hilbert->p, hilbert->NbRows, sizeof(Value *), lex_cmp); /* Remove rays from Hilbert basis */ for (i = 0, j = 0, k = 0; i < hilbert->NbRows && j < rays->NbRows; ++i) { if (Vector_Equal(hilbert->p[i], rays->p[j], C->Dimension)) ++j; else hilbert->p[k++] = hilbert->p[i]; } hilbert->NbRows = k; /* Now remove points that are overconvex combinations of other points */ value_init(tmp); for (i = 0; hilbert->NbRows > 1 && i < hilbert->NbRows; ++i) { Matrix *LP; Vector *obj; int nray = rays->NbRows; int npoint = hilbert->NbRows; enum lp_result result; LP = Matrix_Alloc(dim + 1 + nray + (npoint-1), 2 + nray + (npoint-1)); for (j = 0; j < dim; ++j) { for (k = 0; k < nray; ++k) value_assign(LP->p[j][k+1], rays->p[k][j]); for (k = 0; k < npoint; ++k) { if (k == i) value_oppose(LP->p[j][1+nray+npoint-1], hilbert->p[k][j]); else value_assign(LP->p[j][1+nray+k-(k>i)], hilbert->p[k][j]); } } value_set_si(LP->p[dim][0], 1); for (k = 0; k < nray+npoint-1; ++k) value_set_si(LP->p[dim][1+k], 1); value_set_si(LP->p[dim][LP->NbColumns-1], -1); for (k = 0; k < LP->NbColumns-2; ++k) { value_set_si(LP->p[dim+1+k][0], 1); value_set_si(LP->p[dim+1+k][1+k], 1); } /* Somewhat arbitrary objective function. */ obj = Vector_Alloc(LP->NbColumns-1); value_set_si(obj->p[0], 1); value_set_si(obj->p[obj->Size-1], 1); result = constraints_opt(LP, obj->p, obj->p[0], lp_min, &tmp, options); /* If the LP is not empty, the point can be discarded */ if (result != lp_empty) { hilbert->NbRows--; if (i < hilbert->NbRows) hilbert->p[i] = hilbert->p[hilbert->NbRows]; --i; } Matrix_Free(LP); Vector_Free(obj); } value_clear(tmp); hull = Matrix_Alloc(rays->NbRows + hilbert->NbRows, dim+1); for (i = 0; i < rays->NbRows; ++i) { Vector_Copy(rays->p[i], hull->p[i], dim); value_set_si(hull->p[i][dim], 1); } for (i = 0; i < hilbert->NbRows; ++i) { Vector_Copy(hilbert->p[i], hull->p[rays->NbRows+i], dim); value_set_si(hull->p[rays->NbRows+i][dim], 1); } Matrix_Free(rays); Matrix_Free(hilbert); options->MaxRays = MaxRays; return hull; }