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

compress_parms.c

Go to the documentation of this file.
00001 /** 
00002  * $Id: compress_parms.c,v 1.32 2006/11/03 17:34:26 skimo Exp $
00003  *
00004  * The integer points in a parametric linear subspace of Q^n are generally
00005  * lying on a sub-lattice of Z^n.  
00006  * Functions here use and compute validity lattices, i.e. lattices induced on a
00007  * set of variables by such equalities involving another set of integer
00008  * variables.
00009  * @author B. Meister 12/2003-2006 meister@icps.u-strasbg.fr
00010  * LSIIT -ICPS 
00011  * UMR 7005 CNRS
00012  * Louis Pasteur University (ULP), Strasbourg, France 
00013 */
00014 
00015 #include <stdlib.h>
00016 #include <polylib/polylib.h>
00017 
00018 /** 
00019  * debug flags (2 levels)
00020  */
00021 #define dbgCompParm 0
00022 #define dbgCompParmMore 0
00023 
00024 #define dbgStart(a) if (dbgCompParmMore) { printf(" -- begin "); \
00025                                            printf(#a);        \
00026                                            printf(" --\n"); }   \
00027                                            while(0)
00028 #define dbgEnd(a) if (dbgCompParmMore) { printf(" -- end "); \
00029                                          printf(#a);      \
00030                                          printf(" --\n"); } \
00031                                          while(0)
00032 
00033 
00034 /** 
00035  * Given a full-row-rank nxm matrix M made of m row-vectors), computes the
00036  * basis K (made of n-m column-vectors) of the integer kernel of the rows of M
00037  * so we have: M.K = 0
00038 */
00039 Matrix * int_ker(Matrix * M) {
00040   Matrix *U, *Q, *H, *H2, *K=NULL;
00041   int i, j, rk;
00042 
00043   if (dbgCompParm)
00044     show_matrix(M);
00045   /* eliminate redundant rows : UM = H*/
00046   right_hermite(M, &H, &Q, &U);
00047   for (rk=H->NbRows-1; (rk>=0) && Vector_IsZero(H->p[rk], H->NbColumns); rk--);
00048   rk++;
00049   if (dbgCompParmMore) {
00050     printf("rank = %d\n", rk);
00051   }
00052     
00053   /* there is a non-null kernel if and only if the dimension m of 
00054      the space spanned by the rows 
00055      is inferior to the number n of variables */
00056   if (M->NbColumns <= rk) {
00057     Matrix_Free(H);
00058     Matrix_Free(Q);
00059     Matrix_Free(U);
00060     K = Matrix_Alloc(M->NbColumns, 0);
00061     return K;
00062   }
00063   Matrix_Free(U); 
00064   Matrix_Free(Q);
00065   /* fool left_hermite  by giving NbRows =rank of M*/
00066   H->NbRows=rk;
00067   /* computes MU = [H 0] */
00068   left_hermite(H, &H2, &Q, &U); 
00069    if (dbgCompParmMore) {
00070     printf("-- Int. Kernel -- \n");
00071     show_matrix(M);
00072     printf(" = \n");
00073     show_matrix(H2);
00074     show_matrix(U); 
00075   }
00076   H->NbRows==M->NbRows;
00077   Matrix_Free(H);
00078   /* the Integer Kernel is made of the last n-rk columns of U */
00079   Matrix_subMatrix(U, 0, rk, U->NbRows, U->NbColumns, &K);
00080 
00081   /* clean up */
00082   Matrix_Free(H2);
00083   Matrix_Free(U);
00084   Matrix_Free(Q);
00085   return K;
00086 } /* int_ker */
00087 
00088 
00089 /** 
00090  * Computes the intersection of two linear lattices, whose base vectors are
00091  * respectively represented in A and B.
00092  * If I and/or Lb is set to NULL, then the matrix is allocated. 
00093  * Else, the matrix is assumed to be allocated already. 
00094  * I and Lb are rk x rk, where rk is the rank of A (or B).
00095  * @param A the full-row rank matrix whose column-vectors are the basis for the
00096  * first linear lattice.
00097  * @param B the matrix whose column-vectors are the basis for the second linear
00098  * lattice.
00099  * @param Lb the matrix such that B.Lb = I, where I is the intersection.
00100  * @return their intersection.
00101  */
00102 static void linearInter(Matrix * A, Matrix * B, Matrix ** I, Matrix **Lb) {
00103   Matrix * AB=NULL;
00104   int rk = A->NbRows;
00105   int a = A->NbColumns;
00106   int b = B->NbColumns;
00107   int i,j, z=0;
00108 
00109   Matrix * H, *U, *Q;
00110   /* ensure that the spanning vectors are in the same space */
00111   assert(B->NbRows==rk);
00112   /* 1- build the matrix 
00113    * (A 0 1)
00114    * (0 B 1)
00115    */
00116   AB = Matrix_Alloc(2*rk, a+b+rk);
00117   Matrix_copySubMatrix(A, 0, 0, rk, a, AB, 0, 0);
00118   Matrix_copySubMatrix(B, 0, 0, rk, b, AB, rk, a);
00119   for (i=0; i< rk; i++) {
00120       value_set_si(AB->p[i][a+b+i], 1);
00121       value_set_si(AB->p[i+rk][a+b+i], 1);
00122   }
00123   if (dbgCompParm) {
00124     show_matrix(AB);
00125   }
00126 
00127   /* 2- Compute its left Hermite normal form. AB.U = [H 0] */
00128   left_hermite(AB, &H, &Q, &U);
00129   Matrix_Free(AB);
00130   Matrix_Free(Q);
00131   /* count the number of non-zero colums in H */ 
00132   for (z=H->NbColumns-1; value_zero_p(H->p[H->NbRows-1][z]); z--);
00133   z++;
00134   if (dbgCompParm) {
00135     show_matrix(H);
00136     printf("z=%d\n", z);
00137   }
00138   Matrix_Free(H);
00139   /* if you split U in 9 submatrices, you have: 
00140    * A.U_13 = -U_33
00141    * B.U_23 = -U_33,
00142    * where the nb of cols of U_{*3} equals the nb of zero-cols of H
00143    * U_33 is a (the smallest) combination of col-vectors of A and B at the same
00144    * time: their intersection.
00145   */
00146   Matrix_subMatrix(U, a+b, z, U->NbColumns, U->NbColumns, I);
00147   Matrix_subMatrix(U, a, z, a+b, U->NbColumns, Lb);
00148   if (dbgCompParm) {
00149     show_matrix(U);
00150   }
00151   Matrix_Free(U);
00152 } /* linearInter */
00153 
00154 
00155 /** 
00156  * Given a system of equalities, looks if it has an integer solution in the
00157  * combined space, and if yes, returns one solution.
00158  * <p>pre-condition: the equalities are full-row rank (without the constant
00159  * part)</p>
00160  * @param Eqs the system of equations (as constraints)
00161  * @param I a feasible integer solution if it exists, else NULL. Allocated if
00162  * initially set to NULL, else reused.
00163  */
00164 void Equalities_integerSolution(Matrix * Eqs, Matrix **I) {
00165   Matrix * Hm, *H=NULL, *U, *Q, *M=NULL, *C=NULL, *Hi;
00166   Matrix *Ip;
00167   int i;
00168   Value mod;
00169   unsigned int rk;
00170   if (Eqs==NULL){
00171     if ((*I)!=NULL) Matrix_Free(*I);
00172     I = NULL;
00173     return;
00174   }
00175   /* we use: AI = C = (Ha 0).Q.I = (Ha 0)(I' 0)^T */
00176   /* with I = Qinv.I' = U.I'*/
00177   /* 1- compute I' = Hainv.(-C) */
00178   /* HYP: the equalities are full-row rank */
00179   rk = Eqs->NbRows;
00180   Matrix_subMatrix(Eqs, 0, 1, rk, Eqs->NbColumns-1, &M);
00181   left_hermite(M, &Hm, &Q, &U);
00182   Matrix_Free(M);
00183   Matrix_subMatrix(Hm, 0, 0, rk, rk, &H);
00184   if (dbgCompParmMore) {
00185     show_matrix(Hm);
00186     show_matrix(H);
00187     show_matrix(U);
00188   }
00189   Matrix_Free(Q);
00190   Matrix_Free(Hm);
00191   Matrix_subMatrix(Eqs, 0, Eqs->NbColumns-1, rk, Eqs->NbColumns, &C);
00192   Matrix_oppose(C);
00193   Hi = Matrix_Alloc(rk, rk+1);
00194   MatInverse(H, Hi);
00195   if (dbgCompParmMore) {
00196     show_matrix(C);
00197     show_matrix(Hi);
00198   }
00199   /* put the numerator of Hinv back into H */
00200   Matrix_subMatrix(Hi, 0, 0, rk, rk, &H);
00201   Ip = Matrix_Alloc(Eqs->NbColumns-2, 1);
00202   /* fool Matrix_Product on the size of Ip */
00203   Ip->NbRows = rk;
00204   Matrix_Product(H, C, Ip);
00205   Ip->NbRows = Eqs->NbColumns-2;
00206   Matrix_Free(H);
00207   Matrix_Free(C);
00208   value_init(mod);
00209   for (i=0; i< rk; i++) {
00210     /* if Hinv.C is not integer, return NULL (no solution) */
00211     value_pmodulus(mod, Ip->p[i][0], Hi->p[i][rk]);
00212     if (value_notzero_p(mod)) { 
00213       if ((*I)!=NULL) Matrix_Free(*I);
00214       value_clear(mod);
00215       Matrix_Free(U);
00216       Matrix_Free(Ip);
00217       Matrix_Free(Hi);
00218       I = NULL;
00219       return;
00220     }
00221     else {
00222       value_pdivision(Ip->p[i][0], Ip->p[i][0], Hi->p[i][rk]);
00223     }
00224   }
00225   /* fill the rest of I' with zeros */
00226   for (i=rk; i< Eqs->NbColumns-2; i++) {
00227     value_set_si(Ip->p[i][0], 0);
00228   }
00229   value_clear(mod);
00230   Matrix_Free(Hi);
00231   /* 2 - Compute the particular solution I = U.(I' 0) */
00232   ensureMatrix((*I), Eqs->NbColumns-2, 1);
00233   Matrix_Product(U, Ip, (*I));
00234   Matrix_Free(U);
00235   Matrix_Free(Ip);
00236   if (dbgCompParm) {
00237     show_matrix(*I);
00238   }
00239 }
00240 
00241 
00242 /** 
00243  * Computes the validity lattice of a set of equalities. I.e., the lattice
00244  * induced on the last <tt>b</tt> variables by the equalities involving the
00245  * first <tt>a</tt> integer existential variables.  The submatrix of Eqs that
00246  * concerns only the existential variables (so the first a columns) is assumed
00247  * to be full-row rank.
00248  * @param Eqs the equalities
00249  * @param a the number of existential integer variables, placed as first
00250  * variables
00251  * @param vl the (returned) validity lattice, in homogeneous form. It is
00252  * allocated if initially set to null, or reused if already allocated.
00253  */
00254 void Equalities_validityLattice(Matrix * Eqs, int a, Matrix** vl) {
00255   unsigned int b = Eqs->NbColumns-2-a;
00256   unsigned int r = Eqs->NbRows;
00257   Matrix * A=NULL, * B=NULL, *I = NULL, *Lb=NULL, *sol=NULL;
00258   Matrix *H, *U, *Q;
00259   unsigned int i;
00260 
00261   if (dbgCompParm) {
00262     printf("Computing validity lattice induced by the %d first variables of:"
00263            ,a);
00264     show_matrix(Eqs);
00265   }
00266   if (b==0) {
00267     ensureMatrix((*vl), 1, 1);
00268     value_set_si((*vl)->p[0][0], 1);
00269     return;
00270   }
00271 
00272   /* 1- check that there is an integer solution to the equalities */
00273   /* OPT: could change integerSolution's profile to allocate or not*/
00274   Equalities_integerSolution(Eqs, &sol);
00275   /* if there is no integer solution, there is no validity lattice */
00276   if (sol==NULL) {
00277     if ((*vl)!=NULL) Matrix_Free(*vl);
00278     return;
00279   }
00280   Matrix_subMatrix(Eqs, 0, 1, r, 1+a, &A);
00281   Matrix_subMatrix(Eqs, 0, 1+a, r, 1+a+b, &B);
00282   linearInter(A, B, &I, &Lb);
00283   Matrix_Free(A);
00284   Matrix_Free(B);
00285   Matrix_Free(I);
00286   if (dbgCompParm) {
00287     show_matrix(Lb);
00288   }
00289   
00290   /* 2- The linear part of the validity lattice is the left HNF of Lb */
00291   left_hermite(Lb, &H, &Q, &U);
00292   Matrix_Free(Lb);
00293   Matrix_Free(Q);
00294   Matrix_Free(U);
00295 
00296   /* 3- build the validity lattice */
00297   ensureMatrix((*vl), b+1, b+1);
00298   Matrix_copySubMatrix(H, 0, 0, b, b, (*vl), 0,0);
00299   Matrix_Free(H);
00300   for (i=0; i< b; i++) {
00301     value_assign((*vl)->p[i][b], sol->p[0][a+i]);
00302   }
00303   Matrix_Free(sol);
00304   Vector_Set((*vl)->p[b],0, b);
00305   value_set_si((*vl)->p[b][b], 1);
00306   
00307 } /* validityLattice */
00308 
00309 
00310 /**
00311  * Eliminate the columns corresponding to a list of eliminated parameters.
00312  * @param M the constraints matrix whose columns are to be removed
00313  * @param nbVars an offset to be added to the ranks of the variables to be
00314  * removed
00315  * @param elimParms the list of ranks of the variables to be removed
00316  * @param newM (output) the matrix without the removed columns
00317  */
00318 void Constraints_removeElimCols(Matrix * M, unsigned int nbVars, 
00319                            unsigned int *elimParms, Matrix ** newM) {
00320   unsigned int i, j, k;
00321   if (elimParms[0]==0) {
00322     Matrix_clone(M, newM);
00323     return;
00324   }
00325   if ((*newM)==NULL) {
00326     (*newM) = Matrix_Alloc(M->NbRows, M->NbColumns - elimParms[0]);
00327   }
00328   else {
00329     assert ((*newM)->NbColumns==M->NbColumns - elimParms[0]);
00330   }
00331   for (i=0; i< M->NbRows; i++) {
00332     value_assign((*newM)->p[i][0], M->p[i][0]); /* kind of cstr */
00333     k=0;
00334     Vector_Copy(&(M->p[i][1]), &((*newM)->p[i][1]), nbVars);
00335     for (j=0; j< M->NbColumns-2-nbVars; j++) {
00336       if (j!=elimParms[k+1]) {
00337         value_assign((*newM)->p[i][j-k+nbVars+1], M->p[i][j+nbVars+1]);
00338       }
00339       else {
00340         k++;
00341       }
00342     }
00343     value_assign((*newM)->p[i][(*newM)->NbColumns-1], 
00344                  M->p[i][M->NbColumns-1]); /* cst part */
00345   }
00346 } /* Constraints_removeElimCols */
00347 
00348 
00349 /**
00350  * Eliminates all the equalities in a set of constraints and returns the set of
00351  * constraints defining a full-dimensional polyhedron, such that there is a
00352  * bijection between integer points of the original polyhedron and these of the
00353  * resulting (projected) polyhedron).
00354  * If VL is set to NULL, this funciton allocates it. Else, it assumes that
00355  * (*VL) points to a matrix of the right size.
00356  * <p> The following things are done: 
00357  * <ol>
00358  * <li> remove equalities involving only parameters, and remove as many
00359  *      parameters as there are such equalities. From that, the list of
00360  *      eliminated parameters <i>elimParms</i> is built.
00361  * <li> remove equalities that involve variables. This requires a compression
00362  *      of the parameters and of the other variables that are not eliminated.
00363  *      The affine compresson is represented by matrix VL (for <i>validity
00364  *      lattice</i>) and is such that (N I 1)^T = VL.(N' I' 1), where N', I'
00365  *      are integer (they are the parameters and variables after compression).
00366  *</ol>
00367  *</p>
00368  */
00369 void Constraints_fullDimensionize(Matrix ** M, Matrix ** C, Matrix ** VL, 
00370                                   Matrix ** Eqs, Matrix ** ParmEqs, 
00371                                   unsigned int ** elimVars, 
00372                                   unsigned int ** elimParms,
00373                                   int maxRays) {
00374   unsigned int i, j;
00375   Matrix * A=NULL, *B=NULL;
00376   Matrix * Ineqs=NULL;
00377   unsigned int nbVars = (*M)->NbColumns - (*C)->NbColumns;
00378   unsigned int nbParms;
00379   int nbElimVars;
00380   Matrix * fullDim = NULL;
00381 
00382   /* variables for permutations */
00383   unsigned int * permutation;
00384   Matrix * permutedEqs=NULL, * permutedIneqs=NULL;
00385   
00386   /* 1- Eliminate the equalities involving only parameters. */
00387   (*ParmEqs) = Constraints_removeParmEqs(M, C, 0, elimParms);
00388   /* if the polyehdron is empty, return now. */
00389   if ((*M)->NbColumns==0) return;
00390   /* eliminate the columns corresponding to the eliminated parameters */
00391   if (elimParms[0]!=0) {
00392     Constraints_removeElimCols(*M, nbVars, (*elimParms), &A);
00393     Matrix_Free(*M);
00394     (*M) = A;
00395     Constraints_removeElimCols(*C, 0, (*elimParms), &B);
00396     Matrix_Free(*C);
00397     (*C) = B;
00398     if (dbgCompParm) {
00399       printf("After false parameter elimination: \n");
00400       show_matrix(*M);
00401       show_matrix(*C);
00402     }
00403   }
00404   nbParms = (*C)->NbColumns-2;
00405 
00406   /* 2- Eliminate the equalities involving variables */
00407   /*   a- extract the (remaining) equalities from the poyhedron */
00408   split_constraints((*M), Eqs, &Ineqs);
00409   nbElimVars = (*Eqs)->NbRows;
00410   /*    if the polyhedron is already full-dimensional, return */
00411   if ((*Eqs)->NbRows==0) {
00412     Matrix_identity(nbParms+1, VL);
00413     return;
00414   }
00415   /*   b- choose variables to be eliminated */
00416   permutation = find_a_permutation((*Eqs), nbParms);
00417 
00418   if (dbgCompParm) {
00419     printf("Permuting the vars/parms this way: [ ");
00420     for (i=0; i< (*Eqs)->NbColumns-2; i++) {
00421       printf("%d ", permutation[i]);
00422     }
00423     printf("]\n");
00424   }
00425 
00426   Constraints_permute((*Eqs), permutation, &permutedEqs);
00427   Equalities_validityLattice(permutedEqs, (*Eqs)->NbRows, VL);
00428 
00429   if (dbgCompParm) {
00430     printf("Validity lattice: ");
00431     show_matrix(*VL);
00432   }
00433   Constraints_compressLastVars(permutedEqs, (*VL));
00434   Constraints_permute(Ineqs, permutation, &permutedIneqs);
00435   if (dbgCompParmMore) {
00436     show_matrix(permutedIneqs);
00437     show_matrix(permutedEqs);
00438   }
00439   Matrix_Free(*Eqs);
00440   Matrix_Free(Ineqs);
00441   Constraints_compressLastVars(permutedIneqs, (*VL));
00442   if (dbgCompParm) {
00443     printf("After compression: ");
00444     show_matrix(permutedIneqs);
00445   }
00446   /*   c- eliminate the first variables */
00447   assert(Constraints_eliminateFirstVars(permutedEqs, permutedIneqs));
00448   if (dbgCompParmMore) {
00449     printf("After elimination of the variables: ");
00450     show_matrix(permutedIneqs);
00451   }
00452 
00453   /*   d- get rid of the first (zero) columns, 
00454        which are now useless, and put the parameters back at the end */
00455   fullDim = Matrix_Alloc(permutedIneqs->NbRows,
00456                          permutedIneqs->NbColumns-nbElimVars);
00457   for (i=0; i< permutedIneqs->NbRows; i++) {
00458     value_set_si(fullDim->p[i][0], 1);
00459     for (j=0; j< nbParms; j++) {
00460       value_assign(fullDim->p[i][j+fullDim->NbColumns-nbParms-1], 
00461                    permutedIneqs->p[i][j+nbElimVars+1]);
00462     }
00463     for (j=0; j< permutedIneqs->NbColumns-nbParms-2-nbElimVars; j++) {
00464       value_assign(fullDim->p[i][j+1], 
00465                    permutedIneqs->p[i][nbElimVars+nbParms+j+1]);
00466     }
00467     value_assign(fullDim->p[i][fullDim->NbColumns-1], 
00468                  permutedIneqs->p[i][permutedIneqs->NbColumns-1]);
00469   }
00470   Matrix_Free(permutedIneqs);
00471 
00472 } /* Constraints_fullDimensionize */
00473 
00474 
00475 /**
00476  * Given a matrix that defines a full-dimensional affine lattice, returns the 
00477  * affine sub-lattice spanned in the k first dimensions.
00478  * Useful for instance when you only look for the parameters' validity lattice.
00479  * @param lat the original full-dimensional lattice
00480  * @param subLat the sublattice
00481  */
00482 void Lattice_extractSubLattice(Matrix * lat, unsigned int k, Matrix ** subLat) {
00483   Matrix * H, *Q, *U, *linLat = NULL;
00484   unsigned int i;
00485   dbgStart(Lattice_extractSubLattice);
00486   /* if the dimension is already good, just copy the initial lattice */
00487   if (k==lat->NbRows-1) {
00488     if (*subLat==NULL) {
00489       (*subLat) = Matrix_Copy(lat);
00490     }
00491     else {
00492       Matrix_copySubMatrix(lat, 0, 0, lat->NbRows, lat->NbColumns, (*subLat), 0, 0);
00493     }
00494     return;
00495   }
00496   assert(k<lat->NbRows-1);
00497   /* 1- Make the linear part of the lattice triangular to eliminate terms from 
00498      other dimensions */
00499   Matrix_subMatrix(lat, 0, 0, lat->NbRows, lat->NbColumns-1, &linLat);
00500   /* OPT: any integer column-vector elimination is ok indeed. */
00501   /* OPT: could test if the lattice is already in triangular form. */
00502   left_hermite(linLat, &H, &Q, &U);
00503   if (dbgCompParmMore) {
00504     show_matrix(H);
00505   }
00506   Matrix_Free(Q);
00507   Matrix_Free(U);
00508   Matrix_Free(linLat);
00509   /* if not allocated yet, allocate it */
00510   if (*subLat==NULL) {
00511     (*subLat) = Matrix_Alloc(k+1, k+1);
00512   }
00513   Matrix_copySubMatrix(H, 0, 0, k, k, (*subLat), 0, 0);
00514   Matrix_Free(H);
00515   Matrix_copySubMatrix(lat, 0, lat->NbColumns-1, k, 1, (*subLat), 0, k);
00516   for (i=0; i<k; i++) {
00517     value_set_si((*subLat)->p[k][i], 0);
00518   }
00519   value_set_si((*subLat)->p[k][k], 1);
00520   dbgEnd(Lattice_extractSubLattice);
00521 } /* Lattice_extractSubLattice */
00522 
00523 
00524 /** 
00525  * Computes the overall period of the variables I for (MI) mod |d|, where M is
00526  * a matrix and |d| a vector. Produce a diagonal matrix S = (s_k) where s_k is
00527  * the overall period of i_k 
00528  * @param M the set of affine functions of I (row-vectors)
00529  * @param d the column-vector representing the modulos
00530 */
00531 Matrix * affine_periods(Matrix * M, Matrix * d) {
00532   Matrix * S;
00533   unsigned int i,j;
00534   Value tmp;
00535   Value * periods = (Value *)malloc(sizeof(Value) * M->NbColumns);
00536   value_init(tmp);
00537   for(i=0; i< M->NbColumns; i++) {
00538     value_init(periods[i]);
00539     value_set_si(periods[i], 1);
00540   }
00541   for (i=0; i<M->NbRows; i++) {
00542     for (j=0; j< M->NbColumns; j++) {
00543       value_gcd(tmp, d->p[i][0], M->p[i][j]);
00544       value_divexact(tmp, d->p[i][0], tmp);
00545       value_lcm(periods[j], periods[j], tmp);
00546      }
00547   }
00548   value_clear(tmp);
00549 
00550   /* 2- build S */
00551   S = Matrix_Alloc(M->NbColumns, M->NbColumns);
00552   for (i=0; i< M->NbColumns; i++) 
00553     for (j=0; j< M->NbColumns; j++)
00554       if (i==j) value_assign(S->p[i][j],periods[j]);
00555       else value_set_si(S->p[i][j], 0);
00556 
00557   /* 3- clean up */
00558   for(i=0; i< M->NbColumns; i++) value_clear(periods[i]);
00559   free(periods);
00560   return S;
00561 } /* affine_periods */
00562 
00563 
00564 /** 
00565  * Given an integer matrix B with m rows and integer m-vectors C and d,
00566  * computes the basis of the integer solutions to (BN+C) mod d = 0 (1).
00567  * This is an affine lattice (G): (N 1)^T= G(N' 1)^T, forall N' in Z^b.
00568  * If there is no solution, returns NULL.
00569  * @param B B, a (m x b) matrix
00570  * @param C C, a (m x 1) integer matrix
00571  * @param d d, a (1 x m) integer matrix
00572  * @param imb the affine (b+1)x(b+1) basis of solutions, in the homogeneous
00573  * form. Allocated if initially set to NULL, reused if not.
00574 */
00575 void Equalities_intModBasis(Matrix * B, Matrix * C, Matrix * d, Matrix ** imb) {
00576   int b = B->NbColumns;
00577   /* FIXME: treat the case d=0 as a regular equality B_kN+C_k = 0: */
00578   /* OPT: could keep only equalities for which d>1 */
00579   int nbEqs = B->NbRows;
00580   unsigned int i;
00581 
00582   /* 1- buid the problem DI+BN+C = 0 */
00583   Matrix * eqs = Matrix_Alloc(nbEqs, nbEqs+b+1);
00584   for (i=0; i< nbEqs; i++) {
00585     value_assign(eqs->p[i][i], d->p[0][i]);
00586   }
00587   Matrix_copySubMatrix(B, 0, 0, nbEqs, b, eqs, 0, nbEqs);
00588   Matrix_copySubMatrix(C, 0, 0, nbEqs, 1, eqs, 0, nbEqs+b);
00589 
00590   /* 2- the solution is the validity lattice of the equalities */
00591   Equalities_validityLattice(eqs, nbEqs, imb);
00592   Matrix_Free(eqs);
00593 } /* Equalities_intModBasis */
00594 
00595 
00596 /** kept here for backwards compatiblity. Wrapper to Equalities_intModBasis() */
00597 Matrix * int_mod_basis(Matrix * B, Matrix * C, Matrix * d) {
00598   Matrix * imb = NULL;
00599   Equalities_intModBasis(B, C, d, &imb);
00600   return imb;
00601 } /* int_mod_basis */
00602 
00603 
00604 /**
00605  * Given a parameterized constraints matrix with m equalities, computes the
00606  * compression matrix G such that there is an integer solution in the variables
00607  * space for each value of N', with N = G N' (N are the "nb_parms" parameters)
00608  * @param E a matrix of parametric equalities @param nb_parms the number of
00609  * parameters
00610  * <b>Note: </b>this function is mostly here for backwards
00611  * compatibility. Prefer the use of <tt>Equalities_validityLattice</tt>.
00612 */
00613 Matrix * compress_parms(Matrix * E, int nbParms) {
00614   Matrix * vl=NULL;
00615   Equalities_validityLattice(E, E->NbColumns-2-nbParms, &vl);
00616   return vl;
00617 }/* compress_parms */
00618 
00619 
00620 /** Removes the equalities that involve only parameters, by eliminating some
00621  * parameters in the polyhedron's constraints and in the context.<p> 
00622  * <b>Updates M and Ctxt.</b>
00623  * @param M1 the polyhedron's constraints
00624  * @param Ctxt1 the constraints of the polyhedron's context
00625  * @param renderSpace tells if the returned equalities must be expressed in the
00626  * parameters space (renderSpace=0) or in the combined var/parms space
00627  * (renderSpace = 1)
00628  * @param elimParms the list of parameters that have been removed: an array
00629  * whose 1st element is the number of elements in the list.  (returned)
00630  * @return the system of equalities that involve only parameters.
00631  */
00632 Matrix * Constraints_Remove_parm_eqs(Matrix ** M1, Matrix ** Ctxt1, 
00633                                      int renderSpace, 
00634                                      unsigned int ** elimParms) {
00635   int i, j, k, nbEqsParms =0;
00636   int nbEqsM, nbEqsCtxt, allZeros, nbTautoM = 0, nbTautoCtxt = 0;
00637   Matrix * M = (*M1);
00638   Matrix * Ctxt = (*Ctxt1);
00639   int nbVars = M->NbColumns-Ctxt->NbColumns;
00640   Matrix * Eqs;
00641   Matrix * EqsMTmp;
00642   
00643   /* 1- build the equality matrix(ces) */
00644   nbEqsM = 0;
00645   for (i=0; i< M->NbRows; i++) {
00646     k = First_Non_Zero(M->p[i], M->NbColumns);
00647     /* if it is a tautology, count it as such */
00648     if (k==-1) {
00649       nbTautoM++;
00650     }
00651     else {
00652       /* if it only involves parameters, count it */
00653       if (k>= nbVars+1) nbEqsM++;
00654     }
00655   }
00656 
00657   nbEqsCtxt = 0;
00658   for (i=0; i< Ctxt->NbRows; i++) {
00659     if (value_zero_p(Ctxt->p[i][0])) {
00660       if (First_Non_Zero(Ctxt->p[i], Ctxt->NbColumns)==-1) {
00661         nbTautoCtxt++;
00662       }
00663       else {
00664         nbEqsCtxt ++;
00665       }
00666     }
00667   }
00668   nbEqsParms = nbEqsM + nbEqsCtxt; 
00669 
00670   /* nothing to do in this case */
00671   if (nbEqsParms+nbTautoM+nbTautoCtxt==0) {
00672     (*elimParms) = (unsigned int*) malloc(sizeof(int));
00673     (*elimParms)[0] = 0;
00674     if (renderSpace==0) {
00675       return Matrix_Alloc(0,Ctxt->NbColumns);
00676     }
00677     else {
00678       return Matrix_Alloc(0,M->NbColumns);
00679     }
00680   }
00681   
00682   Eqs= Matrix_Alloc(nbEqsParms, Ctxt->NbColumns);
00683   EqsMTmp= Matrix_Alloc(nbEqsParms, M->NbColumns);
00684   
00685   /* copy equalities from the context */
00686   k = 0;
00687   for (i=0; i< Ctxt->NbRows; i++) {
00688     if (value_zero_p(Ctxt->p[i][0]) 
00689                      && First_Non_Zero(Ctxt->p[i], Ctxt->NbColumns)!=-1) {
00690       Vector_Copy(Ctxt->p[i], Eqs->p[k], Ctxt->NbColumns);
00691       Vector_Copy(Ctxt->p[i]+1, EqsMTmp->p[k]+nbVars+1, 
00692                   Ctxt->NbColumns-1);
00693       k++;
00694     }
00695   }
00696   for (i=0; i< M->NbRows; i++) {
00697     j=First_Non_Zero(M->p[i], M->NbColumns);
00698     /* copy equalities that involve only parameters from M */
00699     if (j>=nbVars+1) {
00700       Vector_Copy(M->p[i]+nbVars+1, Eqs->p[k]+1, Ctxt->NbColumns-1);
00701       Vector_Copy(M->p[i]+nbVars+1, EqsMTmp->p[k]+nbVars+1, 
00702                   Ctxt->NbColumns-1);
00703       /* mark these equalities for removal */
00704       value_set_si(M->p[i][0], 2);
00705       k++;
00706     }
00707     /* mark the all-zero equalities for removal */
00708     if (j==-1) {
00709       value_set_si(M->p[i][0], 2);
00710     }
00711   }
00712 
00713   /* 2- eliminate parameters until all equalities are used or until we find a
00714   contradiction (overconstrained system) */
00715   (*elimParms) = (unsigned int *) malloc((Eqs->NbRows+1) * sizeof(int));
00716   (*elimParms)[0] = 0;
00717   allZeros = 0;
00718   for (i=0; i< Eqs->NbRows; i++) {
00719     /* find a variable that can be eliminated */
00720     k = First_Non_Zero(Eqs->p[i], Eqs->NbColumns);
00721     if (k!=-1) { /* nothing special to do for tautologies */
00722 
00723       /* if there is a contradiction, return empty matrices */
00724       if (k==Eqs->NbColumns-1) {
00725         printf("Contradiction in %dth row of Eqs: ",k);
00726         show_matrix(Eqs);
00727         Matrix_Free(Eqs);
00728         Matrix_Free(EqsMTmp);
00729         (*M1) = Matrix_Alloc(0, M->NbColumns);
00730         Matrix_Free(M);
00731         (*Ctxt1) = Matrix_Alloc(0,Ctxt->NbColumns);
00732         Matrix_Free(Ctxt);
00733         free(*elimParms);
00734         (*elimParms) = (unsigned int *) malloc(sizeof(int));
00735         (*elimParms)[0] = 0;
00736         if (renderSpace==1) {
00737           return Matrix_Alloc(0,(*M1)->NbColumns);
00738         }
00739         else {
00740           return Matrix_Alloc(0,(*Ctxt1)->NbColumns);
00741         }
00742       } 
00743       /* if we have something we can eliminate, do it in 3 places:
00744          Eqs, Ctxt, and M */
00745       else {
00746         k--; /* k is the rank of the variable, now */
00747         (*elimParms)[0]++;
00748         (*elimParms)[(*elimParms[0])]=k;
00749         for (j=0; j< Eqs->NbRows; j++) {
00750           if (i!=j) {
00751             eliminate_var_with_constr(Eqs, i, Eqs, j, k);
00752             eliminate_var_with_constr(EqsMTmp, i, EqsMTmp, j, k+nbVars);
00753           }
00754         }
00755         for (j=0; j< Ctxt->NbRows; j++) {
00756           if (value_notzero_p(Ctxt->p[i][0])) {
00757             eliminate_var_with_constr(Eqs, i, Ctxt, j, k);
00758           }
00759         }
00760         for (j=0; j< M->NbRows; j++) {
00761           if (value_cmp_si(M->p[i][0], 2)) {
00762             eliminate_var_with_constr(EqsMTmp, i, M, j, k+nbVars);
00763           }
00764         }
00765       }
00766     }
00767     /* if (k==-1): count the tautologies in Eqs to remove them later */
00768     else {
00769       allZeros++;
00770     }
00771   }
00772   
00773   /* elimParms may have been overallocated. Now we know how many parms have
00774      been eliminated so we can reallocate the right amount of memory. */
00775   if (!realloc((*elimParms), ((*elimParms)[0]+1)*sizeof(int))) {
00776     fprintf(stderr, "Constraints_Remove_parm_eqs > cannot realloc()");
00777   }
00778 
00779   Matrix_Free(EqsMTmp);
00780 
00781   /* 3- remove the "bad" equalities from the input matrices
00782      and copy the equalities involving only parameters */
00783   EqsMTmp = Matrix_Alloc(M->NbRows-nbEqsM-nbTautoM, M->NbColumns);
00784   k=0;
00785   for (i=0; i< M->NbRows; i++) {
00786     if (value_cmp_si(M->p[i][0], 2)) {
00787       Vector_Copy(M->p[i], EqsMTmp->p[k], M->NbColumns);
00788       k++;
00789     }
00790   }
00791   Matrix_Free(M);
00792   (*M1) = EqsMTmp;
00793   
00794   EqsMTmp = Matrix_Alloc(Ctxt->NbRows-nbEqsCtxt-nbTautoCtxt, Ctxt->NbColumns);
00795   k=0;
00796   for (i=0; i< Ctxt->NbRows; i++) {
00797     if (value_notzero_p(Ctxt->p[i][0])) {
00798       Vector_Copy(Ctxt->p[i], EqsMTmp->p[k], Ctxt->NbColumns);
00799       k++;
00800     }
00801   }
00802   Matrix_Free(Ctxt);
00803   (*Ctxt1) = EqsMTmp;
00804   
00805   if (renderSpace==0) {/* renderSpace=0: equalities in the parameter space */
00806     EqsMTmp = Matrix_Alloc(Eqs->NbRows-allZeros, Eqs->NbColumns);
00807     k=0;
00808     for (i=0; i<Eqs->NbRows; i++) {
00809       if (First_Non_Zero(Eqs->p[i], Eqs->NbColumns)!=-1) {
00810         Vector_Copy(Eqs->p[i], EqsMTmp->p[k], Eqs->NbColumns);
00811         k++;
00812       }
00813     }
00814   }
00815   else {/* renderSpace=1: equalities rendered in the combined space */
00816     EqsMTmp = Matrix_Alloc(Eqs->NbRows-allZeros, (*M1)->NbColumns);
00817     k=0;
00818     for (i=0; i<Eqs->NbRows; i++) {
00819       if (First_Non_Zero(Eqs->p[i], Eqs->NbColumns)!=-1) {
00820         Vector_Copy(Eqs->p[i], &(EqsMTmp->p[k][nbVars]), Eqs->NbColumns);
00821         k++;
00822       }
00823     }
00824   }
00825   Matrix_Free(Eqs);
00826   Eqs = EqsMTmp;
00827 
00828   return Eqs;
00829 } /* Constraints_Remove_parm_eqs */
00830 
00831 
00832 /** Removes equalities involving only parameters, but starting from a
00833  * Polyhedron and its context.
00834  * @param P the polyhedron
00835  * @param C P's context
00836  * @param renderSpace: 0 for the parameter space, =1 for the combined space.
00837  * @maxRays Polylib's usual <i>workspace</i>.
00838  */
00839 Polyhedron * Polyhedron_Remove_parm_eqs(Polyhedron ** P, Polyhedron ** C, 
00840                                         int renderSpace, 
00841                                         unsigned int ** elimParms, 
00842                                         int maxRays) {
00843   Matrix * Eqs;
00844   Polyhedron * Peqs;
00845   Matrix * M = Polyhedron2Constraints((*P));
00846   Matrix * Ct = Polyhedron2Constraints((*C));
00847 
00848   /* if the Minkowski representation is not computed yet, do not compute it in
00849      Constraints2Polyhedron */
00850   if (F_ISSET((*P), POL_VALID | POL_INEQUALITIES) && 
00851       (F_ISSET((*C), POL_VALID | POL_INEQUALITIES))) {
00852     FL_INIT(maxRays, POL_NO_DUAL);
00853   }
00854     
00855   Eqs = Constraints_Remove_parm_eqs(&M, &Ct, renderSpace, elimParms);
00856   Peqs = Constraints2Polyhedron(Eqs, maxRays);
00857   Matrix_Free(Eqs);
00858 
00859   /* particular case: no equality involving only parms is found */
00860   if (Eqs->NbRows==0) {
00861     Matrix_Free(M);
00862     Matrix_Free(Ct);
00863     return Peqs;
00864   }
00865   Polyhedron_Free(*P);
00866   Polyhedron_Free(*C);
00867   (*P) = Constraints2Polyhedron(M, maxRays);
00868   (*C) = Constraints2Polyhedron(Ct, maxRays);
00869   Matrix_Free(M);
00870   Matrix_Free(Ct);
00871   return Peqs;
00872 } /* Polyhedron_Remove_parm_eqs */
00873 
00874 
00875 /**
00876  * Given a matrix with m parameterized equations, compress the nb_parms
00877  * parameters and n-m variables so that m variables are integer, and transform
00878  * the variable space into a n-m space by eliminating the m variables (using
00879  * the equalities) the variables to be eliminated are chosen automatically by
00880  * the function.
00881  * <b>Deprecated.</b> Try to use Constraints_fullDimensionize instead.
00882  * @param M the constraints 
00883  * @param the number of parameters
00884  * @param validityLattice the the integer lattice underlying the integer
00885  * solutions.
00886 */
00887 Matrix * full_dimensionize(Matrix const * M, int nbParms, 
00888                            Matrix ** validityLattice) {
00889   Matrix * Eqs, * Ineqs;
00890   Matrix * permutedEqs, * permutedIneqs;
00891   Matrix * Full_Dim;
00892   Matrix * WVL; /* The Whole Validity Lattice (vars+parms) */
00893   unsigned int i,j;
00894   int nbElimVars;
00895   unsigned int * permutation, * permutationInv;
00896   /* 0- Split the equalities and inequalities from each other */
00897   split_constraints(M, &Eqs, &Ineqs);
00898 
00899   /* 1- if the polyhedron is already full-dimensional, return it */
00900   if (Eqs->NbRows==0) {
00901     Matrix_Free(Eqs);
00902     (*validityLattice) = Identity_Matrix(nbParms+1);
00903     return Ineqs;
00904   }
00905   nbElimVars = Eqs->NbRows;
00906 
00907   /* 2- put the vars to be eliminated at the first positions, 
00908      and compress the other vars/parms
00909      -> [ variables to eliminate / parameters / variables to keep ] */
00910   permutation = find_a_permutation(Eqs, nbParms);
00911   if (dbgCompParm) {
00912     printf("Permuting the vars/parms this way: [ ");
00913     for (i=0; i< Eqs->NbColumns; i++) {
00914       printf("%d ", permutation[i]);
00915     }
00916     printf("]\n");
00917   }
00918   permutedEqs = mpolyhedron_permute(Eqs, permutation);
00919   WVL = compress_parms(permutedEqs, Eqs->NbColumns-2-Eqs->NbRows);
00920   if (dbgCompParm) {
00921     printf("Whole validity lattice: ");
00922     show_matrix(WVL);
00923   }
00924   mpolyhedron_compress_last_vars(permutedEqs, WVL);
00925   permutedIneqs = mpolyhedron_permute(Ineqs, permutation);
00926   if (dbgCompParm) {
00927     show_matrix(permutedEqs);
00928   }
00929   Matrix_Free(Eqs);
00930   Matrix_Free(Ineqs);
00931   mpolyhedron_compress_last_vars(permutedIneqs, WVL);
00932   if (dbgCompParm) {
00933     printf("After compression: ");
00934     show_matrix(permutedIneqs);
00935   }
00936   /* 3- eliminate the first variables */
00937   if (!mpolyhedron_eliminate_first_variables(permutedEqs, permutedIneqs)) {
00938     fprintf(stderr,"full-dimensionize > variable elimination failed. \n"); 
00939     return NULL;
00940   }
00941   if (dbgCompParm) {
00942     printf("After elimination of the variables: ");
00943     show_matrix(permutedIneqs);
00944   }
00945 
00946   /* 4- get rid of the first (zero) columns, 
00947      which are now useless, and put the parameters back at the end */
00948   Full_Dim = Matrix_Alloc(permutedIneqs->NbRows,
00949                           permutedIneqs->NbColumns-nbElimVars);
00950   for (i=0; i< permutedIneqs->NbRows; i++) {
00951     value_set_si(Full_Dim->p[i][0], 1);
00952     for (j=0; j< nbParms; j++) 
00953       value_assign(Full_Dim->p[i][j+Full_Dim->NbColumns-nbParms-1], 
00954                    permutedIneqs->p[i][j+nbElimVars+1]);
00955     for (j=0; j< permutedIneqs->NbColumns-nbParms-2-nbElimVars; j++) 
00956       value_assign(Full_Dim->p[i][j+1], 
00957                    permutedIneqs->p[i][nbElimVars+nbParms+j+1]);
00958     value_assign(Full_Dim->p[i][Full_Dim->NbColumns-1], 
00959                  permutedIneqs->p[i][permutedIneqs->NbColumns-1]);
00960   }
00961   Matrix_Free(permutedIneqs);
00962   
00963   /* 5- Keep only the the validity lattice restricted to the parameters */
00964   *validityLattice = Matrix_Alloc(nbParms+1, nbParms+1);
00965   for (i=0; i< nbParms; i++) {
00966     for (j=0; j< nbParms; j++)
00967       value_assign((*validityLattice)->p[i][j], 
00968                    WVL->p[i][j]);
00969     value_assign((*validityLattice)->p[i][nbParms], 
00970                  WVL->p[i][WVL->NbColumns-1]);
00971   }
00972   for (j=0; j< nbParms; j++) 
00973     value_set_si((*validityLattice)->p[nbParms][j], 0);
00974   value_assign((*validityLattice)->p[nbParms][nbParms], 
00975                WVL->p[WVL->NbColumns-1][WVL->NbColumns-1]);
00976 
00977   /* 6- Clean up */
00978   Matrix_Free(WVL);
00979   return Full_Dim;
00980 } /* full_dimensionize */
00981 
00982 #undef dbgCompParm
00983 #undef dbgCompParmMore

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