Main Page | Class List | File List | Class Members | File Members

matrix_addon.c

Go to the documentation of this file.
00001 /** 
00002  * $Id: matrix_addon.c,v 1.17 2007/03/18 18:49:08 skimo Exp $
00003  * 
00004  * Polylib matrix addons
00005  * Mainly, deals with polyhedra represented as a matrix (implicit form)
00006  * @author Benoit Meister <meister@icps.u-strasbg.fr>
00007  * 
00008  */
00009 
00010 #include <stdlib.h>
00011 #include<polylib/polylib.h>
00012 #include <polylib/matrix_addon.h>
00013 
00014 /** Creates a view of the constraints of a polyhedron as a Matrix * */
00015 Matrix * constraintsView(Polyhedron * P) {
00016   Matrix * view = (Matrix *)malloc(sizeof(Matrix));
00017   view->NbRows = P->NbConstraints;
00018   view->NbColumns = P->Dimension+2;
00019   view->p = P->Constraint;
00020   return view;
00021 }
00022 
00023 /** "Frees" a view of the constraints of a polyhedron */
00024 void constraintsView_Free(Matrix * M) {
00025   free(M);
00026 }
00027 
00028 /** 
00029  * splits a matrix of constraints M into a matrix of equalities Eqs and a
00030  *  matrix of inequalities Ineqs allocs the new matrices. 
00031  * Allocates Eqs and Ineqs.
00032 */
00033 void split_constraints(Matrix const * M, Matrix ** Eqs, Matrix **Ineqs) {
00034   unsigned int i, j, k_eq, k_ineq, nb_eqs=0;
00035 
00036   /* 1- count the number of equations */
00037   for (i=0; i< M->NbRows; i++)     
00038     if (value_zero_p(M->p[i][0])) nb_eqs++;
00039 
00040   /* 2- extract the two matrices of equations */
00041   (*Eqs) = Matrix_Alloc(nb_eqs, M->NbColumns);
00042   (*Ineqs) = Matrix_Alloc(M->NbRows-nb_eqs, M->NbColumns);
00043 
00044   k_eq = k_ineq = 0;
00045   for(i=0; i< M->NbRows; i++) {
00046     if (value_zero_p(M->p[i][0])) 
00047       {
00048         for(j=0; j< M->NbColumns; j++)
00049           value_assign((*Eqs)->p[k_eq][j], M->p[i][j]);
00050         k_eq++;
00051       }
00052     else
00053        {
00054         for(j=0; j< M->NbColumns; j++)
00055           value_assign((*Ineqs)->p[k_ineq][j], M->p[i][j]);
00056         k_ineq++;
00057       }
00058   }
00059 }
00060 
00061 
00062 /* returns the dim-dimensional identity matrix */
00063 Matrix * Identity_Matrix(unsigned int dim) {
00064   Matrix * ret = Matrix_Alloc(dim, dim);
00065   unsigned int i,j;
00066   for (i=0; i< dim; i++) {
00067     for (j=0; j< dim; j++) {
00068       if (i==j) 
00069         { value_set_si(ret->p[i][j], 1); } 
00070       else value_set_si(ret->p[i][j], 0);
00071     }
00072   }
00073   return ret;
00074 } /* Identity_Matrix */
00075 
00076 
00077 /** 
00078  * returns the dim-dimensional identity matrix. 
00079  * If I is set to NULL, allocates it first. 
00080  * Else, assumes an existing, allocated Matrix.
00081 */
00082 void Matrix_identity(unsigned int dim, Matrix ** I) {
00083   int i,j;
00084   if (*I==NULL) {
00085     (*I) = Identity_Matrix(dim);
00086   }
00087   else {
00088     assert((*I)->NbRows>=dim && (*I)->NbColumns>=dim);
00089     for (i=0; i< dim; i++) {
00090       for (j=0; j< dim; j++) {
00091         if (i==j) { 
00092             value_set_si((*I)->p[i][j], 1); 
00093           } 
00094         else {
00095           value_set_si((*I)->p[i][j], 0);
00096         }
00097       }
00098     }
00099   }
00100 } /* Matrix_identity */
00101 
00102 
00103 /** given a n x n integer transformation matrix transf, compute its inverse
00104     M/g, where M is a nxn integer matrix.  g is a common denominator for
00105     elements of (transf^{-1}) */
00106 void mtransformation_inverse(Matrix * transf, Matrix ** inverse, Value * g) {
00107   Value factor;
00108   unsigned int i,j;
00109   Matrix *tmp, *inv;
00110 
00111   value_init(*g);
00112   value_set_si(*g,1);
00113 
00114   /* a - compute the inverse as usual (n x (n+1) matrix) */
00115   tmp = Matrix_Copy(transf);
00116   inv = Matrix_Alloc(transf->NbRows, transf->NbColumns+1);
00117   MatInverse(tmp, inv);
00118   Matrix_Free(tmp);
00119 
00120   /* b - as it is rational, put it to the same denominator*/
00121   (*inverse) = Matrix_Alloc(transf->NbRows, transf->NbRows);
00122   for (i=0; i< inv->NbRows; i++) 
00123     value_lcm(*g, *g, inv->p[i][inv->NbColumns-1]);
00124   for (i=0; i< inv->NbRows; i++) {
00125     value_division(factor, *g, inv->p[i][inv->NbColumns-1]);
00126     for (j=0; j< (*inverse)->NbColumns; j++) 
00127       value_multiply((*inverse)->p[i][j], inv->p[i][j],  factor);
00128   }
00129 
00130   /* c- clean up */
00131   value_clear(factor);
00132   Matrix_Free(inv);
00133 } /* mtransformation_inverse */
00134 
00135 
00136 /** takes a transformation matrix, and expands it to a higher dimension with
00137     the identity matrix regardless of it homogeneousness */
00138 Matrix * mtransformation_expand_left_to_dim(Matrix * M, int new_dim) {
00139   Matrix * ret = Identity_Matrix(new_dim);
00140   int offset = new_dim-M->NbRows;
00141   unsigned int i,j;
00142 
00143   assert(new_dim>=M->NbColumns);
00144   assert(M->NbRows==M->NbColumns);
00145 
00146   for (i=0; i< M->NbRows; i++)
00147     for (j=0; j< M->NbRows; j++)
00148       value_assign(ret->p[offset+i][offset+j], M->p[i][j]);
00149   return ret;
00150 } /* mtransformation_expand_left_to_dim */
00151 
00152 
00153 /** simplify a matrix seen as a polyhedron, by dividing its rows by the gcd of
00154    their elements. */
00155 void mpolyhedron_simplify(Matrix * polyh) {
00156   int i, j;
00157   Value cur_gcd;
00158   value_init(cur_gcd);
00159   for (i=0; i< polyh->NbRows; i++) {
00160     Vector_Gcd(polyh->p[i]+1, polyh->NbColumns-1, &cur_gcd);
00161     printf(" gcd[%d] = ", i); 
00162     value_print(stdout, VALUE_FMT, cur_gcd);printf("\n");
00163     Vector_AntiScale(polyh->p[i]+1, polyh->p[i]+1, cur_gcd, polyh->NbColumns-1);
00164   }
00165   value_clear(cur_gcd);
00166 } /* mpolyhedron_simplify */
00167 
00168 
00169 /** inflates a polyhedron (represented as a matrix) P, so that the apx of its
00170     Ehrhart Polynomial is an upper bound of the Ehrhart polynomial of P
00171     WARNING: this inflation is supposed to be applied on full-dimensional
00172     polyhedra. */
00173 void mpolyhedron_inflate(Matrix * polyh, unsigned int nb_parms) {
00174   unsigned int i,j;
00175   unsigned nb_vars = polyh->NbColumns-nb_parms-2;
00176   Value infl;
00177   value_init(infl);
00178   /* subtract the sum of the negative coefficients of each inequality */
00179   for (i=0; i< polyh->NbRows; i++) {
00180     value_set_si(infl, 0);
00181     for (j=0; j< nb_vars; j++) {
00182       if (value_neg_p(polyh->p[i][1+j]))
00183         value_addto(infl, infl, polyh->p[i][1+j]);
00184     }
00185     /* here, we subtract a negative value */
00186     value_subtract(polyh->p[i][polyh->NbColumns-1], 
00187                    polyh->p[i][polyh->NbColumns-1], infl);
00188   }
00189   value_clear(infl);
00190 } /* mpolyhedron_inflate */
00191 
00192 
00193 /** deflates a polyhedron (represented as a matrix) P, so that the apx of its
00194     Ehrhart Polynomial is a lower bound of the Ehrhart polynomial of P WARNING:
00195     this deflation is supposed to be applied on full-dimensional polyhedra. */
00196 void mpolyhedron_deflate(Matrix * polyh, unsigned int nb_parms) {
00197   unsigned int i,j;
00198   unsigned nb_vars = polyh->NbColumns-nb_parms-2;
00199   Value defl;
00200   value_init(defl);
00201   /* substract the sum of the negative coefficients of each inequality */
00202   for (i=0; i< polyh->NbRows; i++) {
00203     value_set_si(defl, 0);
00204     for (j=0; j< nb_vars; j++) {
00205       if (value_pos_p(polyh->p[i][1+j]))
00206         value_addto(defl, defl, polyh->p[i][1+j]);
00207     }
00208     /* here, we substract a negative value */
00209     value_subtract(polyh->p[i][polyh->NbColumns-1], 
00210                    polyh->p[i][polyh->NbColumns-1], defl);
00211   }
00212   value_clear(defl);
00213 } /* mpolyhedron_deflate */
00214 
00215 
00216 /** use an eliminator row to eliminate a variable in a victim row (without
00217  * changing the sign of the victim row -> important if it is an inequality).
00218  * @param Eliminator the matrix containing the eliminator row
00219  * @param eliminator_row the index of the eliminator row in <tt>Eliminator</tt>
00220  * @param Victim the matrix containing the row to be eliminated
00221  * @param victim_row the row to be eliminated in <tt>Victim</tt>
00222  * @param var_to_elim the variable to be eliminated.
00223  */
00224 void eliminate_var_with_constr(Matrix * Eliminator, 
00225                                unsigned int eliminator_row, Matrix * Victim, 
00226                                unsigned int victim_row, 
00227                                unsigned int var_to_elim) {
00228   Value cur_lcm, mul_a, mul_b;
00229   Value tmp, tmp2;
00230   int k; 
00231 
00232   value_init(cur_lcm); 
00233   value_init(mul_a); 
00234   value_init(mul_b); 
00235   value_init(tmp); 
00236   value_init(tmp2);
00237   /* if the victim coefficient is not zero */
00238   if (value_notzero_p(Victim->p[victim_row][var_to_elim+1])) {
00239     value_lcm(cur_lcm, Eliminator->p[eliminator_row][var_to_elim+1], 
00240               Victim->p[victim_row][var_to_elim+1]);
00241     /* multiplication factors */
00242     value_division(mul_a, cur_lcm, 
00243                    Eliminator->p[eliminator_row][var_to_elim+1]);
00244     value_division(mul_b, cur_lcm, 
00245                    Victim->p[victim_row][var_to_elim+1]);
00246     /* the multiplicator for the vitim row has to be positive */
00247     if (value_pos_p(mul_b)) {
00248       value_oppose(mul_a, mul_a);
00249     }
00250     else {
00251       value_oppose(mul_b, mul_b);
00252     }
00253     value_clear(cur_lcm); 
00254     /* now we have a.mul_a = -(b.mul_b) and mul_a > 0 */
00255     for (k=1; k<Victim->NbColumns; k++) {
00256       value_multiply(tmp, Eliminator->p[eliminator_row][k], mul_a);
00257       value_multiply(tmp2, Victim->p[victim_row][k], mul_b);
00258       value_addto(Victim->p[victim_row][k], tmp, tmp2);
00259     }
00260   }
00261   value_clear(mul_a); 
00262   value_clear(mul_b); 
00263   value_clear(tmp); 
00264   value_clear(tmp2);
00265 }
00266 /* eliminate_var_with_constr */
00267 
00268 
00269 /* STUFF WITH PARTIAL MAPPINGS (Mappings to a subset of the
00270    variables/parameters) : on the first or last variables/parameters */
00271 
00272 /** compress the last vars/pars of the polyhedron M expressed as a polylib
00273     matrix
00274  - adresses the full-rank compressions only
00275  - modfies M */
00276 void mpolyhedron_compress_last_vars(Matrix * M, Matrix * compression) {
00277   unsigned int i, j, k;
00278   unsigned int offset = M->NbColumns - compression->NbRows; 
00279   /* the computations on M will begin on column "offset" */
00280 
00281   Matrix * M_tmp = Matrix_Alloc(1, M->NbColumns);
00282   assert(compression->NbRows==compression->NbColumns);
00283   /* basic matrix multiplication (using a temporary row instead of a whole
00284      temporary matrix), but with a column offset */
00285   for(i=0; i< M->NbRows; i++) {
00286     for (j=0; j< compression->NbRows; j++) {
00287       value_set_si(M_tmp->p[0][j], 0);
00288       for (k=0; k< compression->NbRows; k++) {
00289         value_addmul(M_tmp->p[0][j], M->p[i][k+offset],compression->p[k][j]);
00290       }
00291     }
00292     for (j=0; j< compression->NbRows; j++) 
00293       value_assign(M->p[i][j+offset], M_tmp->p[0][j]);
00294   }
00295   Matrix_Free(M_tmp);
00296 } /* mpolyhedron_compress_last_vars */
00297 
00298 
00299 /** use a set of m equalities Eqs to eliminate m variables in the polyhedron
00300     Ineqs represented as a matrix
00301  eliminates the m first variables
00302  - assumes that Eqs allow to eliminate the m equalities
00303  - modifies Eqs and Ineqs */
00304 unsigned int mpolyhedron_eliminate_first_variables(Matrix * Eqs, 
00305                                                    Matrix *Ineqs) {
00306   unsigned int i, j, k;
00307   /* eliminate one variable (index i) after each other */
00308   for (i=0; i< Eqs->NbRows; i++) {
00309     /* find j, the first (non-marked) row of Eqs with a non-zero coefficient */
00310     for (j=0; j<Eqs->NbRows && (Eqs->p[j][i+1]==0 || 
00311                                 ( !value_cmp_si(Eqs->p[j][0],2) )); 
00312          j++);
00313     /* if no row is found in Eqs that allows to eliminate variable i, return an
00314        error code (0) */
00315     if (j==Eqs->NbRows) return 0;
00316     /* else, eliminate variable i in Eqs and Ineqs with the j^th row of Eqs
00317        (and mark this row so we don't use it again for an elimination) */
00318     for (k=j+1; k<Eqs->NbRows; k++)
00319       eliminate_var_with_constr(Eqs, j, Eqs, k, i);
00320     for (k=0; k< Ineqs->NbRows; k++)
00321       eliminate_var_with_constr(Eqs, j, Ineqs, k, i);
00322     /* mark the row */
00323     value_set_si(Eqs->p[j][0],2);
00324   }
00325   /* un-mark all the rows */
00326   for (i=0; i< Eqs->NbRows; i++) value_set_si(Eqs->p[i][0],0);
00327   return 1;
00328 } /* mpolyhedron_eliminate_first_variables */
00329 
00330 
00331 /** returns a contiguous submatrix of a matrix.
00332  * @param M the input matrix
00333  * @param sr the index of the starting row
00334  * @param sc the index of the starting column
00335  * @param er the index ofthe ending row (excluded)
00336  * @param ec the ined of the ending colummn (excluded)
00337  * @param sub (returned), the submatrix. Allocated if set to NULL, assumed to
00338  * be already allocated else.
00339  */
00340 void Matrix_subMatrix(Matrix * M, unsigned int sr, unsigned int sc, 
00341                           unsigned int er, unsigned int ec, Matrix ** sub) {
00342   int i;
00343   int nbR = er-sr;
00344   int nbC = ec-sc;
00345   assert (er<=M->NbRows && ec<=M->NbColumns);
00346   if ((*sub)==NULL) {
00347     (*sub) = Matrix_Alloc(nbR, nbC);
00348   }
00349   if (nbR==0 || nbC==0) return;
00350   for (i=0; i< nbR; i++) {
00351     Vector_Copy(&(M->p[i+sr][sc]), (*sub)->p[i], nbC);
00352   }
00353 }/* Matrix_subMatrix */
00354 
00355 
00356 /**
00357  * Cloning function. Similar to Matrix_Copy() but allocates the target matrix
00358  * if it is set to NULL.
00359  */
00360 void Matrix_clone(Matrix * M, Matrix ** Cl) {
00361   Matrix_subMatrix(M, 0,0, M->NbRows, M->NbColumns, Cl);
00362 } 
00363 
00364 
00365 /**
00366  * Copies a contiguous submatrix of M1 into M2, at the indicated position.
00367  * M1 and M2 are assumed t be allocated already.
00368  * @param M1 the source matrix
00369  * @param sr1 the starting source row
00370  * @param sc1 the starting source column
00371  * @param nbR the number of rows
00372  * @param nbC the number of columns
00373  * @param M2 the target matrix
00374  * @param sr2 the starting target row
00375  * @param sc2 the starting target column
00376 */
00377 void Matrix_copySubMatrix(Matrix *M1,
00378                           unsigned int sr1, unsigned int sc1,
00379                           unsigned int nbR, unsigned int nbC,
00380                           Matrix * M2,
00381                           unsigned int sr2, unsigned int sc2) {
00382   int i;
00383   for (i=0; i< nbR; i++) {
00384     Vector_Copy(&(M1->p[i+sr1][sc1]), &(M2->p[i+sr2][sc2]), nbC);
00385   }
00386 } /* Matrix_copySubMatrix */
00387 
00388 
00389 /** 
00390  * transforms a matrix M into -M
00391  */
00392 void Matrix_oppose(Matrix * M) {
00393   int i,j;
00394   for (i=0; i<M->NbRows; i++) {
00395     for (j=0; j< M->NbColumns; j++) {
00396       value_oppose(M->p[i][j], M->p[i][j]);
00397     }
00398   }
00399 }

Generated on Thu Sep 4 15:28:58 2008 for polylib by doxygen 1.3.5