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

compress_parms.c

Go to the documentation of this file.
00001 // (C) B. Meister 12/2003-2005
00002 // LSIIT -ICPS 
00003 // UMR 7005 CNRS
00004 // Louis Pasteur University (ULP), Strasbourg, France
00005 
00006 #include <polylib/polylib.h>
00007 
00008 
00009 // custom modulo, following the definition:
00010 // a mod b = c equivalent_to: there exists a bigest nonnegative n in Z such that a = n|b| + c
00011 // (the % C operator is not trusted here)
00012 void b_modulo(Value g, Value a, Value b) {
00013   if (value_zero_p(b)) {value_assign(g,a); return;}
00014   if (value_neg_p(b)) value_oppose(b, b);
00015   if (value_posz_p(a)) {
00016     value_division(g,a,b);
00017     value_multiply(g,g,b);
00018     value_subtract(g, a, g);
00019     return;
00020   }
00021   else {
00022     value_division(g, a, b);
00023     value_decrement(g, g);
00024     value_multiply(g, g, b);
00025     value_subtract(g, a, g);
00026     //  g =a-(a/b-1)*b;
00027     if (value_eq(g,b)) value_set_si(g,0);
00028   }
00029 }
00030 
00031 // given a full-row-rank nxm matrix M(made of m row-vectors), 
00032 // computes the basis K (made of n-m column-vectors) of the integer kernel of the rows of M
00033 // so we have: M.K = 0
00034 // assumed: M is full row rank
00035 Matrix * int_ker(Matrix * M) {
00036   Matrix *U, *Q, *H, *K;
00037   int i, j;
00038 
00039   // there is a non-null kernel if and only if the dimension m of the space spanned by the rows 
00040   // is inferior to the number n of variables 
00041   if (M->NbColumns <= M->NbRows) {
00042     K = Matrix_Alloc(M->NbColumns, 0);
00043     return K;
00044   }
00045   // computes MU = [H 0]
00046   left_hermite(M, &H, &Q, &U);
00047   
00048   // obviously, the Integer Kernel is made of the last n-m columns of U
00049   K=Matrix_Alloc(U->NbRows, U->NbColumns-M->NbRows);
00050   for (i=0; i< U->NbRows; i++)
00051     for(j=0; j< U->NbColumns-M->NbRows; j++) 
00052       value_assign(K->p[i][j], U->p[i][M->NbRows+j]);
00053 
00054   // clean up
00055   Matrix_Free(H);
00056   Matrix_Free(U);
00057   Matrix_Free(Q);
00058   return K;
00059 }
00060 
00061 // Compute the overall period of the variables I for (MI) mod |d|,
00062 // where M is a matrix and |d| a vector
00063 // Produce a diagonal matrix S = (s_k) where s_k is the overall period of i_k
00064 Matrix * affine_periods(Matrix * M, Matrix * d) {
00065   Matrix * S;
00066   unsigned int i,j;
00067   Value tmp;
00068   value_init(tmp);
00069   // 1- compute the overall periods
00070   Value * periods = (Value *)malloc(sizeof(Value) * M->NbColumns);
00071   for(i=0; i< M->NbColumns; i++) {
00072     value_init(periods[i]);
00073     value_set_si(periods[i], 1);
00074   }
00075   for (i=0; i<M->NbRows; i++) {
00076     for (j=0; j< M->NbColumns; j++) {
00077       Gcd(d->p[i][0], M->p[i][j], &tmp);
00078       value_division(tmp, d->p[i][0], tmp);
00079       Lcm3(periods[j], tmp, &(periods[j]));
00080      }
00081   }
00082   value_clear(tmp);
00083 
00084   // 2- build S
00085   S = Matrix_Alloc(M->NbColumns, M->NbColumns);
00086   for (i=0; i< M->NbColumns; i++) 
00087     for (j=0; j< M->NbColumns; j++)
00088       if (i==j) value_assign(S->p[i][j],periods[j]);
00089       else value_set_si(S->p[i][j], 0);
00090 
00091   // 3- clean up
00092   for(i=0; i< M->NbColumns; i++) value_clear(periods[i]);
00093   free(periods);
00094   return S;
00095 }
00096 
00097 // given a matrix B' with m rows and m-vectors C' and d, computes the 
00098 // basis of the integer solutions to (B'N+C') mod d = 0 (1).
00099 // the Ns verifying the system B'N+C' = 0 are solutions of (1)
00100 // K is a basis of the integer kernel of B: its column-vectors link two solutions of (1)
00101 // Moreover, B'_iN mod d is periodic of period (s_ik):
00102 // B'N mod d is periodic of period (s_k) = lcm_i(s_ik)
00103 // The linear part of G is given by the HNF of (K | S), where S is the full-dimensional diagonal matrix (s_k)
00104 // the constant part of G is a particular solution of (1)
00105 // if no integer constant part is found, there is no solution and this function returns NULL.
00106 Matrix * int_mod_basis(Matrix * Bp, Matrix * Cp, Matrix * d) {
00107   int nb_eqs = Bp->NbRows;
00108   unsigned int nb_parms=Bp->NbColumns;
00109   unsigned int i, j;
00110 
00111   //   a/ compute K and S
00112   Matrix * K = int_ker(Bp);
00113   Matrix * S = affine_periods(Bp, d);
00114 
00115   // show_matrix(K);
00116   // show_matrix(S);
00117   
00118   //   b/ compute the linear part of G : HNF(K|S)
00119 
00120   // fill K|S
00121   Matrix * KS = Matrix_Alloc(nb_parms, K->NbColumns+ nb_parms);
00122   for(i=0; i< KS->NbRows; i++) {
00123     for(j=0; j< K->NbColumns; j++) value_assign(KS->p[i][j], K->p[i][j]);
00124     for(j=0; j< S->NbColumns; j++) value_assign(KS->p[i][j+K->NbColumns], S->p[i][j]);
00125   }
00126 
00127   // show_matrix(KS);
00128 
00129   // HNF(K|S)
00130   Matrix * H, * U, * Q;
00131   left_hermite(KS, &H, &U, &Q);
00132   Matrix_Free(KS);
00133   Matrix_Free(U);
00134   Matrix_Free(Q);
00135   
00136   // printf("HNF(K|S) = ");show_matrix(H);
00137 
00138   // put HNF(K|S) in the p x p matrix S (which has already the appropriate size so we economize a Matrix_Alloc)
00139   for (i=0; i< nb_parms; i++) {
00140     for (j=0; j< nb_parms; j++) 
00141       value_assign(S->p[i][j], H->p[i][j]);
00142 
00143 
00144   }
00145   Matrix_Free(H);
00146 
00147   //   c/ compute U_M.N'_0 = N_0: 
00148   Matrix * M = Matrix_Alloc(nb_eqs, nb_parms+nb_eqs);
00149   // N'_0 = M_H^{-1}.(-C' mod d), which must be integer (à verifier encore un peu mais je crois que c'est ok)
00150   // and where H_M = HNF(M) with M = (B' D) : M.U_M = [H_M 0]
00151 
00152   //       copy the B' part
00153   for (i=0; i< nb_eqs; i++) {
00154     for (j=0; j< nb_parms; j++) {
00155       value_assign(M->p[i][j], Bp->p[i][j]);
00156     }
00157   //       copy the D part
00158     for (j=0; j< nb_eqs; j++) {
00159       if (i==j) value_assign(M->p[i][j+nb_parms], d->p[i][0]);
00160       else value_set_si(M->p[i][j+nb_parms], 0);
00161     }
00162   }
00163   
00164   //       compute inv_H_M, the inverse of the HNF H of M = (B' D)
00165   left_hermite(M, &H, &Q, &U);
00166   Matrix * inv_H_M=Matrix_Alloc(nb_eqs, nb_eqs+1);
00167   //again, do a square Matrix from H, using the non-used Matrix Ha
00168   Matrix * Ha = Matrix_Alloc(nb_eqs, nb_eqs);
00169   for(i=0; i< nb_eqs; i++) {
00170     for(j=0; j< nb_eqs; j++) {
00171       value_assign(Ha->p[i][j], H->p[i][j]); 
00172     }
00173   }
00174   MatInverse(Ha, inv_H_M);
00175   Matrix_Free(Ha);
00176   Matrix_Free(H);
00177   Matrix_Free(Q); // only inv_H_M and U_M (alias U) are needed
00178 
00179   //       compute (-C') mod d
00180   Matrix * Cp_mod_d = Matrix_Alloc(nb_eqs, 1);
00181   for (i=0; i< nb_eqs; i++) {
00182     value_oppose(Cp->p[i][0], Cp->p[i][0]);
00183     b_modulo(Cp_mod_d->p[i][0], Cp->p[i][0], d->p[i][0]);
00184     //    Cp_mod_d->p[i][0] = b_modulo(-Cp->p[i][0], d->p[i][0]);
00185   }
00186 
00187   // Compute N'_0 = inv_H_M.((-C') mod d)
00188   //  actually compute (N' \\ 0) sot hat N = U^{-1}.(N' \\ 0)
00189   Matrix * Np_0= Matrix_Alloc(U->NbColumns, 1);
00190   for(i=0; i< nb_eqs; i++) 
00191     {
00192       value_set_si(Np_0->p[i][0], 0);
00193       for(j=0; j< nb_eqs; j++) {
00194         value_addmul(Np_0->p[i][0], inv_H_M->p[i][j], Cp_mod_d->p[j][0]);
00195       }
00196     }
00197   for(i=nb_eqs; i< U->NbColumns; i++) value_set_si(Np_0->p[i][0], 0);
00198   
00199 
00200   // it is still needed to divide the rows of N'_0 by the common denominator of the rows of H_M
00201   // if these rows are not divisible, there is no integer N'_0 so return NULL
00202   for (i=0; i< nb_eqs; i++) {
00203 
00204 #ifdef GNUMP
00205     if (mpz_divisible_p(Np_0->p[i][0], inv_H_M->p[i][nb_eqs])) mpz_divexact(Np_0->p[i][0], Np_0->p[i][0], inv_H_M->p[i][nb_eqs]);
00206 #else
00207     if (!(Np_0->p[i][0]%inv_H_M->p[i][nb_eqs])) Np_0->p[i][0]/=inv_H_M->p[i][nb_eqs];
00208 #endif
00209     else {
00210       Matrix_Free(S);
00211       Matrix_Free(inv_H_M);
00212       Matrix_Free(Np_0);
00213       fprintf(stderr, "int_mod_basis > No particular solution: polyhedron without integer points.\n");
00214       return NULL;
00215     }
00216   }
00217   // show_matrix(Np_0);
00218 
00219   // now compute the actual particular solution N_0 = U_M. N'_0
00220   Matrix * N_0 = Matrix_Alloc(U->NbColumns, 1);
00221   // OPT: seules les nb_eq premières valeurs de N_0 sont utiles en fait.
00222   Matrix_Product(U, Np_0, N_0);
00223   // show_matrix(N_0);
00224   Matrix_Free(Np_0);
00225   Matrix_Free(U);
00226 
00227   // build the whole compression matrix: 
00228   Matrix * G = Matrix_Alloc(S->NbRows+1, S->NbRows+1);
00229   for (i=0; i< S->NbRows; i++) {
00230     for(j=0; j< S->NbRows; j++) 
00231       value_assign(G->p[i][j], S->p[i][j]);
00232     value_assign(G->p[i][S->NbRows], N_0->p[i][0]);
00233   }
00234 
00235   for (j=0; j< S->NbRows; j++) value_set_si(G->p[S->NbRows][j],0);
00236   value_set_si(G->p[S->NbRows][S->NbRows],1);
00237 
00238   // clean up
00239   Matrix_Free(S);
00240   Matrix_Free(N_0);
00241   return G;
00242 } // int_mod_basis
00243 
00244 
00245 // utiliy function: given a matrix of constraints AI+BN+C, extract the part of A, corresponding to the variables.
00246 Matrix * get_linear_part(Matrix const * E, int nb_parms){ 
00247   unsigned int i,j, k, nb_eqs=E->NbRows;
00248   int nb_vars=E->NbColumns - nb_parms -2;
00249   Matrix * A = Matrix_Alloc(nb_eqs, nb_vars);
00250   k=0;
00251   for (i=0; i< E->NbRows; i++) {
00252     for (j=0; j< nb_vars; j++)
00253       value_assign(A->p[i][j],E->p[i][j+1]);
00254   }
00255   return A;
00256 }
00257 
00258 // utiliy function: given a matrix of constraints AI = -BN-C, extract the part of -B, corresponding to the parametrs.
00259 Matrix * get_parameter_part(Matrix const * E, int nb_parms) {
00260   unsigned int i,j, k, nb_eqs=E->NbRows;
00261   int nb_vars=E->NbColumns - nb_parms -2;
00262   Matrix * B = Matrix_Alloc(nb_eqs,nb_parms);
00263   k=0;
00264   // well a priori the minus is unuseful here, but
00265   for(i=0; i< E->NbRows; i++) {
00266     for(j=0; j< nb_parms; j++)
00267       value_assign(B->p[i][j], E->p[i][1+nb_vars+j]);
00268   }
00269   return B;
00270 }
00271 
00272 // utiliy function: given a matrix of constraints (0 or 1)AI= -BN-C, extract the part of -C corresponding to the equalities.
00273 Matrix * get_constant_part(Matrix const * E, int nb_parms) {
00274   unsigned int i,j, k, nb_eqs=E->NbRows;
00275   int nb_vars=E->NbColumns - nb_parms -2;
00276   Matrix * C = Matrix_Alloc(nb_eqs,1);
00277   k=0;
00278   for(i=0; i< E->NbRows; i++) {
00279     value_assign(C->p[i][0], E->p[i][E->NbColumns-1]);
00280   }
00281   return C;
00282 }
00283 
00284 // utility function: given a matrix containing the equations AI+BN+C=0, compute the HNF of A : A = [Ha 0].Q and return :  
00285 // . B'= H^-1.(-B) 
00286 // . C'= H^-1.(-C)
00287 // . U = Q^-1 (-> return value)
00288 // . D, where Ha^-1 = D^-1.H^-1 with H and D integer matrices 
00289 // in fact, as D is diagonal, we return d, a column-vector
00290 Matrix * extract_funny_stuff(Matrix const * E, int nb_parms, Matrix ** Bp, Matrix **Cp, Matrix **d) {
00291 unsigned int i,j, k, nb_eqs=E->NbRows;
00292   int nb_vars=E->NbColumns - nb_parms -2;
00293   Matrix * A, * Ap, * Ha, * U, * Q, * H, *B, *C;
00294 
00295   // particular case: 
00296   if (nb_eqs==0) {
00297     *Bp = Matrix_Alloc(0, E->NbColumns);
00298     *Cp = Matrix_Alloc(0, E->NbColumns);
00299     *d = NULL;
00300     return NULL;
00301   }
00302 
00303   // 1- build A
00304   A = get_linear_part(E, nb_parms);
00305   // show_matrix(A);
00306 
00307   // 2- Compute Ha^-1, where Ha is the left HNF of A
00308   //   a/ Compute H = [Ha 0]
00309   left_hermite(A, &H, &Q, &U);
00310   Matrix_Free(A);
00311   Matrix_Free(Q);
00312   
00313   //   b/ just keep the m x m matrix Ha
00314   Ha = Matrix_Alloc(nb_eqs, nb_eqs);
00315   for (i=0; i< nb_eqs; i++) {
00316     for (j=0; j< nb_eqs; j++) {
00317       value_assign(Ha->p[i][j],H->p[i][j]);
00318     }
00319   }
00320   Matrix_Free(H);
00321 
00322   // show_matrix(Ha);
00323 
00324   //  c/ Invert Ha
00325   Matrix * Ha_pre_inv = Matrix_Alloc(nb_eqs, nb_eqs+1);
00326   if(!MatInverse(Ha, Ha_pre_inv)) { 
00327     fprintf(stderr,"extract_funny_stuff > Matrix Ha is non-invertible.");
00328   }
00329   // store back Ha^-1  in Ha, to spare a MatrixAlloc/MatrixFree
00330   for(i=0; i< nb_eqs; i++) {
00331     for(j=0; j< nb_eqs; j++) {
00332       value_assign(Ha->p[i][j],Ha_pre_inv->p[i][j]);
00333     }
00334   }
00335   // show_matrix(Ha_pre_inv);
00336   // the diagonal elements of D are stored in the last column of Ha_pre_inv (property of MatInverse).
00337   (*d) = Matrix_Alloc(Ha_pre_inv->NbRows, 1);
00338 
00339   for (i=0; i< Ha_pre_inv->NbRows; i++) value_assign((*d)->p[i][0], Ha_pre_inv->p[i][Ha_pre_inv->NbColumns-1]);
00340  
00341   // show_matrix(Ha);
00342   // show_matrix(*d);
00343   Matrix_Free(Ha_pre_inv);
00344 
00345   // 3- Build B'and C'
00346   // compute B'
00347   B = get_parameter_part(E, nb_parms);
00348   (*Bp) = Matrix_Alloc(B->NbRows,B->NbColumns);
00349   // show_matrix(B);
00350   Matrix_Product(Ha, B, (*Bp));
00351   // show_matrix(*Bp);
00352   Matrix_Free(B);
00353   
00354   // compute C'
00355   C = get_constant_part(E, nb_parms);
00356   // show_matrix(C);
00357   (*Cp) = Matrix_Alloc(nb_eqs, 1);
00358   Matrix_Product(Ha, C, (*Cp));
00359   // show_matrix(*Cp);
00360   Matrix_Free(C);
00361 
00362   Matrix_Free(Ha);
00363   return U;
00364 }
00365   
00366 
00367 // given a parameterized constraints matrix with m equalities, computes the compression matrix 
00368 // G such that there is an integer solution in the variables space for each value of N', with 
00369 // N = G N' (N are the "nb_parms" parameters)
00370 Matrix * compress_parms(Matrix * E, int nb_parms) {
00371   unsigned int i,j, k, nb_eqs=0;
00372   int nb_vars=E->NbColumns - nb_parms -2;
00373   Matrix *U, *d, *Bp, *Cp;
00374 
00375   // particular case where there is no equation
00376   if (E->NbRows==0) return Identity_Matrix(nb_parms+1);
00377 
00378   U = extract_funny_stuff(E, nb_parms, &Bp, & Cp, &d); 
00379 
00380   Matrix_Free(U);
00381   // The compression matrix N = G.N' must be such that (B'N+C') mod d = 0 (1)
00382 
00383   // the Ns verifying the system B'N+C' = 0 are solutions of (1)
00384   // K is a basis of the integer kernel of B: its column-vectors link two solutions of (1)
00385   // Moreover, B'_iN mod d is periodic of period (s_ik):
00386   // B'N mod d is periodic of period (s_k) = lcm_i(s_ik)
00387   // The linear part of G is given by the HNF of (K | S), where S is the full-dimensional diagonal matrix (s_k)
00388   // the constant part of G is a particular solution of (1)
00389   // if no integer constant part is found, there is no solution.
00390 
00391   Matrix * G = int_mod_basis(Bp, Cp, d);
00392   Matrix_Free(Bp);
00393   Matrix_Free(Cp);
00394   Matrix_Free(d);
00395   return G;
00396 }
00397 
00398 // given a matrix with m parameterized equations, compress the nb_parms parameters and n-m variables so that m variables are integer,
00399 // and transform the variable space into a n-m space by eliminating the m variables (using the equalities)
00400 // the variables to be eliminated are chosen automatically by the function
00401 Matrix * full_dimensionize(Matrix const * M, int nb_parms, Matrix ** Validity_Lattice) {
00402   Matrix * Eqs, * Ineqs;
00403   Matrix * Permuted_Eqs, * Permuted_Ineqs;
00404   Matrix * Full_Dim;
00405   Matrix * Whole_Validity_Lattice;
00406   unsigned int i,j;
00407   int nb_elim_vars;
00408   int * permutation, * permutation_inv;
00409   // 0- Split the equalities and inequalities from each other
00410   split_constraints(M, &Eqs, &Ineqs);
00411 
00412   // 1- if the polyhedron is already full-dimensional, return it
00413   if (Eqs->NbRows==0) {
00414     Matrix_Free(Eqs);
00415     (*Validity_Lattice) = Identity_Matrix(nb_parms+1);
00416     return Ineqs;
00417   }
00418   nb_elim_vars = Eqs->NbRows;
00419 
00420   // 2- put the vars to be eliminated at the first positions, and compress the other vars/parms
00421   // -> [ variables to eliminate / parameters / variables to keep ]
00422   permutation = find_a_permutation(Eqs, nb_parms);
00423   Permuted_Eqs = mpolyhedron_permute(Eqs, permutation);
00424   Whole_Validity_Lattice = compress_parms(Permuted_Eqs, Eqs->NbColumns-2-Eqs->NbRows);
00425   mpolyhedron_compress_last_vars(Permuted_Eqs, Whole_Validity_Lattice);
00426   Permuted_Ineqs = mpolyhedron_permute(Ineqs, permutation);
00427   Matrix_Free(Eqs);
00428   Matrix_Free(Ineqs);
00429   mpolyhedron_compress_last_vars(Permuted_Ineqs, Whole_Validity_Lattice);
00430   // show_matrix(Whole_Validity_Lattice);
00431 
00432   // 3- eliminate the first variables
00433   if (!mpolyhedron_eliminate_first_variables(Permuted_Eqs, Permuted_Ineqs)) {
00434     fprintf(stderr,"full-dimensionize > variable elimination failed. \n"); 
00435     return NULL;
00436   }
00437   // show_matrix(Permuted_Eqs);
00438   // show_matrix(Permuted_Ineqs);
00439 
00440   // 4- get rid of the first (zero) columns, which are now useless, and put the parameters back at the end
00441   Full_Dim = Matrix_Alloc(Permuted_Ineqs->NbRows, Permuted_Ineqs->NbColumns-nb_elim_vars);
00442   for (i=0; i< Permuted_Ineqs->NbRows; i++) {
00443     value_set_si(Full_Dim->p[i][0], 1);
00444     for (j=0; j< nb_parms; j++) 
00445       value_assign(Full_Dim->p[i][j+Full_Dim->NbColumns-nb_parms-1], Permuted_Ineqs->p[i][j+nb_elim_vars+1]);
00446     for (j=0; j< Permuted_Ineqs->NbColumns-nb_parms-2-nb_elim_vars; j++) 
00447       value_assign(Full_Dim->p[i][j+1], Permuted_Ineqs->p[i][nb_elim_vars+nb_parms+j+1]);
00448     value_assign(Full_Dim->p[i][Full_Dim->NbColumns-1], Permuted_Ineqs->p[i][Permuted_Ineqs->NbColumns-1]);
00449   }
00450   Matrix_Free(Permuted_Ineqs);
00451   // show_matrix(Full_Dim);
00452   // return Full_Dim; 
00453 
00454   // 4- Un-permute (so that the parameters are at the end as usual)
00455   // -> [variables (mixed) / parameters ]
00456   /*permutation_inv = permutation_inverse(permutation, Permuted_Ineqs->NbColumns-1);
00457   free(permutation);
00458   Ineqs = mpolyhedron_permute(Permuted_Ineqs, permutation_inv);
00459   mpolyhedron_simplify(Ineqs);
00460   free(permutation_inv); */
00461   
00462   // 5- Keep only the the validity lattice restricted to the parameters
00463   *Validity_Lattice = Matrix_Alloc(nb_parms+1, nb_parms+1);
00464   for (i=0; i< nb_parms; i++) {
00465     for (j=0; j< nb_parms; j++)
00466       value_assign((*Validity_Lattice)->p[i][j], Whole_Validity_Lattice->p[i][j]);
00467     value_assign((*Validity_Lattice)->p[i][nb_parms], Whole_Validity_Lattice->p[i][Whole_Validity_Lattice->NbColumns-1]);
00468   }
00469   for (j=0; j< nb_parms; j++) value_set_si((*Validity_Lattice)->p[nb_parms][j], 0);
00470   value_assign((*Validity_Lattice)->p[nb_parms][nb_parms], Whole_Validity_Lattice->p[Whole_Validity_Lattice->NbColumns-1][Whole_Validity_Lattice->NbColumns-1]);
00471 
00472   // 6- Clean up 
00473   Matrix_Free(Whole_Validity_Lattice);
00474   return Full_Dim;
00475   //  return Ineqs;
00476 } // full_dimensionize
00477 

Generated on Mon Sep 12 15:15:10 2005 for polylib by doxygen 1.3.5