rename barvinok_maximize to barvinok_bound
[barvinok.git] / util.c
blobd4dbe5093d537b58642c3db7b797ca9785558528
1 #include <stdlib.h>
2 #include <assert.h>
3 #include <barvinok/util.h>
4 #include <barvinok/options.h>
5 #include <polylib/ranking.h>
6 #include "config.h"
8 #define ALLOC(type) (type*)malloc(sizeof(type))
9 #define ALLOCN(type,n) (type*)malloc((n) * sizeof(type))
11 #ifdef __GNUC__
12 #define NALLOC(p,n) p = (typeof(p))malloc((n) * sizeof(*p))
13 #else
14 #define NALLOC(p,n) p = (void *)malloc((n) * sizeof(*p))
15 #endif
17 void manual_count(Polyhedron *P, Value* result)
19 Polyhedron *U = Universe_Polyhedron(0);
20 Enumeration *en = Polyhedron_Enumerate(P,U,1024,NULL);
21 Value *v = compute_poly(en,NULL);
22 value_assign(*result, *v);
23 value_clear(*v);
24 free(v);
25 Enumeration_Free(en);
26 Polyhedron_Free(U);
29 #include <barvinok/evalue.h>
30 #include <barvinok/util.h>
31 #include <barvinok/barvinok.h>
33 /* Return random value between 0 and max-1 inclusive
35 int random_int(int max) {
36 return (int) (((double)(max))*rand()/(RAND_MAX+1.0));
39 Polyhedron *Polyhedron_Read(unsigned MaxRays)
41 int vertices = 0;
42 unsigned NbRows, NbColumns;
43 Matrix *M;
44 Polyhedron *P;
45 char s[128];
47 while (fgets(s, sizeof(s), stdin)) {
48 if (*s == '#')
49 continue;
50 if (strncasecmp(s, "vertices", sizeof("vertices")-1) == 0)
51 vertices = 1;
52 if (sscanf(s, "%u %u", &NbRows, &NbColumns) == 2)
53 break;
55 if (feof(stdin))
56 return NULL;
57 M = Matrix_Alloc(NbRows,NbColumns);
58 Matrix_Read_Input(M);
59 if (vertices)
60 P = Rays2Polyhedron(M, MaxRays);
61 else
62 P = Constraints2Polyhedron(M, MaxRays);
63 Matrix_Free(M);
64 return P;
67 /* Inplace polarization
69 void Polyhedron_Polarize(Polyhedron *P)
71 unsigned NbRows = P->NbConstraints + P->NbRays;
72 int i;
73 Value **q;
75 q = (Value **)malloc(NbRows * sizeof(Value *));
76 assert(q);
77 for (i = 0; i < P->NbRays; ++i)
78 q[i] = P->Ray[i];
79 for (; i < NbRows; ++i)
80 q[i] = P->Constraint[i-P->NbRays];
81 P->NbConstraints = NbRows - P->NbConstraints;
82 P->NbRays = NbRows - P->NbRays;
83 free(P->Constraint);
84 P->Constraint = q;
85 P->Ray = q + P->NbConstraints;
89 * Rather general polar
90 * We can optimize it significantly if we assume that
91 * P includes zero
93 * Also, we calculate the polar as defined in Schrijver
94 * The opposite should probably work as well and would
95 * eliminate the need for multiplying by -1
97 Polyhedron* Polyhedron_Polar(Polyhedron *P, unsigned NbMaxRays)
99 int i;
100 Value mone;
101 unsigned dim = P->Dimension + 2;
102 Matrix *M = Matrix_Alloc(P->NbRays, dim);
104 assert(M);
105 value_init(mone);
106 value_set_si(mone, -1);
107 for (i = 0; i < P->NbRays; ++i) {
108 Vector_Scale(P->Ray[i], M->p[i], mone, dim);
109 value_multiply(M->p[i][0], M->p[i][0], mone);
110 value_multiply(M->p[i][dim-1], M->p[i][dim-1], mone);
112 P = Constraints2Polyhedron(M, NbMaxRays);
113 assert(P);
114 Matrix_Free(M);
115 value_clear(mone);
116 return P;
120 * Returns the supporting cone of P at the vertex with index v
122 Polyhedron* supporting_cone(Polyhedron *P, int v)
124 Matrix *M;
125 Value tmp;
126 int i, n, j;
127 unsigned char *supporting = (unsigned char *)malloc(P->NbConstraints);
128 unsigned dim = P->Dimension + 2;
130 assert(v >=0 && v < P->NbRays);
131 assert(value_pos_p(P->Ray[v][dim-1]));
132 assert(supporting);
134 value_init(tmp);
135 for (i = 0, n = 0; i < P->NbConstraints; ++i) {
136 Inner_Product(P->Constraint[i] + 1, P->Ray[v] + 1, dim - 1, &tmp);
137 if ((supporting[i] = value_zero_p(tmp)))
138 ++n;
140 assert(n >= dim - 2);
141 value_clear(tmp);
142 M = Matrix_Alloc(n, dim);
143 assert(M);
144 for (i = 0, j = 0; i < P->NbConstraints; ++i)
145 if (supporting[i]) {
146 value_set_si(M->p[j][dim-1], 0);
147 Vector_Copy(P->Constraint[i], M->p[j++], dim-1);
149 free(supporting);
150 P = Constraints2Polyhedron(M, P->NbRays+1);
151 assert(P);
152 Matrix_Free(M);
153 return P;
156 #define INT_BITS (sizeof(unsigned) * 8)
158 unsigned *supporting_constraints(Matrix *Constraints, Param_Vertices *v, int *n)
160 Value lcm, tmp, tmp2;
161 unsigned dim = Constraints->NbColumns;
162 unsigned nparam = v->Vertex->NbColumns - 2;
163 unsigned nvar = dim - nparam - 2;
164 int len = (Constraints->NbRows+INT_BITS-1)/INT_BITS;
165 unsigned *supporting = (unsigned *)calloc(len, sizeof(unsigned));
166 int i, j;
167 Vector *row;
168 int ix;
169 unsigned bx;
171 assert(supporting);
172 row = Vector_Alloc(nparam+1);
173 assert(row);
174 value_init(lcm);
175 value_init(tmp);
176 value_init(tmp2);
177 value_set_si(lcm, 1);
178 for (i = 0, *n = 0, ix = 0, bx = MSB; i < Constraints->NbRows; ++i) {
179 Vector_Set(row->p, 0, nparam+1);
180 for (j = 0 ; j < nvar; ++j) {
181 value_set_si(tmp, 1);
182 value_assign(tmp2, Constraints->p[i][j+1]);
183 if (value_ne(lcm, v->Vertex->p[j][nparam+1])) {
184 value_assign(tmp, lcm);
185 value_lcm(lcm, lcm, v->Vertex->p[j][nparam+1]);
186 value_division(tmp, lcm, tmp);
187 value_multiply(tmp2, tmp2, lcm);
188 value_division(tmp2, tmp2, v->Vertex->p[j][nparam+1]);
190 Vector_Combine(row->p, v->Vertex->p[j], row->p,
191 tmp, tmp2, nparam+1);
193 value_set_si(tmp, 1);
194 Vector_Combine(row->p, Constraints->p[i]+1+nvar, row->p, tmp, lcm, nparam+1);
195 for (j = 0; j < nparam+1; ++j)
196 if (value_notzero_p(row->p[j]))
197 break;
198 if (j == nparam + 1) {
199 supporting[ix] |= bx;
200 ++*n;
202 NEXT(ix, bx);
204 assert(*n >= nvar);
205 value_clear(tmp);
206 value_clear(tmp2);
207 value_clear(lcm);
208 Vector_Free(row);
210 return supporting;
213 Polyhedron* supporting_cone_p(Polyhedron *P, Param_Vertices *v)
215 Matrix *M;
216 unsigned dim = P->Dimension + 2;
217 unsigned nparam = v->Vertex->NbColumns - 2;
218 unsigned nvar = dim - nparam - 2;
219 int i, n, j;
220 int ix;
221 unsigned bx;
222 unsigned *supporting;
223 Matrix View;
225 Polyhedron_Matrix_View(P, &View, P->NbConstraints);
226 supporting = supporting_constraints(&View, v, &n);
227 M = Matrix_Alloc(n, nvar+2);
228 assert(M);
229 for (i = 0, j = 0, ix = 0, bx = MSB; i < P->NbConstraints; ++i) {
230 if (supporting[ix] & bx) {
231 value_set_si(M->p[j][nvar+1], 0);
232 Vector_Copy(P->Constraint[i], M->p[j++], nvar+1);
234 NEXT(ix, bx);
236 free(supporting);
237 P = Constraints2Polyhedron(M, P->NbRays+1);
238 assert(P);
239 Matrix_Free(M);
240 return P;
243 Polyhedron* triangulate_cone(Polyhedron *P, unsigned NbMaxCons)
245 struct barvinok_options *options = barvinok_options_new_with_defaults();
246 options->MaxRays = NbMaxCons;
247 P = triangulate_cone_with_options(P, options);
248 barvinok_options_free(options);
249 return P;
252 Polyhedron* triangulate_cone_with_options(Polyhedron *P,
253 struct barvinok_options *options)
255 const static int MAX_TRY=10;
256 int i, j, r, n, t;
257 Value tmp;
258 unsigned dim = P->Dimension;
259 Matrix *M = Matrix_Alloc(P->NbRays+1, dim+3);
260 Matrix *M2, *M3;
261 Polyhedron *L, *R, *T;
262 assert(P->NbEq == 0);
264 L = NULL;
265 R = NULL;
266 value_init(tmp);
268 Vector_Set(M->p[0]+1, 0, dim+1);
269 value_set_si(M->p[0][0], 1);
270 value_set_si(M->p[0][dim+2], 1);
271 Vector_Set(M->p[P->NbRays]+1, 0, dim+2);
272 value_set_si(M->p[P->NbRays][0], 1);
273 value_set_si(M->p[P->NbRays][dim+1], 1);
275 for (i = 0, r = 1; i < P->NbRays; ++i) {
276 if (value_notzero_p(P->Ray[i][dim+1]))
277 continue;
278 Vector_Copy(P->Ray[i], M->p[r], dim+1);
279 value_set_si(M->p[r][dim+2], 0);
280 ++r;
283 M2 = Matrix_Alloc(dim+1, dim+2);
285 t = 0;
286 if (options->try_Delaunay_triangulation) {
287 /* Delaunay triangulation */
288 for (r = 1; r < P->NbRays; ++r) {
289 Inner_Product(M->p[r]+1, M->p[r]+1, dim, &tmp);
290 value_assign(M->p[r][dim+1], tmp);
292 M3 = Matrix_Copy(M);
293 L = Rays2Polyhedron(M3, options->MaxRays);
294 Matrix_Free(M3);
295 ++t;
296 } else {
297 try_again:
298 /* Usually R should still be 0 */
299 Domain_Free(R);
300 Polyhedron_Free(L);
301 for (r = 1; r < P->NbRays; ++r) {
302 value_set_si(M->p[r][dim+1], random_int((t+1)*dim*P->NbRays)+1);
304 M3 = Matrix_Copy(M);
305 L = Rays2Polyhedron(M3, options->MaxRays);
306 Matrix_Free(M3);
307 ++t;
309 assert(t <= MAX_TRY);
311 R = NULL;
312 n = 0;
314 POL_ENSURE_FACETS(L);
315 for (i = 0; i < L->NbConstraints; ++i) {
316 /* Ignore perpendicular facets, i.e., facets with 0 z-coordinate */
317 if (value_negz_p(L->Constraint[i][dim+1]))
318 continue;
319 if (value_notzero_p(L->Constraint[i][dim+2]))
320 continue;
321 for (j = 1, r = 1; j < M->NbRows; ++j) {
322 Inner_Product(M->p[j]+1, L->Constraint[i]+1, dim+1, &tmp);
323 if (value_notzero_p(tmp))
324 continue;
325 if (r > dim)
326 goto try_again;
327 Vector_Copy(M->p[j]+1, M2->p[r]+1, dim);
328 value_set_si(M2->p[r][0], 1);
329 value_set_si(M2->p[r][dim+1], 0);
330 ++r;
332 assert(r == dim+1);
333 Vector_Set(M2->p[0]+1, 0, dim);
334 value_set_si(M2->p[0][0], 1);
335 value_set_si(M2->p[0][dim+1], 1);
336 T = Rays2Polyhedron(M2, P->NbConstraints+1);
337 T->next = R;
338 R = T;
339 ++n;
341 Matrix_Free(M2);
343 Polyhedron_Free(L);
344 value_clear(tmp);
345 Matrix_Free(M);
347 return R;
350 void check_triangulization(Polyhedron *P, Polyhedron *T)
352 Polyhedron *C, *D, *E, *F, *G, *U;
353 for (C = T; C; C = C->next) {
354 if (C == T)
355 U = C;
356 else
357 U = DomainConvex(DomainUnion(U, C, 100), 100);
358 for (D = C->next; D; D = D->next) {
359 F = C->next;
360 G = D->next;
361 C->next = NULL;
362 D->next = NULL;
363 E = DomainIntersection(C, D, 600);
364 assert(E->NbRays == 0 || E->NbEq >= 1);
365 Polyhedron_Free(E);
366 C->next = F;
367 D->next = G;
370 assert(PolyhedronIncludes(U, P));
371 assert(PolyhedronIncludes(P, U));
374 /* Computes x, y and g such that g = gcd(a,b) and a*x+b*y = g */
375 void Extended_Euclid(Value a, Value b, Value *x, Value *y, Value *g)
377 Value c, d, e, f, tmp;
379 value_init(c);
380 value_init(d);
381 value_init(e);
382 value_init(f);
383 value_init(tmp);
384 value_absolute(c, a);
385 value_absolute(d, b);
386 value_set_si(e, 1);
387 value_set_si(f, 0);
388 while(value_pos_p(d)) {
389 value_division(tmp, c, d);
390 value_multiply(tmp, tmp, f);
391 value_subtract(e, e, tmp);
392 value_division(tmp, c, d);
393 value_multiply(tmp, tmp, d);
394 value_subtract(c, c, tmp);
395 value_swap(c, d);
396 value_swap(e, f);
398 value_assign(*g, c);
399 if (value_zero_p(a))
400 value_set_si(*x, 0);
401 else if (value_pos_p(a))
402 value_assign(*x, e);
403 else value_oppose(*x, e);
404 if (value_zero_p(b))
405 value_set_si(*y, 0);
406 else {
407 value_multiply(tmp, a, *x);
408 value_subtract(tmp, c, tmp);
409 value_division(*y, tmp, b);
411 value_clear(c);
412 value_clear(d);
413 value_clear(e);
414 value_clear(f);
415 value_clear(tmp);
418 static int unimodular_complete_1(Matrix *m)
420 Value g, b, c, old, tmp;
421 unsigned i, j;
422 int ok;
424 value_init(b);
425 value_init(c);
426 value_init(g);
427 value_init(old);
428 value_init(tmp);
429 value_assign(g, m->p[0][0]);
430 for (i = 1; value_zero_p(g) && i < m->NbColumns; ++i) {
431 for (j = 0; j < m->NbColumns; ++j) {
432 if (j == i-1)
433 value_set_si(m->p[i][j], 1);
434 else
435 value_set_si(m->p[i][j], 0);
437 value_assign(g, m->p[0][i]);
439 for (; i < m->NbColumns; ++i) {
440 value_assign(old, g);
441 Extended_Euclid(old, m->p[0][i], &c, &b, &g);
442 value_oppose(b, b);
443 for (j = 0; j < m->NbColumns; ++j) {
444 if (j < i) {
445 value_multiply(tmp, m->p[0][j], b);
446 value_division(m->p[i][j], tmp, old);
447 } else if (j == i)
448 value_assign(m->p[i][j], c);
449 else
450 value_set_si(m->p[i][j], 0);
453 ok = value_one_p(g);
454 value_clear(b);
455 value_clear(c);
456 value_clear(g);
457 value_clear(old);
458 value_clear(tmp);
459 return ok;
462 int unimodular_complete(Matrix *M, int row)
464 int r;
465 int ok = 1;
466 Matrix *H, *Q, *U;
468 if (row == 1)
469 return unimodular_complete_1(M);
471 left_hermite(M, &H, &Q, &U);
472 Matrix_Free(U);
473 for (r = 0; ok && r < row; ++r)
474 if (value_notone_p(H->p[r][r]))
475 ok = 0;
476 Matrix_Free(H);
477 for (r = row; r < M->NbRows; ++r)
478 Vector_Copy(Q->p[r], M->p[r], M->NbColumns);
479 Matrix_Free(Q);
480 return ok;
484 * left_hermite may leave positive entries below the main diagonal in H.
485 * This function postprocesses the output of left_hermite to make
486 * the non-zero entries below the main diagonal negative.
488 void neg_left_hermite(Matrix *A, Matrix **H_p, Matrix **Q_p, Matrix **U_p)
490 int row, col, i, j;
491 Matrix *H, *U, *Q;
493 left_hermite(A, &H, &Q, &U);
494 *H_p = H;
495 *Q_p = Q;
496 *U_p = U;
498 for (row = 0, col = 0; col < H->NbColumns; ++col, ++row) {
499 while (value_zero_p(H->p[row][col]))
500 ++row;
501 for (i = 0; i < col; ++i) {
502 if (value_negz_p(H->p[row][i]))
503 continue;
505 /* subtract column col from column i in H and U */
506 for (j = 0; j < H->NbRows; ++j)
507 value_subtract(H->p[j][i], H->p[j][i], H->p[j][col]);
508 for (j = 0; j < U->NbRows; ++j)
509 value_subtract(U->p[j][i], U->p[j][i], U->p[j][col]);
511 /* add row i to row col in Q */
512 for (j = 0; j < Q->NbColumns; ++j)
513 value_addto(Q->p[col][j], Q->p[col][j], Q->p[i][j]);
519 * Returns a full-dimensional polyhedron with the same number
520 * of integer points as P
522 Polyhedron *remove_equalities(Polyhedron *P, unsigned MaxRays)
524 Polyhedron *Q = Polyhedron_Copy(P);
525 unsigned dim = P->Dimension;
526 Matrix *m1, *m2;
527 int i;
529 if (Q->NbEq == 0)
530 return Q;
532 Q = DomainConstraintSimplify(Q, MaxRays);
533 if (emptyQ2(Q))
534 return Q;
536 m1 = Matrix_Alloc(dim, dim);
537 for (i = 0; i < Q->NbEq; ++i)
538 Vector_Copy(Q->Constraint[i]+1, m1->p[i], dim);
540 /* m1 may not be unimodular, but we won't be throwing anything away */
541 unimodular_complete(m1, Q->NbEq);
543 m2 = Matrix_Alloc(dim+1-Q->NbEq, dim+1);
544 for (i = Q->NbEq; i < dim; ++i)
545 Vector_Copy(m1->p[i], m2->p[i-Q->NbEq], dim);
546 value_set_si(m2->p[dim-Q->NbEq][dim], 1);
547 Matrix_Free(m1);
549 P = Polyhedron_Image(Q, m2, MaxRays);
550 Matrix_Free(m2);
551 Polyhedron_Free(Q);
553 return P;
557 * Returns a full-dimensional polyhedron with the same number
558 * of integer points as P
559 * nvar specifies the number of variables
560 * The remaining dimensions are assumed to be parameters
561 * Destroys P
562 * factor is NbEq x (nparam+2) matrix, containing stride constraints
563 * on the parameters; column nparam is the constant;
564 * column nparam+1 is the stride
566 * if factor is NULL, only remove equalities that don't affect
567 * the number of points
569 Polyhedron *remove_equalities_p(Polyhedron *P, unsigned nvar, Matrix **factor,
570 unsigned MaxRays)
572 Value g;
573 Polyhedron *Q;
574 unsigned dim = P->Dimension;
575 Matrix *m1, *m2, *f;
576 int i, j;
578 if (P->NbEq == 0)
579 return P;
581 m1 = Matrix_Alloc(nvar, nvar);
582 P = DomainConstraintSimplify(P, MaxRays);
583 if (factor) {
584 f = Matrix_Alloc(P->NbEq, dim-nvar+2);
585 *factor = f;
587 value_init(g);
588 for (i = 0, j = 0; i < P->NbEq; ++i) {
589 if (First_Non_Zero(P->Constraint[i]+1, nvar) == -1)
590 continue;
592 Vector_Gcd(P->Constraint[i]+1, nvar, &g);
593 if (!factor && value_notone_p(g))
594 continue;
596 if (factor) {
597 Vector_Copy(P->Constraint[i]+1+nvar, f->p[j], dim-nvar+1);
598 value_assign(f->p[j][dim-nvar+1], g);
601 Vector_Copy(P->Constraint[i]+1, m1->p[j], nvar);
603 ++j;
605 value_clear(g);
607 unimodular_complete(m1, j);
609 m2 = Matrix_Alloc(dim+1-j, dim+1);
610 for (i = 0; i < nvar-j ; ++i)
611 Vector_Copy(m1->p[i+j], m2->p[i], nvar);
612 Matrix_Free(m1);
613 for (i = nvar-j; i <= dim-j; ++i)
614 value_set_si(m2->p[i][i+j], 1);
616 Q = Polyhedron_Image(P, m2, MaxRays);
617 Matrix_Free(m2);
618 Polyhedron_Free(P);
620 return Q;
623 void Line_Length(Polyhedron *P, Value *len)
625 Value tmp, pos, neg;
626 int p = 0, n = 0;
627 int i;
629 assert(P->Dimension == 1);
631 value_init(tmp);
632 value_init(pos);
633 value_init(neg);
635 for (i = 0; i < P->NbConstraints; ++i) {
636 value_oppose(tmp, P->Constraint[i][2]);
637 if (value_pos_p(P->Constraint[i][1])) {
638 mpz_cdiv_q(tmp, tmp, P->Constraint[i][1]);
639 if (!p || value_gt(tmp, pos))
640 value_assign(pos, tmp);
641 p = 1;
642 } else if (value_neg_p(P->Constraint[i][1])) {
643 mpz_fdiv_q(tmp, tmp, P->Constraint[i][1]);
644 if (!n || value_lt(tmp, neg))
645 value_assign(neg, tmp);
646 n = 1;
648 if (n && p) {
649 value_subtract(tmp, neg, pos);
650 value_increment(*len, tmp);
651 } else
652 value_set_si(*len, -1);
655 value_clear(tmp);
656 value_clear(pos);
657 value_clear(neg);
661 * Factors the polyhedron P into polyhedra Q_i such that
662 * the number of integer points in P is equal to the product
663 * of the number of integer points in the individual Q_i
665 * If no factors can be found, NULL is returned.
666 * Otherwise, a linked list of the factors is returned.
668 * If there are factors and if T is not NULL, then a matrix will be
669 * returned through T expressing the old variables in terms of the
670 * new variables as they appear in the sequence of factors.
672 * The algorithm works by first computing the Hermite normal form
673 * and then grouping columns linked by one or more constraints together,
674 * where a constraints "links" two or more columns if the constraint
675 * has nonzero coefficients in the columns.
677 Polyhedron* Polyhedron_Factor(Polyhedron *P, unsigned nparam, Matrix **T,
678 unsigned NbMaxRays)
680 int i, j, k;
681 Matrix *M, *H, *Q, *U;
682 int *pos; /* for each column: row position of pivot */
683 int *group; /* group to which a column belongs */
684 int *cnt; /* number of columns in the group */
685 int *rowgroup; /* group to which a constraint belongs */
686 int nvar = P->Dimension - nparam;
687 Polyhedron *F = NULL;
689 if (nvar <= 1)
690 return NULL;
692 NALLOC(pos, nvar);
693 NALLOC(group, nvar);
694 NALLOC(cnt, nvar);
695 NALLOC(rowgroup, P->NbConstraints);
697 M = Matrix_Alloc(P->NbConstraints, nvar);
698 for (i = 0; i < P->NbConstraints; ++i)
699 Vector_Copy(P->Constraint[i]+1, M->p[i], nvar);
700 left_hermite(M, &H, &Q, &U);
701 Matrix_Free(M);
702 Matrix_Free(Q);
704 for (i = 0; i < P->NbConstraints; ++i)
705 rowgroup[i] = -1;
706 for (i = 0, j = 0; i < H->NbColumns; ++i) {
707 for ( ; j < H->NbRows; ++j)
708 if (value_notzero_p(H->p[j][i]))
709 break;
710 assert (j < H->NbRows);
711 pos[i] = j;
713 for (i = 0; i < nvar; ++i) {
714 group[i] = i;
715 cnt[i] = 1;
717 for (i = 0; i < H->NbColumns && cnt[0] < nvar; ++i) {
718 if (rowgroup[pos[i]] == -1)
719 rowgroup[pos[i]] = i;
720 for (j = pos[i]+1; j < H->NbRows; ++j) {
721 if (value_zero_p(H->p[j][i]))
722 continue;
723 if (rowgroup[j] != -1)
724 continue;
725 rowgroup[j] = group[i];
726 for (k = i+1; k < H->NbColumns && j >= pos[k]; ++k) {
727 int g = group[k];
728 while (cnt[g] == 0)
729 g = group[g];
730 group[k] = g;
731 if (group[k] != group[i] && value_notzero_p(H->p[j][k])) {
732 assert(cnt[group[k]] != 0);
733 assert(cnt[group[i]] != 0);
734 if (group[i] < group[k]) {
735 cnt[group[i]] += cnt[group[k]];
736 cnt[group[k]] = 0;
737 group[k] = group[i];
738 } else {
739 cnt[group[k]] += cnt[group[i]];
740 cnt[group[i]] = 0;
741 group[i] = group[k];
748 if (cnt[0] != nvar) {
749 /* Extract out pure context constraints separately */
750 Polyhedron **next = &F;
751 int tot_d = 0;
752 if (T)
753 *T = Matrix_Alloc(nvar, nvar);
754 for (i = nparam ? -1 : 0; i < nvar; ++i) {
755 int d;
757 if (i == -1) {
758 for (j = 0, k = 0; j < P->NbConstraints; ++j)
759 if (rowgroup[j] == -1) {
760 if (First_Non_Zero(P->Constraint[j]+1+nvar,
761 nparam) == -1)
762 rowgroup[j] = -2;
763 else
764 ++k;
766 if (k == 0)
767 continue;
768 d = 0;
769 } else {
770 if (cnt[i] == 0)
771 continue;
772 d = cnt[i];
773 for (j = 0, k = 0; j < P->NbConstraints; ++j)
774 if (rowgroup[j] >= 0 && group[rowgroup[j]] == i) {
775 rowgroup[j] = i;
776 ++k;
780 if (T)
781 for (j = 0; j < nvar; ++j) {
782 int l, m;
783 for (l = 0, m = 0; m < d; ++l) {
784 if (group[l] != i)
785 continue;
786 value_assign((*T)->p[j][tot_d+m++], U->p[j][l]);
790 M = Matrix_Alloc(k, d+nparam+2);
791 for (j = 0, k = 0; j < P->NbConstraints; ++j) {
792 int l, m;
793 if (rowgroup[j] != i)
794 continue;
795 value_assign(M->p[k][0], P->Constraint[j][0]);
796 for (l = 0, m = 0; m < d; ++l) {
797 if (group[l] != i)
798 continue;
799 value_assign(M->p[k][1+m++], H->p[j][l]);
801 Vector_Copy(P->Constraint[j]+1+nvar, M->p[k]+1+m, nparam+1);
802 ++k;
804 *next = Constraints2Polyhedron(M, NbMaxRays);
805 next = &(*next)->next;
806 Matrix_Free(M);
807 tot_d += d;
810 Matrix_Free(U);
811 Matrix_Free(H);
812 free(pos);
813 free(group);
814 free(cnt);
815 free(rowgroup);
816 return F;
820 * Project on final dim dimensions
822 Polyhedron* Polyhedron_Project(Polyhedron *P, int dim)
824 int i;
825 int remove = P->Dimension - dim;
826 Matrix *T;
827 Polyhedron *I;
829 if (P->Dimension == dim)
830 return Polyhedron_Copy(P);
832 T = Matrix_Alloc(dim+1, P->Dimension+1);
833 for (i = 0; i < dim+1; ++i)
834 value_set_si(T->p[i][i+remove], 1);
835 I = Polyhedron_Image(P, T, P->NbConstraints);
836 Matrix_Free(T);
837 return I;
840 /* Constructs a new constraint that ensures that
841 * the first constraint is (strictly) smaller than
842 * the second.
844 static void smaller_constraint(Value *a, Value *b, Value *c, int pos, int shift,
845 int len, int strict, Value *tmp)
847 value_oppose(*tmp, b[pos+1]);
848 value_set_si(c[0], 1);
849 Vector_Combine(a+1+shift, b+1+shift, c+1, *tmp, a[pos+1], len-shift-1);
850 if (strict)
851 value_decrement(c[len-shift-1], c[len-shift-1]);
852 ConstraintSimplify(c, c, len-shift, tmp);
856 /* For each pair of lower and upper bounds on the first variable,
857 * calls fn with the set of constraints on the remaining variables
858 * where these bounds are active, i.e., (stricly) larger/smaller than
859 * the other lower/upper bounds, the lower and upper bound and the
860 * call back data.
862 * If the first variable is equal to an affine combination of the
863 * other variables then fn is called with both lower and upper
864 * pointing to the corresponding equality.
866 * If there is no lower (or upper) bound, then NULL is passed
867 * as the corresponding bound.
869 void for_each_lower_upper_bound(Polyhedron *P,
870 for_each_lower_upper_bound_init init,
871 for_each_lower_upper_bound_fn fn,
872 void *cb_data)
874 unsigned dim = P->Dimension;
875 Matrix *M;
876 int *pos;
877 int i, j, p, n, z;
878 int k, l, k2, l2, q;
879 Value g;
881 if (value_zero_p(P->Constraint[0][0]) &&
882 value_notzero_p(P->Constraint[0][1])) {
883 M = Matrix_Alloc(P->NbConstraints-1, dim-1+2);
884 for (i = 1; i < P->NbConstraints; ++i) {
885 value_assign(M->p[i-1][0], P->Constraint[i][0]);
886 Vector_Copy(P->Constraint[i]+2, M->p[i-1]+1, dim);
888 if (init)
889 init(1, cb_data);
890 fn(M, P->Constraint[0], P->Constraint[0], cb_data);
891 Matrix_Free(M);
892 return;
895 value_init(g);
896 pos = ALLOCN(int, P->NbConstraints);
898 for (i = 0, z = 0; i < P->NbConstraints; ++i)
899 if (value_zero_p(P->Constraint[i][1]))
900 pos[P->NbConstraints-1 - z++] = i;
901 /* put those with positive coefficients first; number: p */
902 for (i = 0, p = 0, n = P->NbConstraints-z-1; i < P->NbConstraints; ++i)
903 if (value_pos_p(P->Constraint[i][1]))
904 pos[p++] = i;
905 else if (value_neg_p(P->Constraint[i][1]))
906 pos[n--] = i;
907 n = P->NbConstraints-z-p;
909 if (init)
910 init(p*n, cb_data);
912 M = Matrix_Alloc((p ? p-1 : 0) + (n ? n-1 : 0) + z + 1, dim-1+2);
913 for (i = 0; i < z; ++i) {
914 value_assign(M->p[i][0], P->Constraint[pos[P->NbConstraints-1 - i]][0]);
915 Vector_Copy(P->Constraint[pos[P->NbConstraints-1 - i]]+2,
916 M->p[i]+1, dim);
918 for (k = p ? 0 : -1; k < p; ++k) {
919 for (k2 = 0; k2 < p; ++k2) {
920 if (k2 == k)
921 continue;
922 q = 1 + z + k2 - (k2 > k);
923 smaller_constraint(
924 P->Constraint[pos[k]],
925 P->Constraint[pos[k2]],
926 M->p[q], 0, 1, dim+2, k2 > k, &g);
928 for (l = n ? p : p-1; l < p+n; ++l) {
929 Value *lower;
930 Value *upper;
931 for (l2 = p; l2 < p+n; ++l2) {
932 if (l2 == l)
933 continue;
934 q = 1 + z + l2-1 - (l2 > l);
935 smaller_constraint(
936 P->Constraint[pos[l2]],
937 P->Constraint[pos[l]],
938 M->p[q], 0, 1, dim+2, l2 > l, &g);
940 if (p && n)
941 smaller_constraint(P->Constraint[pos[k]],
942 P->Constraint[pos[l]],
943 M->p[z], 0, 1, dim+2, 0, &g);
944 lower = p ? P->Constraint[pos[k]] : NULL;
945 upper = n ? P->Constraint[pos[l]] : NULL;
946 fn(M, lower, upper, cb_data);
949 Matrix_Free(M);
951 free(pos);
952 value_clear(g);
955 struct section { Polyhedron * D; evalue E; };
957 struct PLL_data {
958 int nd;
959 unsigned MaxRays;
960 Polyhedron *C;
961 evalue mone;
962 struct section *s;
965 static void PLL_init(unsigned n, void *cb_data)
967 struct PLL_data *data = (struct PLL_data *)cb_data;
969 data->s = ALLOCN(struct section, n);
972 static void PLL_cb(Matrix *M, Value *lower, Value *upper, void *cb_data)
974 struct PLL_data *data = (struct PLL_data *)cb_data;
975 unsigned dim = M->NbColumns-1;
976 Matrix *M2;
977 Polyhedron *T;
978 evalue *L, *U;
980 assert(lower);
981 assert(upper);
983 M2 = Matrix_Copy(M);
984 T = Constraints2Polyhedron(M2, data->MaxRays);
985 Matrix_Free(M2);
986 data->s[data->nd].D = DomainIntersection(T, data->C, data->MaxRays);
987 Domain_Free(T);
989 POL_ENSURE_VERTICES(data->s[data->nd].D);
990 if (emptyQ(data->s[data->nd].D)) {
991 Polyhedron_Free(data->s[data->nd].D);
992 return;
994 L = bv_ceil3(lower+1+1, dim-1+1, lower[0+1], data->s[data->nd].D);
995 U = bv_ceil3(upper+1+1, dim-1+1, upper[0+1], data->s[data->nd].D);
996 eadd(L, U);
997 eadd(&data->mone, U);
998 emul(&data->mone, U);
999 data->s[data->nd].E = *U;
1000 evalue_free(L);
1001 free(U);
1002 ++data->nd;
1005 static evalue *ParamLine_Length_mod(Polyhedron *P, Polyhedron *C, unsigned MaxRays)
1007 unsigned dim = P->Dimension;
1008 unsigned nvar = dim - C->Dimension;
1009 struct PLL_data data;
1010 evalue *F;
1011 int k;
1013 assert(nvar == 1);
1015 value_init(data.mone.d);
1016 evalue_set_si(&data.mone, -1, 1);
1018 data.nd = 0;
1019 data.MaxRays = MaxRays;
1020 data.C = C;
1021 for_each_lower_upper_bound(P, PLL_init, PLL_cb, &data);
1023 F = ALLOC(evalue);
1024 value_init(F->d);
1025 value_set_si(F->d, 0);
1026 F->x.p = new_enode(partition, 2*data.nd, dim-nvar);
1027 for (k = 0; k < data.nd; ++k) {
1028 EVALUE_SET_DOMAIN(F->x.p->arr[2*k], data.s[k].D);
1029 value_clear(F->x.p->arr[2*k+1].d);
1030 F->x.p->arr[2*k+1] = data.s[k].E;
1032 free(data.s);
1034 free_evalue_refs(&data.mone);
1036 return F;
1039 evalue* ParamLine_Length(Polyhedron *P, Polyhedron *C,
1040 struct barvinok_options *options)
1042 evalue* tmp;
1043 tmp = ParamLine_Length_mod(P, C, options->MaxRays);
1044 if (options->lookup_table) {
1045 evalue_mod2table(tmp, C->Dimension);
1046 reduce_evalue(tmp);
1048 return tmp;
1051 Bool isIdentity(Matrix *M)
1053 unsigned i, j;
1054 if (M->NbRows != M->NbColumns)
1055 return False;
1057 for (i = 0;i < M->NbRows; i ++)
1058 for (j = 0; j < M->NbColumns; j ++)
1059 if (i == j) {
1060 if(value_notone_p(M->p[i][j]))
1061 return False;
1062 } else {
1063 if(value_notzero_p(M->p[i][j]))
1064 return False;
1066 return True;
1069 void Param_Polyhedron_Print(FILE* DST, Param_Polyhedron *PP, char **param_names)
1071 Param_Domain *P;
1072 Param_Vertices *V;
1074 for(P=PP->D;P;P=P->next) {
1076 /* prints current val. dom. */
1077 fprintf(DST, "---------------------------------------\n");
1078 fprintf(DST, "Domain :\n");
1079 Print_Domain(DST, P->Domain, param_names);
1081 /* scan the vertices */
1082 fprintf(DST, "Vertices :\n");
1083 FORALL_PVertex_in_ParamPolyhedron(V,P,PP) {
1085 /* prints each vertex */
1086 Print_Vertex(DST, V->Vertex, param_names);
1087 fprintf(DST, "\n");
1089 END_FORALL_PVertex_in_ParamPolyhedron;
1093 void Enumeration_Print(FILE *Dst, Enumeration *en, const char * const *params)
1095 for (; en; en = en->next) {
1096 Print_Domain(Dst, en->ValidityDomain, params);
1097 print_evalue(Dst, &en->EP, params);
1101 void Enumeration_Free(Enumeration *en)
1103 Enumeration *ee;
1105 while( en )
1107 free_evalue_refs( &(en->EP) );
1108 Domain_Free( en->ValidityDomain );
1109 ee = en ->next;
1110 free( en );
1111 en = ee;
1115 void Enumeration_mod2table(Enumeration *en, unsigned nparam)
1117 for (; en; en = en->next) {
1118 evalue_mod2table(&en->EP, nparam);
1119 reduce_evalue(&en->EP);
1123 size_t Enumeration_size(Enumeration *en)
1125 size_t s = 0;
1127 for (; en; en = en->next) {
1128 s += domain_size(en->ValidityDomain);
1129 s += evalue_size(&en->EP);
1131 return s;
1134 void Free_ParamNames(char **params, int m)
1136 while (--m >= 0)
1137 free(params[m]);
1138 free(params);
1141 /* Check whether every set in D2 is included in some set of D1 */
1142 int DomainIncludes(Polyhedron *D1, Polyhedron *D2)
1144 for ( ; D2; D2 = D2->next) {
1145 Polyhedron *P1;
1146 for (P1 = D1; P1; P1 = P1->next)
1147 if (PolyhedronIncludes(P1, D2))
1148 break;
1149 if (!P1)
1150 return 0;
1152 return 1;
1155 int line_minmax(Polyhedron *I, Value *min, Value *max)
1157 int i;
1159 if (I->NbEq >= 1) {
1160 value_oppose(I->Constraint[0][2], I->Constraint[0][2]);
1161 /* There should never be a remainder here */
1162 if (value_pos_p(I->Constraint[0][1]))
1163 mpz_fdiv_q(*min, I->Constraint[0][2], I->Constraint[0][1]);
1164 else
1165 mpz_fdiv_q(*min, I->Constraint[0][2], I->Constraint[0][1]);
1166 value_assign(*max, *min);
1167 } else for (i = 0; i < I->NbConstraints; ++i) {
1168 if (value_zero_p(I->Constraint[i][1])) {
1169 Polyhedron_Free(I);
1170 return 0;
1173 value_oppose(I->Constraint[i][2], I->Constraint[i][2]);
1174 if (value_pos_p(I->Constraint[i][1]))
1175 mpz_cdiv_q(*min, I->Constraint[i][2], I->Constraint[i][1]);
1176 else
1177 mpz_fdiv_q(*max, I->Constraint[i][2], I->Constraint[i][1]);
1179 Polyhedron_Free(I);
1180 return 1;
1183 /**
1185 PROCEDURES TO COMPUTE ENUMERATION. recursive procedure, recurse for
1186 each imbriquation
1188 @param pos index position of current loop index (1..hdim-1)
1189 @param P loop domain
1190 @param context context values for fixed indices
1191 @param exist number of existential variables
1192 @return the number of integer points in this
1193 polyhedron
1196 void count_points_e (int pos, Polyhedron *P, int exist, int nparam,
1197 Value *context, Value *res)
1199 Value LB, UB, k, c;
1201 if (emptyQ(P)) {
1202 value_set_si(*res, 0);
1203 return;
1206 if (!exist) {
1207 count_points(pos, P, context, res);
1208 return;
1211 value_init(LB); value_init(UB); value_init(k);
1212 value_set_si(LB,0);
1213 value_set_si(UB,0);
1215 if (lower_upper_bounds(pos,P,context,&LB,&UB) !=0) {
1216 /* Problem if UB or LB is INFINITY */
1217 value_clear(LB); value_clear(UB); value_clear(k);
1218 if (pos > P->Dimension - nparam - exist)
1219 value_set_si(*res, 1);
1220 else
1221 value_set_si(*res, -1);
1222 return;
1225 #ifdef EDEBUG1
1226 if (!P->next) {
1227 int i;
1228 for (value_assign(k,LB); value_le(k,UB); value_increment(k,k)) {
1229 fprintf(stderr, "(");
1230 for (i=1; i<pos; i++) {
1231 value_print(stderr,P_VALUE_FMT,context[i]);
1232 fprintf(stderr,",");
1234 value_print(stderr,P_VALUE_FMT,k);
1235 fprintf(stderr,")\n");
1238 #endif
1240 value_set_si(context[pos],0);
1241 if (value_lt(UB,LB)) {
1242 value_clear(LB); value_clear(UB); value_clear(k);
1243 value_set_si(*res, 0);
1244 return;
1246 if (!P->next) {
1247 if (exist)
1248 value_set_si(*res, 1);
1249 else {
1250 value_subtract(k,UB,LB);
1251 value_add_int(k,k,1);
1252 value_assign(*res, k);
1254 value_clear(LB); value_clear(UB); value_clear(k);
1255 return;
1258 /*-----------------------------------------------------------------*/
1259 /* Optimization idea */
1260 /* If inner loops are not a function of k (the current index) */
1261 /* i.e. if P->Constraint[i][pos]==0 for all P following this and */
1262 /* for all i, */
1263 /* Then CNT = (UB-LB+1)*count_points(pos+1, P->next, context) */
1264 /* (skip the for loop) */
1265 /*-----------------------------------------------------------------*/
1267 value_init(c);
1268 value_set_si(*res, 0);
1269 for (value_assign(k,LB);value_le(k,UB);value_increment(k,k)) {
1270 /* Insert k in context */
1271 value_assign(context[pos],k);
1272 count_points_e(pos+1, P->next, exist, nparam, context, &c);
1273 if(value_notmone_p(c))
1274 value_addto(*res, *res, c);
1275 else {
1276 value_set_si(*res, -1);
1277 break;
1279 if (pos > P->Dimension - nparam - exist &&
1280 value_pos_p(*res))
1281 break;
1283 value_clear(c);
1285 #ifdef EDEBUG11
1286 fprintf(stderr,"%d\n",CNT);
1287 #endif
1289 /* Reset context */
1290 value_set_si(context[pos],0);
1291 value_clear(LB); value_clear(UB); value_clear(k);
1292 return;
1293 } /* count_points_e */
1295 int DomainContains(Polyhedron *P, Value *list_args, int len,
1296 unsigned MaxRays, int set)
1298 int i;
1299 Value m;
1301 if (P->Dimension == len)
1302 return in_domain(P, list_args);
1304 assert(set); // assume list_args is large enough
1305 assert((P->Dimension - len) % 2 == 0);
1306 value_init(m);
1307 for (i = 0; i < P->Dimension - len; i += 2) {
1308 int j, k;
1309 for (j = 0 ; j < P->NbEq; ++j)
1310 if (value_notzero_p(P->Constraint[j][1+len+i]))
1311 break;
1312 assert(j < P->NbEq);
1313 value_absolute(m, P->Constraint[j][1+len+i]);
1314 k = First_Non_Zero(P->Constraint[j]+1, len);
1315 assert(k != -1);
1316 assert(First_Non_Zero(P->Constraint[j]+1+k+1, len - k - 1) == -1);
1317 mpz_fdiv_q(list_args[len+i], list_args[k], m);
1318 mpz_fdiv_r(list_args[len+i+1], list_args[k], m);
1320 value_clear(m);
1322 return in_domain(P, list_args);
1325 Polyhedron *DomainConcat(Polyhedron *head, Polyhedron *tail)
1327 Polyhedron *S;
1328 if (!head)
1329 return tail;
1330 for (S = head; S->next; S = S->next)
1332 S->next = tail;
1333 return head;
1336 evalue *barvinok_lexsmaller_ev(Polyhedron *P, Polyhedron *D, unsigned dim,
1337 Polyhedron *C, unsigned MaxRays)
1339 evalue *ranking;
1340 Polyhedron *RC, *RD, *Q;
1341 unsigned nparam = dim + C->Dimension;
1342 unsigned exist;
1343 Polyhedron *CA;
1345 RC = LexSmaller(P, D, dim, C, MaxRays);
1346 RD = RC->next;
1347 RC->next = NULL;
1349 exist = RD->Dimension - nparam - dim;
1350 CA = align_context(RC, RD->Dimension, MaxRays);
1351 Q = DomainIntersection(RD, CA, MaxRays);
1352 Polyhedron_Free(CA);
1353 Domain_Free(RD);
1354 Polyhedron_Free(RC);
1355 RD = Q;
1357 for (Q = RD; Q; Q = Q->next) {
1358 evalue *t;
1359 Polyhedron *next = Q->next;
1360 Q->next = 0;
1362 t = barvinok_enumerate_e(Q, exist, nparam, MaxRays);
1364 if (Q == RD)
1365 ranking = t;
1366 else {
1367 eadd(t, ranking);
1368 evalue_free(t);
1371 Q->next = next;
1374 Domain_Free(RD);
1376 return ranking;
1379 Enumeration *barvinok_lexsmaller(Polyhedron *P, Polyhedron *D, unsigned dim,
1380 Polyhedron *C, unsigned MaxRays)
1382 evalue *EP = barvinok_lexsmaller_ev(P, D, dim, C, MaxRays);
1384 return partition2enumeration(EP);
1387 /* "align" matrix to have nrows by inserting
1388 * the necessary number of rows and an equal number of columns in front
1390 Matrix *align_matrix(Matrix *M, int nrows)
1392 int i;
1393 int newrows = nrows - M->NbRows;
1394 Matrix *M2 = Matrix_Alloc(nrows, newrows + M->NbColumns);
1395 for (i = 0; i < newrows; ++i)
1396 value_set_si(M2->p[i][i], 1);
1397 for (i = 0; i < M->NbRows; ++i)
1398 Vector_Copy(M->p[i], M2->p[newrows+i]+newrows, M->NbColumns);
1399 return M2;
1402 static void print_varlist(FILE *out, int n, char **names)
1404 int i;
1405 fprintf(out, "[");
1406 for (i = 0; i < n; ++i) {
1407 if (i)
1408 fprintf(out, ",");
1409 fprintf(out, "%s", names[i]);
1411 fprintf(out, "]");
1414 static void print_term(FILE *out, Value v, int pos, int dim, int nparam,
1415 char **iter_names, char **param_names, int *first)
1417 if (value_zero_p(v)) {
1418 if (first && *first && pos >= dim + nparam)
1419 fprintf(out, "0");
1420 return;
1423 if (first) {
1424 if (!*first && value_pos_p(v))
1425 fprintf(out, "+");
1426 *first = 0;
1428 if (pos < dim + nparam) {
1429 if (value_mone_p(v))
1430 fprintf(out, "-");
1431 else if (!value_one_p(v))
1432 value_print(out, VALUE_FMT, v);
1433 if (pos < dim)
1434 fprintf(out, "%s", iter_names[pos]);
1435 else
1436 fprintf(out, "%s", param_names[pos-dim]);
1437 } else
1438 value_print(out, VALUE_FMT, v);
1441 char **util_generate_names(int n, const char *prefix)
1443 int i;
1444 int len = (prefix ? strlen(prefix) : 0) + 10;
1445 char **names = ALLOCN(char*, n);
1446 if (!names) {
1447 fprintf(stderr, "ERROR: memory overflow.\n");
1448 exit(1);
1450 for (i = 0; i < n; ++i) {
1451 names[i] = ALLOCN(char, len);
1452 if (!names[i]) {
1453 fprintf(stderr, "ERROR: memory overflow.\n");
1454 exit(1);
1456 if (!prefix)
1457 snprintf(names[i], len, "%d", i);
1458 else
1459 snprintf(names[i], len, "%s%d", prefix, i);
1462 return names;
1465 void util_free_names(int n, char **names)
1467 int i;
1468 for (i = 0; i < n; ++i)
1469 free(names[i]);
1470 free(names);
1473 void Polyhedron_pprint(FILE *out, Polyhedron *P, int dim, int nparam,
1474 char **iter_names, char **param_names)
1476 int i, j;
1477 Value tmp;
1479 assert(dim + nparam == P->Dimension);
1481 value_init(tmp);
1483 fprintf(out, "{ ");
1484 if (nparam) {
1485 print_varlist(out, nparam, param_names);
1486 fprintf(out, " -> ");
1488 print_varlist(out, dim, iter_names);
1489 fprintf(out, " : ");
1491 if (emptyQ2(P))
1492 fprintf(out, "FALSE");
1493 else for (i = 0; i < P->NbConstraints; ++i) {
1494 int first = 1;
1495 int v = First_Non_Zero(P->Constraint[i]+1, P->Dimension);
1496 if (v == -1 && value_pos_p(P->Constraint[i][0]))
1497 continue;
1498 if (i)
1499 fprintf(out, " && ");
1500 if (v == -1 && value_notzero_p(P->Constraint[i][1+P->Dimension]))
1501 fprintf(out, "FALSE");
1502 else if (value_pos_p(P->Constraint[i][v+1])) {
1503 print_term(out, P->Constraint[i][v+1], v, dim, nparam,
1504 iter_names, param_names, NULL);
1505 if (value_zero_p(P->Constraint[i][0]))
1506 fprintf(out, " = ");
1507 else
1508 fprintf(out, " >= ");
1509 for (j = v+1; j <= dim+nparam; ++j) {
1510 value_oppose(tmp, P->Constraint[i][1+j]);
1511 print_term(out, tmp, j, dim, nparam,
1512 iter_names, param_names, &first);
1514 } else {
1515 value_oppose(tmp, P->Constraint[i][1+v]);
1516 print_term(out, tmp, v, dim, nparam,
1517 iter_names, param_names, NULL);
1518 fprintf(out, " <= ");
1519 for (j = v+1; j <= dim+nparam; ++j)
1520 print_term(out, P->Constraint[i][1+j], j, dim, nparam,
1521 iter_names, param_names, &first);
1525 fprintf(out, " }\n");
1527 value_clear(tmp);
1530 /* Construct a cone over P with P placed at x_d = 1, with
1531 * x_d the coordinate of an extra dimension
1533 * It's probably a mistake to depend so much on the internal
1534 * representation. We should probably simply compute the
1535 * vertices/facets first.
1537 Polyhedron *Cone_over_Polyhedron(Polyhedron *P)
1539 unsigned NbConstraints = 0;
1540 unsigned NbRays = 0;
1541 Polyhedron *C;
1542 int i;
1544 if (POL_HAS(P, POL_INEQUALITIES))
1545 NbConstraints = P->NbConstraints + 1;
1546 if (POL_HAS(P, POL_POINTS))
1547 NbRays = P->NbRays + 1;
1549 C = Polyhedron_Alloc(P->Dimension+1, NbConstraints, NbRays);
1550 if (POL_HAS(P, POL_INEQUALITIES)) {
1551 C->NbEq = P->NbEq;
1552 for (i = 0; i < P->NbConstraints; ++i)
1553 Vector_Copy(P->Constraint[i], C->Constraint[i], P->Dimension+2);
1554 /* n >= 0 */
1555 value_set_si(C->Constraint[P->NbConstraints][0], 1);
1556 value_set_si(C->Constraint[P->NbConstraints][1+P->Dimension], 1);
1558 if (POL_HAS(P, POL_POINTS)) {
1559 C->NbBid = P->NbBid;
1560 for (i = 0; i < P->NbRays; ++i)
1561 Vector_Copy(P->Ray[i], C->Ray[i], P->Dimension+2);
1562 /* vertex 0 */
1563 value_set_si(C->Ray[P->NbRays][0], 1);
1564 value_set_si(C->Ray[P->NbRays][1+C->Dimension], 1);
1566 POL_SET(C, POL_VALID);
1567 if (POL_HAS(P, POL_INEQUALITIES))
1568 POL_SET(C, POL_INEQUALITIES);
1569 if (POL_HAS(P, POL_POINTS))
1570 POL_SET(C, POL_POINTS);
1571 if (POL_HAS(P, POL_VERTICES))
1572 POL_SET(C, POL_VERTICES);
1573 return C;
1576 /* Returns a (dim+nparam+1)x((dim-n)+nparam+1) matrix
1577 * mapping the transformed subspace back to the original space.
1578 * n is the number of equalities involving the variables
1579 * (i.e., not purely the parameters).
1580 * The remaining n coordinates in the transformed space would
1581 * have constant (parametric) values and are therefore not
1582 * included in the variables of the new space.
1584 Matrix *compress_variables(Matrix *Equalities, unsigned nparam)
1586 unsigned dim = (Equalities->NbColumns-2) - nparam;
1587 Matrix *M, *H, *Q, *U, *C, *ratH, *invH, *Ul, *T1, *T2, *T;
1588 Value mone;
1589 int n, i, j;
1590 int ok;
1592 for (n = 0; n < Equalities->NbRows; ++n)
1593 if (First_Non_Zero(Equalities->p[n]+1, dim) == -1)
1594 break;
1595 if (n == 0)
1596 return Identity(dim+nparam+1);
1597 value_init(mone);
1598 value_set_si(mone, -1);
1599 M = Matrix_Alloc(n, dim);
1600 C = Matrix_Alloc(n+1, nparam+1);
1601 for (i = 0; i < n; ++i) {
1602 Vector_Copy(Equalities->p[i]+1, M->p[i], dim);
1603 Vector_Scale(Equalities->p[i]+1+dim, C->p[i], mone, nparam+1);
1605 value_set_si(C->p[n][nparam], 1);
1606 left_hermite(M, &H, &Q, &U);
1607 Matrix_Free(M);
1608 Matrix_Free(Q);
1609 value_clear(mone);
1611 ratH = Matrix_Alloc(n+1, n+1);
1612 invH = Matrix_Alloc(n+1, n+1);
1613 for (i = 0; i < n; ++i)
1614 Vector_Copy(H->p[i], ratH->p[i], n);
1615 value_set_si(ratH->p[n][n], 1);
1616 ok = Matrix_Inverse(ratH, invH);
1617 assert(ok);
1618 Matrix_Free(H);
1619 Matrix_Free(ratH);
1620 T1 = Matrix_Alloc(n+1, nparam+1);
1621 Matrix_Product(invH, C, T1);
1622 Matrix_Free(C);
1623 Matrix_Free(invH);
1624 if (value_notone_p(T1->p[n][nparam])) {
1625 for (i = 0; i < n; ++i) {
1626 if (!mpz_divisible_p(T1->p[i][nparam], T1->p[n][nparam])) {
1627 Matrix_Free(T1);
1628 Matrix_Free(U);
1629 return NULL;
1631 /* compress_params should have taken care of this */
1632 for (j = 0; j < nparam; ++j)
1633 assert(mpz_divisible_p(T1->p[i][j], T1->p[n][nparam]));
1634 Vector_AntiScale(T1->p[i], T1->p[i], T1->p[n][nparam], nparam+1);
1636 value_set_si(T1->p[n][nparam], 1);
1638 Ul = Matrix_Alloc(dim+1, n+1);
1639 for (i = 0; i < dim; ++i)
1640 Vector_Copy(U->p[i], Ul->p[i], n);
1641 value_set_si(Ul->p[dim][n], 1);
1642 T2 = Matrix_Alloc(dim+1, nparam+1);
1643 Matrix_Product(Ul, T1, T2);
1644 Matrix_Free(Ul);
1645 Matrix_Free(T1);
1647 T = Matrix_Alloc(dim+nparam+1, (dim-n)+nparam+1);
1648 for (i = 0; i < dim; ++i) {
1649 Vector_Copy(U->p[i]+n, T->p[i], dim-n);
1650 Vector_Copy(T2->p[i], T->p[i]+dim-n, nparam+1);
1652 for (i = 0; i < nparam+1; ++i)
1653 value_set_si(T->p[dim+i][(dim-n)+i], 1);
1654 assert(value_one_p(T2->p[dim][nparam]));
1655 Matrix_Free(U);
1656 Matrix_Free(T2);
1658 return T;
1661 /* Computes the left inverse of an affine embedding M and, if Eq is not NULL,
1662 * the equalities that define the affine subspace onto which M maps
1663 * its argument.
1665 Matrix *left_inverse(Matrix *M, Matrix **Eq)
1667 int i, ok;
1668 Matrix *L, *H, *Q, *U, *ratH, *invH, *Ut, *inv;
1669 Vector *t;
1671 if (M->NbColumns == 1) {
1672 inv = Matrix_Alloc(1, M->NbRows);
1673 value_set_si(inv->p[0][M->NbRows-1], 1);
1674 if (Eq) {
1675 *Eq = Matrix_Alloc(M->NbRows-1, 1+(M->NbRows-1)+1);
1676 for (i = 0; i < M->NbRows-1; ++i) {
1677 value_oppose((*Eq)->p[i][1+i], M->p[M->NbRows-1][0]);
1678 value_assign((*Eq)->p[i][1+(M->NbRows-1)], M->p[i][0]);
1681 return inv;
1683 if (Eq)
1684 *Eq = NULL;
1685 L = Matrix_Alloc(M->NbRows-1, M->NbColumns-1);
1686 for (i = 0; i < L->NbRows; ++i)
1687 Vector_Copy(M->p[i], L->p[i], L->NbColumns);
1688 right_hermite(L, &H, &U, &Q);
1689 Matrix_Free(L);
1690 Matrix_Free(Q);
1691 t = Vector_Alloc(U->NbColumns);
1692 for (i = 0; i < U->NbColumns; ++i)
1693 value_oppose(t->p[i], M->p[i][M->NbColumns-1]);
1694 if (Eq) {
1695 *Eq = Matrix_Alloc(H->NbRows - H->NbColumns, 2 + U->NbColumns);
1696 for (i = 0; i < H->NbRows - H->NbColumns; ++i) {
1697 Vector_Copy(U->p[H->NbColumns+i], (*Eq)->p[i]+1, U->NbColumns);
1698 Inner_Product(U->p[H->NbColumns+i], t->p, U->NbColumns,
1699 (*Eq)->p[i]+1+U->NbColumns);
1702 ratH = Matrix_Alloc(H->NbColumns+1, H->NbColumns+1);
1703 invH = Matrix_Alloc(H->NbColumns+1, H->NbColumns+1);
1704 for (i = 0; i < H->NbColumns; ++i)
1705 Vector_Copy(H->p[i], ratH->p[i], H->NbColumns);
1706 value_set_si(ratH->p[ratH->NbRows-1][ratH->NbColumns-1], 1);
1707 Matrix_Free(H);
1708 ok = Matrix_Inverse(ratH, invH);
1709 assert(ok);
1710 Matrix_Free(ratH);
1711 Ut = Matrix_Alloc(invH->NbRows, U->NbColumns+1);
1712 for (i = 0; i < Ut->NbRows-1; ++i) {
1713 Vector_Copy(U->p[i], Ut->p[i], U->NbColumns);
1714 Inner_Product(U->p[i], t->p, U->NbColumns, &Ut->p[i][Ut->NbColumns-1]);
1716 Matrix_Free(U);
1717 Vector_Free(t);
1718 value_set_si(Ut->p[Ut->NbRows-1][Ut->NbColumns-1], 1);
1719 inv = Matrix_Alloc(invH->NbRows, Ut->NbColumns);
1720 Matrix_Product(invH, Ut, inv);
1721 Matrix_Free(Ut);
1722 Matrix_Free(invH);
1723 return inv;
1726 /* Check whether all rays are revlex positive in the parameters
1728 int Polyhedron_has_revlex_positive_rays(Polyhedron *P, unsigned nparam)
1730 int r;
1731 for (r = 0; r < P->NbRays; ++r) {
1732 int i;
1733 if (value_notzero_p(P->Ray[r][P->Dimension+1]))
1734 continue;
1735 for (i = P->Dimension-1; i >= P->Dimension-nparam; --i) {
1736 if (value_neg_p(P->Ray[r][i+1]))
1737 return 0;
1738 if (value_pos_p(P->Ray[r][i+1]))
1739 break;
1741 /* A ray independent of the parameters */
1742 if (i < P->Dimension-nparam)
1743 return 0;
1745 return 1;
1748 static Polyhedron *Recession_Cone(Polyhedron *P, unsigned nparam, unsigned MaxRays)
1750 int i;
1751 unsigned nvar = P->Dimension - nparam;
1752 Matrix *M = Matrix_Alloc(P->NbConstraints, 1 + nvar + 1);
1753 Polyhedron *R;
1754 for (i = 0; i < P->NbConstraints; ++i)
1755 Vector_Copy(P->Constraint[i], M->p[i], 1+nvar);
1756 R = Constraints2Polyhedron(M, MaxRays);
1757 Matrix_Free(M);
1758 return R;
1761 int Polyhedron_is_unbounded(Polyhedron *P, unsigned nparam, unsigned MaxRays)
1763 int i;
1764 int is_unbounded;
1765 Polyhedron *R = Recession_Cone(P, nparam, MaxRays);
1766 POL_ENSURE_VERTICES(R);
1767 if (R->NbBid == 0)
1768 for (i = 0; i < R->NbRays; ++i)
1769 if (value_zero_p(R->Ray[i][1+R->Dimension]))
1770 break;
1771 is_unbounded = R->NbBid > 0 || i < R->NbRays;
1772 Polyhedron_Free(R);
1773 return is_unbounded;
1776 static void SwapColumns(Value **V, int n, int i, int j)
1778 int r;
1780 for (r = 0; r < n; ++r)
1781 value_swap(V[r][i], V[r][j]);
1784 void Polyhedron_ExchangeColumns(Polyhedron *P, int Column1, int Column2)
1786 SwapColumns(P->Constraint, P->NbConstraints, Column1, Column2);
1787 SwapColumns(P->Ray, P->NbRays, Column1, Column2);
1788 if (P->NbEq) {
1789 Matrix M;
1790 Polyhedron_Matrix_View(P, &M, P->NbConstraints);
1791 Gauss(&M, P->NbEq, P->Dimension+1);
1795 /* perform transposition inline; assumes M is a square matrix */
1796 void Matrix_Transposition(Matrix *M)
1798 int i, j;
1800 assert(M->NbRows == M->NbColumns);
1801 for (i = 0; i < M->NbRows; ++i)
1802 for (j = i+1; j < M->NbColumns; ++j)
1803 value_swap(M->p[i][j], M->p[j][i]);
1806 /* Matrix "view" of first rows rows */
1807 void Polyhedron_Matrix_View(Polyhedron *P, Matrix *M, unsigned rows)
1809 M->NbRows = rows;
1810 M->NbColumns = P->Dimension+2;
1811 M->p_Init = P->p_Init;
1812 M->p = P->Constraint;