polylib 5.22.8
compress_parms.c
Go to the documentation of this file.
1/**
2 * $Id: compress_parms.c,v 1.32 2006/11/03 17:34:26 skimo Exp $
3 *
4 * The integer points in a parametric linear subspace of Q^n are generally
5 * lying on a sub-lattice of Z^n.
6 * Functions here use and compute validity lattices, i.e. lattices induced on a
7 * set of variables by such equalities involving another set of integer
8 * variables.
9 * @author B. Meister 12/2003-2006 meister@icps.u-strasbg.fr
10 * LSIIT -ICPS
11 * UMR 7005 CNRS
12 * Louis Pasteur University (ULP), Strasbourg, France
13 */
14
15#include <polylib/polylib.h>
16#include <stdlib.h>
17
18/**
19 * debug flags (2 levels)
20 */
21#define dbgCompParm 0
22#define dbgCompParmMore 0
23
24#define dbgStart(a) \
25 if (dbgCompParmMore) { \
26 printf(" -- begin "); \
27 printf(#a); \
28 printf(" --\n"); \
29 } \
30 while (0)
31#define dbgEnd(a) \
32 if (dbgCompParmMore) { \
33 printf(" -- end "); \
34 printf(#a); \
35 printf(" --\n"); \
36 } \
37 while (0)
38
39/**
40 * Given a full-row-rank nxm matrix M made of m row-vectors), computes the
41 * basis K (made of n-m column-vectors) of the integer kernel of the rows of M
42 * so we have: M.K = 0
43 */
45 Matrix *U, *Q, *H, *H2, *K = NULL;
46 int i, j, rk;
47
48 if (dbgCompParm)
49 show_matrix(M);
50 /* eliminate redundant rows : UM = H*/
51 right_hermite(M, &H, &Q, &U);
52 for (rk = H->NbRows - 1; (rk >= 0) && Vector_IsZero(H->p[rk], H->NbColumns);
53 rk--)
54 ;
55 rk++;
56 if (dbgCompParmMore) {
57 printf("rank = %d\n", rk);
58 }
59
60 /* there is a non-null kernel if and only if the dimension m of
61 the space spanned by the rows
62 is inferior to the number n of variables */
63 if (M->NbColumns <= rk) {
64 Matrix_Free(H);
65 Matrix_Free(Q);
66 Matrix_Free(U);
67 K = Matrix_Alloc(M->NbColumns, 0);
68 return K;
69 }
70 Matrix_Free(U);
71 Matrix_Free(Q);
72 /* fool left_hermite by giving NbRows =rank of M*/
73 H->NbRows = rk;
74 /* computes MU = [H 0] */
75 left_hermite(H, &H2, &Q, &U);
76 if (dbgCompParmMore) {
77 printf("-- Int. Kernel -- \n");
78 show_matrix(M);
79 printf(" = \n");
80 show_matrix(H2);
81 show_matrix(U);
82 }
83 H->NbRows = M->NbRows;
84 Matrix_Free(H);
85 /* the Integer Kernel is made of the last n-rk columns of U */
86 Matrix_subMatrix(U, 0, rk, U->NbRows, U->NbColumns, &K);
87
88 /* clean up */
89 Matrix_Free(H2);
90 Matrix_Free(U);
91 Matrix_Free(Q);
92 return K;
93} /* int_ker */
94
95/**
96 * Computes the intersection of two linear lattices, whose base vectors are
97 * respectively represented in A and B.
98 * If I and/or Lb is set to NULL, then the matrix is allocated.
99 * Else, the matrix is assumed to be allocated already.
100 * I and Lb are rk x rk, where rk is the rank of A (or B).
101 * @param A the full-row rank matrix whose column-vectors are the basis for the
102 * first linear lattice.
103 * @param B the matrix whose column-vectors are the basis for the second linear
104 * lattice.
105 * @param Lb the matrix such that B.Lb = I, where I is the intersection.
106 * @return their intersection.
107 */
108static void linearInter(Matrix *A, Matrix *B, Matrix **I, Matrix **Lb) {
109 Matrix *AB = NULL;
110 int rk = A->NbRows;
111 int a = A->NbColumns;
112 int b = B->NbColumns;
113 int i, j, z = 0;
114
115 Matrix *H, *U, *Q;
116 /* ensure that the spanning vectors are in the same space */
117 assert(B->NbRows == rk);
118 /* 1- build the matrix
119 * (A 0 1)
120 * (0 B 1)
121 */
122 AB = Matrix_Alloc(2 * rk, a + b + rk);
123 Matrix_copySubMatrix(A, 0, 0, rk, a, AB, 0, 0);
124 Matrix_copySubMatrix(B, 0, 0, rk, b, AB, rk, a);
125 for (i = 0; i < rk; i++) {
126 value_set_si(AB->p[i][a + b + i], 1);
127 value_set_si(AB->p[i + rk][a + b + i], 1);
128 }
129 if (dbgCompParm) {
130 show_matrix(AB);
131 }
132
133 /* 2- Compute its left Hermite normal form. AB.U = [H 0] */
134 left_hermite(AB, &H, &Q, &U);
135 Matrix_Free(AB);
136 Matrix_Free(Q);
137 /* count the number of non-zero colums in H */
138 for (z = H->NbColumns - 1; value_zero_p(H->p[H->NbRows - 1][z]); z--)
139 ;
140 z++;
141 if (dbgCompParm) {
142 show_matrix(H);
143 printf("z=%d\n", z);
144 }
145 Matrix_Free(H);
146 /* if you split U in 9 submatrices, you have:
147 * A.U_13 = -U_33
148 * B.U_23 = -U_33,
149 * where the nb of cols of U_{*3} equals the nb of zero-cols of H
150 * U_33 is a (the smallest) combination of col-vectors of A and B at the same
151 * time: their intersection.
152 */
153 Matrix_subMatrix(U, a + b, z, U->NbColumns, U->NbColumns, I);
154 Matrix_subMatrix(U, a, z, a + b, U->NbColumns, Lb);
155 if (dbgCompParm) {
156 show_matrix(U);
157 }
158 Matrix_Free(U);
159} /* linearInter */
160
161/**
162 * Given a system of equalities, looks if it has an integer solution in the
163 * combined space, and if yes, returns one solution.
164 * <p>pre-condition: the equalities are full-row rank (without the constant
165 * part)</p>
166 * @param Eqs the system of equations (as constraints)
167 * @param I a feasible integer solution if it exists, else NULL. Allocated if
168 * initially set to NULL, else reused.
169 */
171 Matrix *Hm, *H = NULL, *U, *Q, *M = NULL, *C = NULL, *Hi;
172 Matrix *Ip;
173 int i;
174 Value mod;
175 unsigned int rk;
176 if (Eqs == NULL) {
177 if ((*I) != NULL)
178 Matrix_Free(*I);
179 I = NULL;
180 return;
181 }
182 /* we use: AI = C = (Ha 0).Q.I = (Ha 0)(I' 0)^T */
183 /* with I = Qinv.I' = U.I'*/
184 /* 1- compute I' = Hainv.(-C) */
185 /* HYP: the equalities are full-row rank */
186 rk = Eqs->NbRows;
187 Matrix_subMatrix(Eqs, 0, 1, rk, Eqs->NbColumns - 1, &M);
188 left_hermite(M, &Hm, &Q, &U);
189 Matrix_Free(M);
190 Matrix_subMatrix(Hm, 0, 0, rk, rk, &H);
191 if (dbgCompParmMore) {
192 show_matrix(Hm);
193 show_matrix(H);
194 show_matrix(U);
195 }
196 Matrix_Free(Q);
197 Matrix_Free(Hm);
198 Matrix_subMatrix(Eqs, 0, Eqs->NbColumns - 1, rk, Eqs->NbColumns, &C);
199 Matrix_oppose(C);
200 Hi = Matrix_Alloc(rk, rk + 1);
201 MatInverse(H, Hi);
202 if (dbgCompParmMore) {
203 show_matrix(C);
204 show_matrix(Hi);
205 }
206 /* put the numerator of Hinv back into H */
207 Matrix_subMatrix(Hi, 0, 0, rk, rk, &H);
208 Ip = Matrix_Alloc(Eqs->NbColumns - 2, 1);
209 /* fool Matrix_Product on the size of Ip */
210 Ip->NbRows = rk;
211 Matrix_Product(H, C, Ip);
212 Ip->NbRows = Eqs->NbColumns - 2;
213 Matrix_Free(H);
214 Matrix_Free(C);
215 value_init(mod);
216 for (i = 0; i < rk; i++) {
217 /* if Hinv.C is not integer, return NULL (no solution) */
218 value_pmodulus(mod, Ip->p[i][0], Hi->p[i][rk]);
219 if (value_notzero_p(mod)) {
220 if ((*I) != NULL)
221 Matrix_Free(*I);
222 value_clear(mod);
223 Matrix_Free(U);
224 Matrix_Free(Ip);
225 Matrix_Free(Hi);
226 I = NULL;
227 return;
228 } else {
229 value_pdivision(Ip->p[i][0], Ip->p[i][0], Hi->p[i][rk]);
230 }
231 }
232 /* fill the rest of I' with zeros */
233 for (i = rk; i < Eqs->NbColumns - 2; i++) {
234 value_set_si(Ip->p[i][0], 0);
235 }
236 value_clear(mod);
237 Matrix_Free(Hi);
238 /* 2 - Compute the particular solution I = U.(I' 0) */
239 ensureMatrix((*I), Eqs->NbColumns - 2, 1);
240 Matrix_Product(U, Ip, (*I));
241 Matrix_Free(U);
242 Matrix_Free(Ip);
243 if (dbgCompParm) {
244 show_matrix(*I);
245 }
246}
247
248/**
249 * Computes the validity lattice of a set of equalities. I.e., the lattice
250 * induced on the last <tt>b</tt> variables by the equalities involving the
251 * first <tt>a</tt> integer existential variables. The submatrix of Eqs that
252 * concerns only the existential variables (so the first a columns) is assumed
253 * to be full-row rank.
254 * @param Eqs the equalities
255 * @param a the number of existential integer variables, placed as first
256 * variables
257 * @param vl the (returned) validity lattice, in homogeneous form. It is
258 * allocated if initially set to null, or reused if already allocated.
259 */
260void Equalities_validityLattice(Matrix *Eqs, int a, Matrix **vl) {
261 unsigned int b = Eqs->NbColumns - 2 - a;
262 unsigned int r = Eqs->NbRows;
263 Matrix *A = NULL, *B = NULL, *I = NULL, *Lb = NULL, *sol = NULL;
264 Matrix *H, *U, *Q;
265 unsigned int i;
266
267 if (dbgCompParm) {
268 printf("Computing validity lattice induced by the %d first variables of:",
269 a);
270 show_matrix(Eqs);
271 }
272 if (b == 0) {
273 ensureMatrix((*vl), 1, 1);
274 value_set_si((*vl)->p[0][0], 1);
275 return;
276 }
277
278 /* 1- check that there is an integer solution to the equalities */
279 /* OPT: could change integerSolution's profile to allocate or not*/
281 /* if there is no integer solution, there is no validity lattice */
282 if (sol == NULL) {
283 if ((*vl) != NULL)
284 Matrix_Free(*vl);
285 return;
286 }
287 Matrix_subMatrix(Eqs, 0, 1, r, 1 + a, &A);
288 Matrix_subMatrix(Eqs, 0, 1 + a, r, 1 + a + b, &B);
289 linearInter(A, B, &I, &Lb);
290 Matrix_Free(A);
291 Matrix_Free(B);
292 Matrix_Free(I);
293 if (dbgCompParm) {
294 show_matrix(Lb);
295 }
296
297 /* 2- The linear part of the validity lattice is the left HNF of Lb */
298 left_hermite(Lb, &H, &Q, &U);
299 Matrix_Free(Lb);
300 Matrix_Free(Q);
301 Matrix_Free(U);
302
303 /* 3- build the validity lattice */
304 ensureMatrix((*vl), b + 1, b + 1);
305 Matrix_copySubMatrix(H, 0, 0, b, b, (*vl), 0, 0);
306 Matrix_Free(H);
307 for (i = 0; i < b; i++) {
308 value_assign((*vl)->p[i][b], sol->p[0][a + i]);
309 }
310 Matrix_Free(sol);
311 Vector_Set((*vl)->p[b], 0, b);
312 value_set_si((*vl)->p[b][b], 1);
313
314} /* validityLattice */
315
316/**
317 * Eliminate the columns corresponding to a list of eliminated parameters.
318 * @param M the constraints matrix whose columns are to be removed
319 * @param nbVars an offset to be added to the ranks of the variables to be
320 * removed
321 * @param elimParms the list of ranks of the variables to be removed
322 * @param newM (output) the matrix without the removed columns
323 */
324void Constraints_removeElimCols(Matrix *M, unsigned int nbVars,
325 unsigned int *elimParms, Matrix **newM) {
326 unsigned int i, j, k;
327 if (elimParms[0] == 0) {
328 Matrix_clone(M, newM);
329 return;
330 }
331 if ((*newM) == NULL) {
332 (*newM) = Matrix_Alloc(M->NbRows, M->NbColumns - elimParms[0]);
333 } else {
334 assert((*newM)->NbColumns == M->NbColumns - elimParms[0]);
335 }
336 for (i = 0; i < M->NbRows; i++) {
337 value_assign((*newM)->p[i][0], M->p[i][0]); /* kind of cstr */
338 k = 0;
339 Vector_Copy(&(M->p[i][1]), &((*newM)->p[i][1]), nbVars);
340 for (j = 0; j < M->NbColumns - 2 - nbVars; j++) {
341 if (j != elimParms[k + 1]) {
342 value_assign((*newM)->p[i][j - k + nbVars + 1],
343 M->p[i][j + nbVars + 1]);
344 } else {
345 k++;
346 }
347 }
348 value_assign((*newM)->p[i][(*newM)->NbColumns - 1],
349 M->p[i][M->NbColumns - 1]); /* cst part */
350 }
351} /* Constraints_removeElimCols */
352
353/**
354 * Eliminates all the equalities in a set of constraints and returns the set of
355 * constraints defining a full-dimensional polyhedron, such that there is a
356 * bijection between integer points of the original polyhedron and these of the
357 * resulting (projected) polyhedron).
358 * If VL is set to NULL, this funciton allocates it. Else, it assumes that
359 * (*VL) points to a matrix of the right size.
360 * <p> The following things are done:
361 * <ol>
362 * <li> remove equalities involving only parameters, and remove as many
363 * parameters as there are such equalities. From that, the list of
364 * eliminated parameters <i>elimParms</i> is built.
365 * <li> remove equalities that involve variables. This requires a compression
366 * of the parameters and of the other variables that are not eliminated.
367 * The affine compresson is represented by matrix VL (for <i>validity
368 * lattice</i>) and is such that (N I 1)^T = VL.(N' I' 1), where N', I'
369 * are integer (they are the parameters and variables after compression).
370 *</ol>
371 *</p>
372 */
374 Matrix **Eqs, Matrix **ParmEqs,
375 unsigned int **elimVars,
376 unsigned int **elimParms, int maxRays) {
377 unsigned int i, j;
378 Matrix *A = NULL, *B = NULL;
379 Matrix *Ineqs = NULL;
380 unsigned int nbVars = (*M)->NbColumns - (*C)->NbColumns;
381 unsigned int nbParms;
382 int nbElimVars;
383 Matrix *fullDim = NULL;
384
385 /* variables for permutations */
386 unsigned int *permutation;
387 Matrix *permutedEqs = NULL, *permutedIneqs = NULL;
388
389 /* 1- Eliminate the equalities involving only parameters. */
390 (*ParmEqs) = Constraints_removeParmEqs(M, C, 0, elimParms);
391 /* if the polyehdron is empty, return now. */
392 if ((*M)->NbColumns == 0)
393 return;
394 /* eliminate the columns corresponding to the eliminated parameters */
395 if (elimParms[0] != 0) {
396 Constraints_removeElimCols(*M, nbVars, (*elimParms), &A);
397 Matrix_Free(*M);
398 (*M) = A;
399 Constraints_removeElimCols(*C, 0, (*elimParms), &B);
400 Matrix_Free(*C);
401 (*C) = B;
402 if (dbgCompParm) {
403 printf("After false parameter elimination: \n");
404 show_matrix(*M);
405 show_matrix(*C);
406 }
407 }
408 nbParms = (*C)->NbColumns - 2;
409
410 /* 2- Eliminate the equalities involving variables */
411 /* a- extract the (remaining) equalities from the poyhedron */
412 split_constraints((*M), Eqs, &Ineqs);
413 nbElimVars = (*Eqs)->NbRows;
414 /* if the polyhedron is already full-dimensional, return */
415 if ((*Eqs)->NbRows == 0) {
416 Matrix_identity(nbParms + 1, VL);
417 return;
418 }
419 /* b- choose variables to be eliminated */
420 permutation = find_a_permutation((*Eqs), nbParms);
421
422 if (dbgCompParm) {
423 printf("Permuting the vars/parms this way: [ ");
424 for (i = 0; i < (*Eqs)->NbColumns - 2; i++) {
425 printf("%d ", permutation[i]);
426 }
427 printf("]\n");
428 }
429
430 Constraints_permute((*Eqs), permutation, &permutedEqs);
431 Equalities_validityLattice(permutedEqs, (*Eqs)->NbRows, VL);
432
433 if (dbgCompParm) {
434 printf("Validity lattice: ");
435 show_matrix(*VL);
436 }
437 Constraints_compressLastVars(permutedEqs, (*VL));
438 Constraints_permute(Ineqs, permutation, &permutedIneqs);
439 if (dbgCompParmMore) {
440 show_matrix(permutedIneqs);
441 show_matrix(permutedEqs);
442 }
443 Matrix_Free(*Eqs);
444 Matrix_Free(Ineqs);
445 Constraints_compressLastVars(permutedIneqs, (*VL));
446 if (dbgCompParm) {
447 printf("After compression: ");
448 show_matrix(permutedIneqs);
449 }
450 /* c- eliminate the first variables */
451 assert(Constraints_eliminateFirstVars(permutedEqs, permutedIneqs));
452 if (dbgCompParmMore) {
453 printf("After elimination of the variables: ");
454 show_matrix(permutedIneqs);
455 }
456
457 /* d- get rid of the first (zero) columns,
458 which are now useless, and put the parameters back at the end */
459 fullDim = Matrix_Alloc(permutedIneqs->NbRows,
460 permutedIneqs->NbColumns - nbElimVars);
461 for (i = 0; i < permutedIneqs->NbRows; i++) {
462 value_set_si(fullDim->p[i][0], 1);
463 for (j = 0; j < nbParms; j++) {
464 value_assign(fullDim->p[i][j + fullDim->NbColumns - nbParms - 1],
465 permutedIneqs->p[i][j + nbElimVars + 1]);
466 }
467 for (j = 0; j < permutedIneqs->NbColumns - nbParms - 2 - nbElimVars; j++) {
468 value_assign(fullDim->p[i][j + 1],
469 permutedIneqs->p[i][nbElimVars + nbParms + j + 1]);
470 }
471 value_assign(fullDim->p[i][fullDim->NbColumns - 1],
472 permutedIneqs->p[i][permutedIneqs->NbColumns - 1]);
473 }
474 Matrix_Free(permutedIneqs);
475
476} /* Constraints_fullDimensionize */
477
478/**
479 * Given a matrix that defines a full-dimensional affine lattice, returns the
480 * affine sub-lattice spanned in the k first dimensions.
481 * Useful for instance when you only look for the parameters' validity lattice.
482 * @param lat the original full-dimensional lattice
483 * @param subLat the sublattice
484 */
485void Lattice_extractSubLattice(Matrix *lat, unsigned int k, Matrix **subLat) {
486 Matrix *H, *Q, *U, *linLat = NULL;
487 unsigned int i;
489 /* if the dimension is already good, just copy the initial lattice */
490 if (k == lat->NbRows - 1) {
491 if (*subLat == NULL) {
492 (*subLat) = Matrix_Copy(lat);
493 } else {
494 Matrix_copySubMatrix(lat, 0, 0, lat->NbRows, lat->NbColumns, (*subLat), 0,
495 0);
496 }
497 return;
498 }
499 assert(k < lat->NbRows - 1);
500 /* 1- Make the linear part of the lattice triangular to eliminate terms from
501 other dimensions */
502 Matrix_subMatrix(lat, 0, 0, lat->NbRows, lat->NbColumns - 1, &linLat);
503 /* OPT: any integer column-vector elimination is ok indeed. */
504 /* OPT: could test if the lattice is already in triangular form. */
505 left_hermite(linLat, &H, &Q, &U);
506 if (dbgCompParmMore) {
507 show_matrix(H);
508 }
509 Matrix_Free(Q);
510 Matrix_Free(U);
511 Matrix_Free(linLat);
512 /* if not allocated yet, allocate it */
513 if (*subLat == NULL) {
514 (*subLat) = Matrix_Alloc(k + 1, k + 1);
515 }
516 Matrix_copySubMatrix(H, 0, 0, k, k, (*subLat), 0, 0);
517 Matrix_Free(H);
518 Matrix_copySubMatrix(lat, 0, lat->NbColumns - 1, k, 1, (*subLat), 0, k);
519 for (i = 0; i < k; i++) {
520 value_set_si((*subLat)->p[k][i], 0);
521 }
522 value_set_si((*subLat)->p[k][k], 1);
524} /* Lattice_extractSubLattice */
525
526/**
527 * Computes the overall period of the variables I for (MI) mod |d|, where M is
528 * a matrix and |d| a vector. Produce a diagonal matrix S = (s_k) where s_k is
529 * the overall period of i_k
530 * @param M the set of affine functions of I (row-vectors)
531 * @param d the column-vector representing the modulos
532 */
534 Matrix *S;
535 unsigned int i, j;
536 Value tmp;
537 Value *periods = (Value *)malloc(sizeof(Value) * M->NbColumns);
538 value_init(tmp);
539 for (i = 0; i < M->NbColumns; i++) {
540 value_init(periods[i]);
541 value_set_si(periods[i], 1);
542 }
543 for (i = 0; i < M->NbRows; i++) {
544 for (j = 0; j < M->NbColumns; j++) {
545 value_gcd(tmp, d->p[i][0], M->p[i][j]);
546 value_divexact(tmp, d->p[i][0], tmp);
547 value_lcm(periods[j], periods[j], tmp);
548 }
549 }
550 value_clear(tmp);
551
552 /* 2- build S */
553 S = Matrix_Alloc(M->NbColumns, M->NbColumns);
554 for (i = 0; i < M->NbColumns; i++)
555 for (j = 0; j < M->NbColumns; j++)
556 if (i == j)
557 value_assign(S->p[i][j], periods[j]);
558 else
559 value_set_si(S->p[i][j], 0);
560
561 /* 3- clean up */
562 for (i = 0; i < M->NbColumns; i++)
563 value_clear(periods[i]);
564 free(periods);
565 return S;
566} /* affine_periods */
567
568/**
569 * Given an integer matrix B with m rows and integer m-vectors C and d,
570 * computes the basis of the integer solutions to (BN+C) mod d = 0 (1).
571 * This is an affine lattice (G): (N 1)^T= G(N' 1)^T, forall N' in Z^b.
572 * If there is no solution, returns NULL.
573 * @param B B, a (m x b) matrix
574 * @param C C, a (m x 1) integer matrix
575 * @param d d, a (1 x m) integer matrix
576 * @param imb the affine (b+1)x(b+1) basis of solutions, in the homogeneous
577 * form. Allocated if initially set to NULL, reused if not.
578 */
580 int b = B->NbColumns;
581 /* FIXME: treat the case d=0 as a regular equality B_kN+C_k = 0: */
582 /* OPT: could keep only equalities for which d>1 */
583 int nbEqs = B->NbRows;
584 unsigned int i;
585
586 /* 1- buid the problem DI+BN+C = 0 */
587 Matrix *eqs = Matrix_Alloc(nbEqs, nbEqs + b + 1);
588 for (i = 0; i < nbEqs; i++) {
589 value_assign(eqs->p[i][i], d->p[0][i]);
590 }
591 Matrix_copySubMatrix(B, 0, 0, nbEqs, b, eqs, 0, nbEqs);
592 Matrix_copySubMatrix(C, 0, 0, nbEqs, 1, eqs, 0, nbEqs + b);
593
594 /* 2- the solution is the validity lattice of the equalities */
595 Equalities_validityLattice(eqs, nbEqs, imb);
596 Matrix_Free(eqs);
597} /* Equalities_intModBasis */
598
599/** kept here for backwards compatiblity. Wrapper to Equalities_intModBasis() */
601 Matrix *imb = NULL;
602 Equalities_intModBasis(B, C, d, &imb);
603 return imb;
604} /* int_mod_basis */
605
606/**
607 * Given a parameterized constraints matrix with m equalities, computes the
608 * compression matrix G such that there is an integer solution in the variables
609 * space for each value of N', with N = G N' (N are the "nb_parms" parameters)
610 * @param E a matrix of parametric equalities @param nb_parms the number of
611 * parameters
612 * <b>Note: </b>this function is mostly here for backwards
613 * compatibility. Prefer the use of <tt>Equalities_validityLattice</tt>.
614 */
615Matrix *compress_parms(Matrix *E, int nbParms) {
616 Matrix *vl = NULL;
617 Equalities_validityLattice(E, E->NbColumns - 2 - nbParms, &vl);
618 return vl;
619} /* compress_parms */
620
621/** Removes the equalities that involve only parameters, by eliminating some
622 * parameters in the polyhedron's constraints and in the context.<p>
623 * <b>Updates M and Ctxt.</b>
624 * @param M1 the polyhedron's constraints
625 * @param Ctxt1 the constraints of the polyhedron's context
626 * @param renderSpace tells if the returned equalities must be expressed in the
627 * parameters space (renderSpace=0) or in the combined var/parms space
628 * (renderSpace = 1)
629 * @param elimParms the list of parameters that have been removed: an array
630 * whose 1st element is the number of elements in the list. (returned)
631 * @return the system of equalities that involve only parameters.
632 */
634 int renderSpace, unsigned int **elimParms) {
635 int i, j, k, nbEqsParms = 0;
636 int nbEqsM, nbEqsCtxt, allZeros, nbTautoM = 0, nbTautoCtxt = 0;
637 Matrix *M = (*M1);
638 Matrix *Ctxt = (*Ctxt1);
639 int nbVars = M->NbColumns - Ctxt->NbColumns;
640 Matrix *Eqs;
641 Matrix *EqsMTmp;
642
643 /* 1- build the equality matrix(ces) */
644 nbEqsM = 0;
645 for (i = 0; i < M->NbRows; i++) {
646 k = First_Non_Zero(M->p[i], M->NbColumns);
647 /* if it is a tautology, count it as such */
648 if (k == -1) {
649 nbTautoM++;
650 } else {
651 /* if it only involves parameters, count it */
652 if (k >= nbVars + 1)
653 nbEqsM++;
654 }
655 }
656
657 nbEqsCtxt = 0;
658 for (i = 0; i < Ctxt->NbRows; i++) {
659 if (value_zero_p(Ctxt->p[i][0])) {
660 if (First_Non_Zero(Ctxt->p[i], Ctxt->NbColumns) == -1) {
661 nbTautoCtxt++;
662 } else {
663 nbEqsCtxt++;
664 }
665 }
666 }
667 nbEqsParms = nbEqsM + nbEqsCtxt;
668
669 /* nothing to do in this case */
670 if (nbEqsParms + nbTautoM + nbTautoCtxt == 0) {
671 (*elimParms) = (unsigned int *)malloc(sizeof(int));
672 (*elimParms)[0] = 0;
673 if (renderSpace == 0) {
674 return Matrix_Alloc(0, Ctxt->NbColumns);
675 } else {
676 return Matrix_Alloc(0, M->NbColumns);
677 }
678 }
679
680 Eqs = Matrix_Alloc(nbEqsParms, Ctxt->NbColumns);
681 EqsMTmp = Matrix_Alloc(nbEqsParms, M->NbColumns);
682
683 /* copy equalities from the context */
684 k = 0;
685 for (i = 0; i < Ctxt->NbRows; i++) {
686 if (value_zero_p(Ctxt->p[i][0]) &&
687 First_Non_Zero(Ctxt->p[i], Ctxt->NbColumns) != -1) {
688 Vector_Copy(Ctxt->p[i], Eqs->p[k], Ctxt->NbColumns);
689 Vector_Copy(Ctxt->p[i] + 1, EqsMTmp->p[k] + nbVars + 1,
690 Ctxt->NbColumns - 1);
691 k++;
692 }
693 }
694 for (i = 0; i < M->NbRows; i++) {
695 j = First_Non_Zero(M->p[i], M->NbColumns);
696 /* copy equalities that involve only parameters from M */
697 if (j >= nbVars + 1) {
698 Vector_Copy(M->p[i] + nbVars + 1, Eqs->p[k] + 1, Ctxt->NbColumns - 1);
699 Vector_Copy(M->p[i] + nbVars + 1, EqsMTmp->p[k] + nbVars + 1,
700 Ctxt->NbColumns - 1);
701 /* mark these equalities for removal */
702 value_set_si(M->p[i][0], 2);
703 k++;
704 }
705 /* mark the all-zero equalities for removal */
706 if (j == -1) {
707 value_set_si(M->p[i][0], 2);
708 }
709 }
710
711 /* 2- eliminate parameters until all equalities are used or until we find a
712 contradiction (overconstrained system) */
713 (*elimParms) = (unsigned int *)malloc((Eqs->NbRows + 1) * sizeof(int));
714 (*elimParms)[0] = 0;
715 allZeros = 0;
716 for (i = 0; i < Eqs->NbRows; i++) {
717 /* find a variable that can be eliminated */
718 k = First_Non_Zero(Eqs->p[i], Eqs->NbColumns);
719 if (k != -1) { /* nothing special to do for tautologies */
720
721 /* if there is a contradiction, return empty matrices */
722 if (k == Eqs->NbColumns - 1) {
723 printf("Contradiction in %dth row of Eqs: ", k);
724 show_matrix(Eqs);
725 Matrix_Free(Eqs);
726 Matrix_Free(EqsMTmp);
727 (*M1) = Matrix_Alloc(0, M->NbColumns);
728 Matrix_Free(M);
729 (*Ctxt1) = Matrix_Alloc(0, Ctxt->NbColumns);
730 Matrix_Free(Ctxt);
731 free(*elimParms);
732 (*elimParms) = (unsigned int *)malloc(sizeof(int));
733 (*elimParms)[0] = 0;
734 if (renderSpace == 1) {
735 return Matrix_Alloc(0, (*M1)->NbColumns);
736 } else {
737 return Matrix_Alloc(0, (*Ctxt1)->NbColumns);
738 }
739 }
740 /* if we have something we can eliminate, do it in 3 places:
741 Eqs, Ctxt, and M */
742 else {
743 k--; /* k is the rank of the variable, now */
744 (*elimParms)[0]++;
745 (*elimParms)[(*elimParms[0])] = k;
746 for (j = 0; j < Eqs->NbRows; j++) {
747 if (i != j) {
748 eliminate_var_with_constr(Eqs, i, Eqs, j, k);
749 eliminate_var_with_constr(EqsMTmp, i, EqsMTmp, j, k + nbVars);
750 }
751 }
752 for (j = 0; j < Ctxt->NbRows; j++) {
753 if (value_notzero_p(Ctxt->p[i][0])) {
754 eliminate_var_with_constr(Eqs, i, Ctxt, j, k);
755 }
756 }
757 for (j = 0; j < M->NbRows; j++) {
758 if (value_cmp_si(M->p[i][0], 2)) {
759 eliminate_var_with_constr(EqsMTmp, i, M, j, k + nbVars);
760 }
761 }
762 }
763 }
764 /* if (k==-1): count the tautologies in Eqs to remove them later */
765 else {
766 allZeros++;
767 }
768 }
769
770 /* elimParms may have been overallocated. Now we know how many parms have
771 been eliminated so we can reallocate the right amount of memory. */
772 if (!realloc((*elimParms), ((*elimParms)[0] + 1) * sizeof(int))) {
773 fprintf(stderr, "Constraints_Remove_parm_eqs > cannot realloc()");
774 }
775
776 Matrix_Free(EqsMTmp);
777
778 /* 3- remove the "bad" equalities from the input matrices
779 and copy the equalities involving only parameters */
780 EqsMTmp = Matrix_Alloc(M->NbRows - nbEqsM - nbTautoM, M->NbColumns);
781 k = 0;
782 for (i = 0; i < M->NbRows; i++) {
783 if (value_cmp_si(M->p[i][0], 2)) {
784 Vector_Copy(M->p[i], EqsMTmp->p[k], M->NbColumns);
785 k++;
786 }
787 }
788 Matrix_Free(M);
789 (*M1) = EqsMTmp;
790
791 EqsMTmp =
792 Matrix_Alloc(Ctxt->NbRows - nbEqsCtxt - nbTautoCtxt, Ctxt->NbColumns);
793 k = 0;
794 for (i = 0; i < Ctxt->NbRows; i++) {
795 if (value_notzero_p(Ctxt->p[i][0])) {
796 Vector_Copy(Ctxt->p[i], EqsMTmp->p[k], Ctxt->NbColumns);
797 k++;
798 }
799 }
800 Matrix_Free(Ctxt);
801 (*Ctxt1) = EqsMTmp;
802
803 if (renderSpace == 0) { /* renderSpace=0: equalities in the parameter space */
804 EqsMTmp = Matrix_Alloc(Eqs->NbRows - allZeros, Eqs->NbColumns);
805 k = 0;
806 for (i = 0; i < Eqs->NbRows; i++) {
807 if (First_Non_Zero(Eqs->p[i], Eqs->NbColumns) != -1) {
808 Vector_Copy(Eqs->p[i], EqsMTmp->p[k], Eqs->NbColumns);
809 k++;
810 }
811 }
812 } else { /* renderSpace=1: equalities rendered in the combined space */
813 EqsMTmp = Matrix_Alloc(Eqs->NbRows - allZeros, (*M1)->NbColumns);
814 k = 0;
815 for (i = 0; i < Eqs->NbRows; i++) {
816 if (First_Non_Zero(Eqs->p[i], Eqs->NbColumns) != -1) {
817 Vector_Copy(Eqs->p[i], &(EqsMTmp->p[k][nbVars]), Eqs->NbColumns);
818 k++;
819 }
820 }
821 }
822 Matrix_Free(Eqs);
823 Eqs = EqsMTmp;
824
825 return Eqs;
826} /* Constraints_Remove_parm_eqs */
827
828/** Removes equalities involving only parameters, but starting from a
829 * Polyhedron and its context.
830 * @param P the polyhedron
831 * @param C P's context
832 * @param renderSpace: 0 for the parameter space, =1 for the combined space.
833 * @maxRays Polylib's usual <i>workspace</i>.
834 */
836 int renderSpace,
837 unsigned int **elimParms, int maxRays) {
838 Matrix *Eqs;
839 Polyhedron *Peqs;
841 Matrix *Ct = Polyhedron2Constraints((*C));
842
843 /* if the Minkowski representation is not computed yet, do not compute it in
844 Constraints2Polyhedron */
845 if (F_ISSET((*P), POL_VALID | POL_INEQUALITIES) &&
848 }
849
850 Eqs = Constraints_Remove_parm_eqs(&M, &Ct, renderSpace, elimParms);
851 Peqs = Constraints2Polyhedron(Eqs, maxRays);
852 Matrix_Free(Eqs);
853
854 /* particular case: no equality involving only parms is found */
855 if (Eqs->NbRows == 0) {
856 Matrix_Free(M);
857 Matrix_Free(Ct);
858 return Peqs;
859 }
860 Polyhedron_Free(*P);
861 Polyhedron_Free(*C);
864 Matrix_Free(M);
865 Matrix_Free(Ct);
866 return Peqs;
867} /* Polyhedron_Remove_parm_eqs */
868
869/**
870 * Given a matrix with m parameterized equations, compress the nb_parms
871 * parameters and n-m variables so that m variables are integer, and transform
872 * the variable space into a n-m space by eliminating the m variables (using
873 * the equalities) the variables to be eliminated are chosen automatically by
874 * the function.
875 * <b>Deprecated.</b> Try to use Constraints_fullDimensionize instead.
876 * @param M the constraints
877 * @param the number of parameters
878 * @param validityLattice the the integer lattice underlying the integer
879 * solutions.
880 */
881Matrix *full_dimensionize(Matrix const *M, int nbParms,
882 Matrix **validityLattice) {
883 Matrix *Eqs, *Ineqs;
884 Matrix *permutedEqs, *permutedIneqs;
885 Matrix *Full_Dim;
886 Matrix *WVL; /* The Whole Validity Lattice (vars+parms) */
887 unsigned int i, j;
888 int nbElimVars;
889 unsigned int *permutation, *permutationInv;
890 /* 0- Split the equalities and inequalities from each other */
891 split_constraints(M, &Eqs, &Ineqs);
892
893 /* 1- if the polyhedron is already full-dimensional, return it */
894 if (Eqs->NbRows == 0) {
895 Matrix_Free(Eqs);
896 (*validityLattice) = Identity_Matrix(nbParms + 1);
897 return Ineqs;
898 }
899 nbElimVars = Eqs->NbRows;
900
901 /* 2- put the vars to be eliminated at the first positions,
902 and compress the other vars/parms
903 -> [ variables to eliminate / parameters / variables to keep ] */
904 permutation = find_a_permutation(Eqs, nbParms);
905 if (dbgCompParm) {
906 printf("Permuting the vars/parms this way: [ ");
907 for (i = 0; i < Eqs->NbColumns; i++) {
908 printf("%d ", permutation[i]);
909 }
910 printf("]\n");
911 }
912 permutedEqs = mpolyhedron_permute(Eqs, permutation);
913 WVL = compress_parms(permutedEqs, Eqs->NbColumns - 2 - Eqs->NbRows);
914 /* Check for empty WVL (no solution) */
915 if (!WVL) {
916 fprintf(stderr, "full_dimensionize > parameters compression failed.\n");
917 Matrix_Free(Eqs);
918 Matrix_Free(Ineqs);
919 (*validityLattice) = Identity_Matrix(nbParms + 1);
920 return NULL;
921 }
922 if (dbgCompParm) {
923 printf("Whole validity lattice: ");
924 show_matrix(WVL);
925 }
926 mpolyhedron_compress_last_vars(permutedEqs, WVL);
927 permutedIneqs = mpolyhedron_permute(Ineqs, permutation);
928 if (dbgCompParm) {
929 show_matrix(permutedEqs);
930 }
931 Matrix_Free(Eqs);
932 Matrix_Free(Ineqs);
933 mpolyhedron_compress_last_vars(permutedIneqs, WVL);
934 if (dbgCompParm) {
935 printf("After compression: ");
936 show_matrix(permutedIneqs);
937 }
938 /* 3- eliminate the first variables */
939 if (!mpolyhedron_eliminate_first_variables(permutedEqs, permutedIneqs)) {
940 fprintf(stderr, "full_dimensionize > variable elimination failed.\n");
941 Matrix_Free(permutedIneqs);
942 (*validityLattice) = Identity_Matrix(nbParms + 1);
943 return NULL;
944 }
945 if (dbgCompParm) {
946 printf("After elimination of the variables: ");
947 show_matrix(permutedIneqs);
948 }
949
950 /* 4- get rid of the first (zero) columns,
951 which are now useless, and put the parameters back at the end */
952 Full_Dim = Matrix_Alloc(permutedIneqs->NbRows,
953 permutedIneqs->NbColumns - nbElimVars);
954 for (i = 0; i < permutedIneqs->NbRows; i++) {
955 value_set_si(Full_Dim->p[i][0], 1);
956 for (j = 0; j < nbParms; j++)
957 value_assign(Full_Dim->p[i][j + Full_Dim->NbColumns - nbParms - 1],
958 permutedIneqs->p[i][j + nbElimVars + 1]);
959 for (j = 0; j < permutedIneqs->NbColumns - nbParms - 2 - nbElimVars; j++)
960 value_assign(Full_Dim->p[i][j + 1],
961 permutedIneqs->p[i][nbElimVars + nbParms + j + 1]);
962 value_assign(Full_Dim->p[i][Full_Dim->NbColumns - 1],
963 permutedIneqs->p[i][permutedIneqs->NbColumns - 1]);
964 }
965 Matrix_Free(permutedIneqs);
966
967 /* 5- Keep only the the validity lattice restricted to the parameters */
968 *validityLattice = Matrix_Alloc(nbParms + 1, nbParms + 1);
969 for (i = 0; i < nbParms; i++) {
970 for (j = 0; j < nbParms; j++)
971 value_assign((*validityLattice)->p[i][j], WVL->p[i][j]);
972 value_assign((*validityLattice)->p[i][nbParms],
973 WVL->p[i][WVL->NbColumns - 1]);
974 }
975 for (j = 0; j < nbParms; j++)
976 value_set_si((*validityLattice)->p[nbParms][j], 0);
977 value_assign((*validityLattice)->p[nbParms][nbParms],
978 WVL->p[WVL->NbColumns - 1][WVL->NbColumns - 1]);
979
980 /* 6- Clean up */
981 Matrix_Free(WVL);
982 return Full_Dim;
983} /* full_dimensionize */
984
985#undef dbgCompParm
986#undef dbgCompParmMore
Matrix * Matrix_Copy(Matrix const *Src)
Definition: Matop.c:98
#define value_pdivision(ref, val1, val2)
Definition: arithmetique.h:550
#define value_notzero_p(val)
Definition: arithmetique.h:576
#define value_divexact(ref, val1, val2)
Definition: arithmetique.h:548
#define value_gcd(ref, val1, val2)
Definition: arithmetique.h:556
#define value_pmodulus(ref, val1, val2)
Definition: arithmetique.h:551
#define value_lcm(ref, val1, val2)
Definition: arithmetique.h:557
#define value_zero_p(val)
Definition: arithmetique.h:575
#define value_assign(v1, v2)
Definition: arithmetique.h:482
#define value_set_si(val, i)
Definition: arithmetique.h:483
#define value_clear(val)
Definition: arithmetique.h:485
#define value_init(val)
Definition: arithmetique.h:481
#define value_cmp_si(val, n)
Definition: arithmetique.h:581
#define assert(ex)
Definition: assert.h:16
void Equalities_integerSolution(Matrix *Eqs, Matrix **I)
Given a system of equalities, looks if it has an integer solution in the combined space,...
Matrix * compress_parms(Matrix *E, int nbParms)
Given a parameterized constraints matrix with m equalities, computes the compression matrix G such th...
void Equalities_intModBasis(Matrix *B, Matrix *C, Matrix *d, Matrix **imb)
Given an integer matrix B with m rows and integer m-vectors C and d, computes the basis of the intege...
void Constraints_fullDimensionize(Matrix **M, Matrix **C, Matrix **VL, Matrix **Eqs, Matrix **ParmEqs, unsigned int **elimVars, unsigned int **elimParms, int maxRays)
Eliminates all the equalities in a set of constraints and returns the set of constraints defining a f...
#define dbgCompParm
debug flags (2 levels)
void Constraints_removeElimCols(Matrix *M, unsigned int nbVars, unsigned int *elimParms, Matrix **newM)
Eliminate the columns corresponding to a list of eliminated parameters.
Matrix * Constraints_Remove_parm_eqs(Matrix **M1, Matrix **Ctxt1, int renderSpace, unsigned int **elimParms)
Removes the equalities that involve only parameters, by eliminating some parameters in the polyhedron...
#define dbgEnd(a)
#define dbgStart(a)
Matrix * affine_periods(Matrix *M, Matrix *d)
Computes the overall period of the variables I for (MI) mod |d|, where M is a matrix and |d| a vector...
static void linearInter(Matrix *A, Matrix *B, Matrix **I, Matrix **Lb)
Computes the intersection of two linear lattices, whose base vectors are respectively represented in ...
Matrix * int_mod_basis(Matrix *B, Matrix *C, Matrix *d)
kept here for backwards compatiblity.
#define dbgCompParmMore
void Lattice_extractSubLattice(Matrix *lat, unsigned int k, Matrix **subLat)
Given a matrix that defines a full-dimensional affine lattice, returns the affine sub-lattice spanned...
Matrix * int_ker(Matrix *M)
Given a full-row-rank nxm matrix M made of m row-vectors), computes the basis K (made of n-m column-v...
void Equalities_validityLattice(Matrix *Eqs, int a, Matrix **vl)
Computes the validity lattice of a set of equalities.
Polyhedron * Polyhedron_Remove_parm_eqs(Polyhedron **P, Polyhedron **C, int renderSpace, unsigned int **elimParms, int maxRays)
Removes equalities involving only parameters, but starting from a Polyhedron and its context.
Matrix * full_dimensionize(Matrix const *M, int nbParms, Matrix **validityLattice)
Given a matrix with m parameterized equations, compress the nb_parms parameters and n-m variables so ...
#define Constraints_removeParmEqs(a, b, c, d)
void Matrix_Product(Matrix *Mat1, Matrix *Mat2, Matrix *Mat3)
Definition: matrix.c:872
Matrix * Matrix_Alloc(unsigned NbRows, unsigned NbColumns)
Definition: matrix.c:24
void right_hermite(Matrix *A, Matrix **Hp, Matrix **Up, Matrix **Qp)
Definition: matrix.c:449
void left_hermite(Matrix *A, Matrix **Hp, Matrix **Qp, Matrix **Up)
Definition: matrix.c:517
int MatInverse(Matrix *Mat, Matrix *MatInv)
Definition: matrix.c:603
void Matrix_Free(Matrix *Mat)
Definition: matrix.c:71
void split_constraints(Matrix const *M, Matrix **Eqs, Matrix **Ineqs)
splits a matrix of constraints M into a matrix of equalities Eqs and a matrix of inequalities Ineqs a...
Definition: matrix_addon.c:31
void Matrix_subMatrix(Matrix *M, unsigned int sr, unsigned int sc, unsigned int er, unsigned int ec, Matrix **sub)
returns a contiguous submatrix of a matrix.
Definition: matrix_addon.c:325
void mpolyhedron_compress_last_vars(Matrix *M, Matrix *compression)
compress the last vars/pars of the polyhedron M expressed as a polylib matrix
Definition: matrix_addon.c:261
Matrix * Identity_Matrix(unsigned int dim)
Definition: matrix_addon.c:58
unsigned int mpolyhedron_eliminate_first_variables(Matrix *Eqs, Matrix *Ineqs)
use a set of m equalities Eqs to eliminate m variables in the polyhedron Ineqs represented as a matri...
Definition: matrix_addon.c:288
void eliminate_var_with_constr(Matrix *Eliminator, unsigned int eliminator_row, Matrix *Victim, unsigned int victim_row, unsigned int var_to_elim)
use an eliminator row to eliminate a variable in a victim row (without changing the sign of the victi...
Definition: matrix_addon.c:213
void Matrix_clone(Matrix *M, Matrix **Cl)
Cloning function.
Definition: matrix_addon.c:345
void Matrix_copySubMatrix(Matrix *M1, unsigned int sr1, unsigned int sc1, unsigned int nbR, unsigned int nbC, Matrix *M2, unsigned int sr2, unsigned int sc2)
Copies a contiguous submatrix of M1 into M2, at the indicated position.
Definition: matrix_addon.c:361
void Matrix_oppose(Matrix *M)
transforms a matrix M into -M
Definition: matrix_addon.c:373
void Matrix_identity(unsigned int dim, Matrix **I)
returns the dim-dimensional identity matrix.
Definition: matrix_addon.c:77
#define show_matrix(M)
Polylib matrix addons Mainly, deals with polyhedra represented in implicit form (set of constraints).
Definition: matrix_addon.h:15
#define ensureMatrix(M, r, c)
Allocates a matrix if it is null, or else asserts that it has at least a certain size.
Definition: matrix_addon.h:28
#define Constraints_compressLastVars(a, b)
Definition: matrix_addon.h:83
#define Constraints_eliminateFirstVars(a, b)
Definition: matrix_addon.h:90
void Constraints_permute(Matrix *C, unsigned int *perm, Matrix **Cp)
permutes the variables of the constraints of a polyhedron
Matrix * mpolyhedron_permute(Matrix *polyh, unsigned int *permutation)
permutes the variables of the constraints of a polyhedron
unsigned int * find_a_permutation(Matrix *Eqs, unsigned int nb_parms)
finds a valid permutation : for a set of m equations, find m variables that will be put at the beginn...
void Polyhedron_Free(Polyhedron *Pol)
Definition: polyhedron.c:1614
Matrix * Polyhedron2Constraints(Polyhedron *Pol)
Definition: polyhedron.c:2060
Polyhedron * Constraints2Polyhedron(Matrix *Constraints, unsigned NbMaxRays)
Given a matrix of constraints ('Constraints'), construct and return a polyhedron.
Definition: polyhedron.c:1905
Definition: types.h:75
unsigned NbRows
Definition: types.h:76
Value ** p
Definition: types.h:77
unsigned NbColumns
Definition: types.h:76
#define maxRays
#define POL_INEQUALITIES
Definition: types.h:100
#define F_ISSET(p, f)
Definition: types.h:91
#define POL_VALID
Definition: types.h:107
#define POL_NO_DUAL
Definition: types.h:66
#define FL_INIT(l, f)
Definition: types.h:83
void Vector_Set(Value *p, int n, unsigned length)
Definition: vector.c:227
int Vector_IsZero(Value *v, unsigned length)
Definition: vector.c:724
void Vector_Copy(Value *p1, Value *p2, unsigned length)
Definition: vector.c:256
int First_Non_Zero(Value *p, unsigned length)
Definition: vector.c:120