use test scripts for performing tests
[barvinok.git] / util.c
blob3a205498864d1913e4af9203eadaaf495c41b549
1 #include <stdlib.h>
2 #include <string.h>
3 #include <assert.h>
4 #include <isl/ctx.h>
5 #include <isl/val.h>
6 #include <isl/val_gmp.h>
7 #include <isl/space.h>
8 #include <isl/set.h>
9 #include <isl_set_polylib.h>
10 #include <barvinok/polylib.h>
11 #include <barvinok/options.h>
12 #include <polylib/ranking.h>
13 #include "lattice_point.h"
15 #define ALLOC(type) (type*)malloc(sizeof(type))
16 #define ALLOCN(type,n) (type*)malloc((n) * sizeof(type))
18 #ifdef __GNUC__
19 #define NALLOC(p,n) p = (typeof(p))malloc((n) * sizeof(*p))
20 #else
21 #define NALLOC(p,n) p = (void *)malloc((n) * sizeof(*p))
22 #endif
24 void manual_count(Polyhedron *P, Value* result)
26 isl_ctx *ctx = isl_ctx_alloc();
27 isl_space *space;
28 isl_set *set;
29 isl_val *v;
30 int nvar = P->Dimension;
32 space = isl_space_set_alloc(ctx, 0, nvar);
33 set = isl_set_new_from_polylib(P, space);
35 v = isl_set_count_val(set);
36 isl_val_get_num_gmp(v, *result);
37 isl_val_free(v);
39 isl_set_free(set);
40 isl_ctx_free(ctx);
42 assert(v);
45 #include <barvinok/evalue.h>
46 #include <barvinok/util.h>
47 #include <barvinok/barvinok.h>
49 /* Return random value between 0 and max-1 inclusive
51 int random_int(int max) {
52 return (int) (((double)(max))*rand()/(RAND_MAX+1.0));
55 Polyhedron *Polyhedron_Read(unsigned MaxRays)
57 int vertices = 0;
58 unsigned NbRows, NbColumns;
59 Matrix *M;
60 Polyhedron *P;
61 char s[128];
63 while (fgets(s, sizeof(s), stdin)) {
64 if (*s == '#')
65 continue;
66 if (strncasecmp(s, "vertices", sizeof("vertices")-1) == 0)
67 vertices = 1;
68 if (sscanf(s, "%u %u", &NbRows, &NbColumns) == 2)
69 break;
71 if (feof(stdin))
72 return NULL;
73 M = Matrix_Alloc(NbRows,NbColumns);
74 Matrix_Read_Input(M);
75 if (vertices)
76 P = Rays2Polyhedron(M, MaxRays);
77 else
78 P = Constraints2Polyhedron(M, MaxRays);
79 Matrix_Free(M);
80 return P;
83 /* Inplace polarization
85 void Polyhedron_Polarize(Polyhedron *P)
87 unsigned NbRows;
88 int i;
89 Value **q;
91 POL_ENSURE_FACETS(P);
92 POL_ENSURE_VERTICES(P);
93 NbRows = P->NbConstraints + P->NbRays;
94 q = (Value **)malloc(NbRows * sizeof(Value *));
95 assert(q);
96 for (i = 0; i < P->NbRays; ++i)
97 q[i] = P->Ray[i];
98 for (; i < NbRows; ++i)
99 q[i] = P->Constraint[i-P->NbRays];
100 P->NbConstraints = NbRows - P->NbConstraints;
101 P->NbRays = NbRows - P->NbRays;
102 free(P->Constraint);
103 P->Constraint = q;
104 P->Ray = q + P->NbConstraints;
108 * Rather general polar
109 * We can optimize it significantly if we assume that
110 * P includes zero
112 * Also, we calculate the polar as defined in Schrijver
113 * The opposite should probably work as well and would
114 * eliminate the need for multiplying by -1
116 Polyhedron* Polyhedron_Polar(Polyhedron *P, unsigned NbMaxRays)
118 int i;
119 Value mone;
120 unsigned dim = P->Dimension + 2;
121 Matrix *M = Matrix_Alloc(P->NbRays, dim);
123 assert(M);
124 value_init(mone);
125 value_set_si(mone, -1);
126 for (i = 0; i < P->NbRays; ++i) {
127 Vector_Scale(P->Ray[i], M->p[i], mone, dim);
128 value_multiply(M->p[i][0], M->p[i][0], mone);
129 value_multiply(M->p[i][dim-1], M->p[i][dim-1], mone);
131 P = Constraints2Polyhedron(M, NbMaxRays);
132 assert(P);
133 Matrix_Free(M);
134 value_clear(mone);
135 return P;
139 * Returns the supporting cone of P at the vertex with index v
141 Polyhedron* supporting_cone(Polyhedron *P, int v)
143 Matrix *M;
144 Value tmp;
145 int i, n, j;
146 unsigned char *supporting = (unsigned char *)malloc(P->NbConstraints);
147 unsigned dim = P->Dimension + 2;
149 assert(v >=0 && v < P->NbRays);
150 assert(value_pos_p(P->Ray[v][dim-1]));
151 assert(supporting);
153 value_init(tmp);
154 for (i = 0, n = 0; i < P->NbConstraints; ++i) {
155 Inner_Product(P->Constraint[i] + 1, P->Ray[v] + 1, dim - 1, &tmp);
156 if ((supporting[i] = value_zero_p(tmp)))
157 ++n;
159 assert(n >= dim - 2);
160 value_clear(tmp);
161 M = Matrix_Alloc(n, dim);
162 assert(M);
163 for (i = 0, j = 0; i < P->NbConstraints; ++i)
164 if (supporting[i]) {
165 value_set_si(M->p[j][dim-1], 0);
166 Vector_Copy(P->Constraint[i], M->p[j++], dim-1);
168 free(supporting);
169 P = Constraints2Polyhedron(M, P->NbRays+1);
170 assert(P);
171 Matrix_Free(M);
172 return P;
175 #define INT_BITS (sizeof(unsigned) * 8)
177 unsigned *supporting_constraints(Matrix *Constraints, Param_Vertices *v, int *n)
179 Value lcm, tmp, tmp2;
180 unsigned dim = Constraints->NbColumns;
181 unsigned nparam = v->Vertex->NbColumns - 2;
182 unsigned nvar = dim - nparam - 2;
183 int len = (Constraints->NbRows+INT_BITS-1)/INT_BITS;
184 unsigned *supporting = (unsigned *)calloc(len, sizeof(unsigned));
185 int i, j;
186 Vector *row;
187 int ix;
188 unsigned bx;
190 assert(supporting);
191 row = Vector_Alloc(nparam+1);
192 assert(row);
193 value_init(lcm);
194 value_init(tmp);
195 value_init(tmp2);
196 value_set_si(lcm, 1);
197 for (i = 0, *n = 0, ix = 0, bx = MSB; i < Constraints->NbRows; ++i) {
198 Vector_Set(row->p, 0, nparam+1);
199 for (j = 0 ; j < nvar; ++j) {
200 value_set_si(tmp, 1);
201 value_assign(tmp2, Constraints->p[i][j+1]);
202 if (value_ne(lcm, v->Vertex->p[j][nparam+1])) {
203 value_assign(tmp, lcm);
204 value_lcm(lcm, lcm, v->Vertex->p[j][nparam+1]);
205 value_division(tmp, lcm, tmp);
206 value_multiply(tmp2, tmp2, lcm);
207 value_division(tmp2, tmp2, v->Vertex->p[j][nparam+1]);
209 Vector_Combine(row->p, v->Vertex->p[j], row->p,
210 tmp, tmp2, nparam+1);
212 value_set_si(tmp, 1);
213 Vector_Combine(row->p, Constraints->p[i]+1+nvar, row->p, tmp, lcm, nparam+1);
214 for (j = 0; j < nparam+1; ++j)
215 if (value_notzero_p(row->p[j]))
216 break;
217 if (j == nparam + 1) {
218 supporting[ix] |= bx;
219 ++*n;
221 NEXT(ix, bx);
223 assert(*n >= nvar);
224 value_clear(tmp);
225 value_clear(tmp2);
226 value_clear(lcm);
227 Vector_Free(row);
229 return supporting;
232 Polyhedron* supporting_cone_p(Polyhedron *P, Param_Vertices *v)
234 Matrix *M;
235 unsigned dim = P->Dimension + 2;
236 unsigned nparam = v->Vertex->NbColumns - 2;
237 unsigned nvar = dim - nparam - 2;
238 int i, n, j;
239 int ix;
240 unsigned bx;
241 unsigned *supporting;
242 Matrix View;
244 Polyhedron_Matrix_View(P, &View, P->NbConstraints);
245 supporting = supporting_constraints(&View, v, &n);
246 M = Matrix_Alloc(n, nvar+2);
247 assert(M);
248 for (i = 0, j = 0, ix = 0, bx = MSB; i < P->NbConstraints; ++i) {
249 if (supporting[ix] & bx) {
250 value_set_si(M->p[j][nvar+1], 0);
251 Vector_Copy(P->Constraint[i], M->p[j++], nvar+1);
253 NEXT(ix, bx);
255 free(supporting);
256 P = Constraints2Polyhedron(M, P->NbRays+1);
257 assert(P);
258 Matrix_Free(M);
259 return P;
262 Polyhedron* triangulate_cone(Polyhedron *P, unsigned NbMaxCons)
264 struct barvinok_options *options = barvinok_options_new_with_defaults();
265 options->MaxRays = NbMaxCons;
266 P = triangulate_cone_with_options(P, options);
267 barvinok_options_free(options);
268 return P;
271 Polyhedron* triangulate_cone_with_options(Polyhedron *P,
272 struct barvinok_options *options)
274 const static int MAX_TRY=10;
275 int i, j, r, n, t;
276 Value tmp;
277 unsigned dim = P->Dimension;
278 Matrix *M = Matrix_Alloc(P->NbRays+1, dim+3);
279 Matrix *M2, *M3;
280 Polyhedron *L, *R, *T;
281 assert(P->NbEq == 0);
283 L = NULL;
284 R = NULL;
285 value_init(tmp);
287 Vector_Set(M->p[0]+1, 0, dim+1);
288 value_set_si(M->p[0][0], 1);
289 value_set_si(M->p[0][dim+2], 1);
290 Vector_Set(M->p[P->NbRays]+1, 0, dim+2);
291 value_set_si(M->p[P->NbRays][0], 1);
292 value_set_si(M->p[P->NbRays][dim+1], 1);
294 for (i = 0, r = 1; i < P->NbRays; ++i) {
295 if (value_notzero_p(P->Ray[i][dim+1]))
296 continue;
297 Vector_Copy(P->Ray[i], M->p[r], dim+1);
298 value_set_si(M->p[r][dim+2], 0);
299 ++r;
302 M2 = Matrix_Alloc(dim+1, dim+2);
304 t = 0;
305 if (options->try_Delaunay_triangulation) {
306 /* Delaunay triangulation */
307 for (r = 1; r < P->NbRays; ++r) {
308 Inner_Product(M->p[r]+1, M->p[r]+1, dim, &tmp);
309 value_assign(M->p[r][dim+1], tmp);
311 M3 = Matrix_Copy(M);
312 L = Rays2Polyhedron(M3, options->MaxRays);
313 Matrix_Free(M3);
314 ++t;
315 } else {
316 try_again:
317 /* Usually R should still be 0 */
318 Domain_Free(R);
319 Polyhedron_Free(L);
320 for (r = 1; r < P->NbRays; ++r) {
321 value_set_si(M->p[r][dim+1], random_int((t+1)*dim*P->NbRays)+1);
323 M3 = Matrix_Copy(M);
324 L = Rays2Polyhedron(M3, options->MaxRays);
325 Matrix_Free(M3);
326 ++t;
328 assert(t <= MAX_TRY);
330 R = NULL;
331 n = 0;
333 POL_ENSURE_FACETS(L);
334 for (i = 0; i < L->NbConstraints; ++i) {
335 /* Ignore perpendicular facets, i.e., facets with 0 z-coordinate */
336 if (value_negz_p(L->Constraint[i][dim+1]))
337 continue;
338 if (value_notzero_p(L->Constraint[i][dim+2]))
339 continue;
340 for (j = 1, r = 1; j < M->NbRows; ++j) {
341 Inner_Product(M->p[j]+1, L->Constraint[i]+1, dim+1, &tmp);
342 if (value_notzero_p(tmp))
343 continue;
344 if (r > dim)
345 goto try_again;
346 Vector_Copy(M->p[j]+1, M2->p[r]+1, dim);
347 value_set_si(M2->p[r][0], 1);
348 value_set_si(M2->p[r][dim+1], 0);
349 ++r;
351 assert(r == dim+1);
352 Vector_Set(M2->p[0]+1, 0, dim);
353 value_set_si(M2->p[0][0], 1);
354 value_set_si(M2->p[0][dim+1], 1);
355 T = Rays2Polyhedron(M2, P->NbConstraints+1);
356 T->next = R;
357 R = T;
358 ++n;
360 Matrix_Free(M2);
362 Polyhedron_Free(L);
363 value_clear(tmp);
364 Matrix_Free(M);
366 return R;
369 void check_triangulization(Polyhedron *P, Polyhedron *T)
371 Polyhedron *C, *D, *E, *F, *G, *U;
372 for (C = T; C; C = C->next) {
373 if (C == T)
374 U = C;
375 else
376 U = DomainConvex(DomainUnion(U, C, 100), 100);
377 for (D = C->next; D; D = D->next) {
378 F = C->next;
379 G = D->next;
380 C->next = NULL;
381 D->next = NULL;
382 E = DomainIntersection(C, D, 600);
383 assert(E->NbRays == 0 || E->NbEq >= 1);
384 Polyhedron_Free(E);
385 C->next = F;
386 D->next = G;
389 assert(PolyhedronIncludes(U, P));
390 assert(PolyhedronIncludes(P, U));
393 /* Computes x, y and g such that g = gcd(a,b) and a*x+b*y = g */
394 void Extended_Euclid(Value a, Value b, Value *x, Value *y, Value *g)
396 Value c, d, e, f, tmp;
398 value_init(c);
399 value_init(d);
400 value_init(e);
401 value_init(f);
402 value_init(tmp);
403 value_absolute(c, a);
404 value_absolute(d, b);
405 value_set_si(e, 1);
406 value_set_si(f, 0);
407 while(value_pos_p(d)) {
408 value_division(tmp, c, d);
409 value_multiply(tmp, tmp, f);
410 value_subtract(e, e, tmp);
411 value_division(tmp, c, d);
412 value_multiply(tmp, tmp, d);
413 value_subtract(c, c, tmp);
414 value_swap(c, d);
415 value_swap(e, f);
417 value_assign(*g, c);
418 if (value_zero_p(a))
419 value_set_si(*x, 0);
420 else if (value_pos_p(a))
421 value_assign(*x, e);
422 else value_oppose(*x, e);
423 if (value_zero_p(b))
424 value_set_si(*y, 0);
425 else {
426 value_multiply(tmp, a, *x);
427 value_subtract(tmp, c, tmp);
428 value_division(*y, tmp, b);
430 value_clear(c);
431 value_clear(d);
432 value_clear(e);
433 value_clear(f);
434 value_clear(tmp);
437 static int unimodular_complete_1(Matrix *m)
439 Value g, b, c, old, tmp;
440 unsigned i, j;
441 int ok;
443 value_init(b);
444 value_init(c);
445 value_init(g);
446 value_init(old);
447 value_init(tmp);
448 value_assign(g, m->p[0][0]);
449 for (i = 1; value_zero_p(g) && i < m->NbColumns; ++i) {
450 for (j = 0; j < m->NbColumns; ++j) {
451 if (j == i-1)
452 value_set_si(m->p[i][j], 1);
453 else
454 value_set_si(m->p[i][j], 0);
456 value_assign(g, m->p[0][i]);
458 for (; i < m->NbColumns; ++i) {
459 value_assign(old, g);
460 Extended_Euclid(old, m->p[0][i], &c, &b, &g);
461 value_oppose(b, b);
462 for (j = 0; j < m->NbColumns; ++j) {
463 if (j < i) {
464 value_multiply(tmp, m->p[0][j], b);
465 value_division(m->p[i][j], tmp, old);
466 } else if (j == i)
467 value_assign(m->p[i][j], c);
468 else
469 value_set_si(m->p[i][j], 0);
472 ok = value_one_p(g);
473 value_clear(b);
474 value_clear(c);
475 value_clear(g);
476 value_clear(old);
477 value_clear(tmp);
478 return ok;
481 int unimodular_complete(Matrix *M, int row)
483 int r;
484 int ok = 1;
485 Matrix *H, *Q, *U;
487 if (row == 1)
488 return unimodular_complete_1(M);
490 left_hermite(M, &H, &Q, &U);
491 Matrix_Free(U);
492 for (r = 0; ok && r < row; ++r)
493 if (value_notone_p(H->p[r][r]))
494 ok = 0;
495 Matrix_Free(H);
496 for (r = row; r < M->NbRows; ++r)
497 Vector_Copy(Q->p[r], M->p[r], M->NbColumns);
498 Matrix_Free(Q);
499 return ok;
503 * left_hermite may leave positive entries below the main diagonal in H.
504 * This function postprocesses the output of left_hermite to make
505 * the non-zero entries below the main diagonal negative.
507 void neg_left_hermite(Matrix *A, Matrix **H_p, Matrix **Q_p, Matrix **U_p)
509 int row, col, i, j;
510 Matrix *H, *U, *Q;
512 left_hermite(A, &H, &Q, &U);
513 *H_p = H;
514 *Q_p = Q;
515 *U_p = U;
517 for (row = 0, col = 0; col < H->NbColumns; ++col, ++row) {
518 while (value_zero_p(H->p[row][col]))
519 ++row;
520 for (i = 0; i < col; ++i) {
521 if (value_negz_p(H->p[row][i]))
522 continue;
524 /* subtract column col from column i in H and U */
525 for (j = 0; j < H->NbRows; ++j)
526 value_subtract(H->p[j][i], H->p[j][i], H->p[j][col]);
527 for (j = 0; j < U->NbRows; ++j)
528 value_subtract(U->p[j][i], U->p[j][i], U->p[j][col]);
530 /* add row i to row col in Q */
531 for (j = 0; j < Q->NbColumns; ++j)
532 value_addto(Q->p[col][j], Q->p[col][j], Q->p[i][j]);
538 * Returns a full-dimensional polyhedron with the same number
539 * of integer points as P
541 Polyhedron *remove_equalities(Polyhedron *P, unsigned MaxRays)
543 Matrix M;
544 Matrix *T;
545 Polyhedron *Q = Polyhedron_Copy(P);
547 if (Q->NbEq == 0)
548 return Q;
550 Q = DomainConstraintSimplify(Q, MaxRays);
551 if (emptyQ2(Q))
552 return Q;
554 Polyhedron_Matrix_View(Q, &M, Q->NbEq);
555 T = compress_variables(&M, 0);
557 if (!T)
558 P = NULL;
559 else {
560 P = Polyhedron_Preimage(Q, T, MaxRays);
561 Matrix_Free(T);
564 Polyhedron_Free(Q);
566 return P;
570 * Returns a full-dimensional polyhedron with the same number
571 * of integer points as P
572 * nvar specifies the number of variables
573 * The remaining dimensions are assumed to be parameters
574 * Destroys P
575 * factor is NbEq x (nparam+2) matrix, containing stride constraints
576 * on the parameters; column nparam is the constant;
577 * column nparam+1 is the stride
579 * if factor is NULL, only remove equalities that don't affect
580 * the number of points
582 Polyhedron *remove_equalities_p(Polyhedron *P, unsigned nvar, Matrix **factor,
583 unsigned MaxRays)
585 Value g;
586 Polyhedron *Q;
587 unsigned dim = P->Dimension;
588 Matrix *m1, *m2, *f;
589 int i, j;
591 if (P->NbEq == 0)
592 return P;
594 m1 = Matrix_Alloc(nvar, nvar);
595 P = DomainConstraintSimplify(P, MaxRays);
596 if (factor) {
597 f = Matrix_Alloc(P->NbEq, dim-nvar+2);
598 *factor = f;
600 value_init(g);
601 for (i = 0, j = 0; i < P->NbEq; ++i) {
602 if (First_Non_Zero(P->Constraint[i]+1, nvar) == -1)
603 continue;
605 Vector_Gcd(P->Constraint[i]+1, nvar, &g);
606 if (!factor && value_notone_p(g))
607 continue;
609 if (factor) {
610 Vector_Copy(P->Constraint[i]+1+nvar, f->p[j], dim-nvar+1);
611 value_assign(f->p[j][dim-nvar+1], g);
614 Vector_Copy(P->Constraint[i]+1, m1->p[j], nvar);
616 ++j;
618 value_clear(g);
620 unimodular_complete(m1, j);
622 m2 = Matrix_Alloc(dim+1-j, dim+1);
623 for (i = 0; i < nvar-j ; ++i)
624 Vector_Copy(m1->p[i+j], m2->p[i], nvar);
625 Matrix_Free(m1);
626 for (i = nvar-j; i <= dim-j; ++i)
627 value_set_si(m2->p[i][i+j], 1);
629 Q = Polyhedron_Image(P, m2, MaxRays);
630 Matrix_Free(m2);
631 Polyhedron_Free(P);
633 return Q;
636 void Line_Length(Polyhedron *P, Value *len)
638 Value tmp, pos, neg;
639 int p = 0, n = 0;
640 int i;
642 assert(P->Dimension == 1);
644 if (P->NbEq > 0) {
645 if (mpz_divisible_p(P->Constraint[0][2], P->Constraint[0][1]))
646 value_set_si(*len, 1);
647 else
648 value_set_si(*len, 0);
649 return;
652 value_init(tmp);
653 value_init(pos);
654 value_init(neg);
656 for (i = 0; i < P->NbConstraints; ++i) {
657 value_oppose(tmp, P->Constraint[i][2]);
658 if (value_pos_p(P->Constraint[i][1])) {
659 mpz_cdiv_q(tmp, tmp, P->Constraint[i][1]);
660 if (!p || value_gt(tmp, pos))
661 value_assign(pos, tmp);
662 p = 1;
663 } else if (value_neg_p(P->Constraint[i][1])) {
664 mpz_fdiv_q(tmp, tmp, P->Constraint[i][1]);
665 if (!n || value_lt(tmp, neg))
666 value_assign(neg, tmp);
667 n = 1;
669 if (n && p) {
670 value_subtract(tmp, neg, pos);
671 value_increment(*len, tmp);
672 } else
673 value_set_si(*len, -1);
676 value_clear(tmp);
677 value_clear(pos);
678 value_clear(neg);
681 /* Update group[k] to the group column k belongs to.
682 * When merging two groups, only the group of the current
683 * group leader is changed. Here we change the group of
684 * the other members to also point to the group that the
685 * old group leader now points to.
687 static void update_group(int *group, int *cnt, int k)
689 int g = group[k];
690 while (cnt[g] == 0)
691 g = group[g];
692 group[k] = g;
696 * Factors the polyhedron P into polyhedra Q_i such that
697 * the number of integer points in P is equal to the product
698 * of the number of integer points in the individual Q_i
700 * If no factors can be found, NULL is returned.
701 * Otherwise, a linked list of the factors is returned.
703 * If there are factors and if T is not NULL, then a matrix will be
704 * returned through T expressing the old variables in terms of the
705 * new variables as they appear in the sequence of factors.
707 * The algorithm works by first computing the Hermite normal form
708 * and then grouping columns linked by one or more constraints together,
709 * where a constraints "links" two or more columns if the constraint
710 * has nonzero coefficients in the columns.
712 Polyhedron* Polyhedron_Factor(Polyhedron *P, unsigned nparam, Matrix **T,
713 unsigned NbMaxRays)
715 int i, j, k;
716 Matrix *M, *H, *Q, *U;
717 int *pos; /* for each column: row position of pivot */
718 int *group; /* group to which a column belongs */
719 int *cnt; /* number of columns in the group */
720 int *rowgroup; /* group to which a constraint belongs */
721 int nvar = P->Dimension - nparam;
722 Polyhedron *F = NULL;
724 if (nvar <= 1)
725 return NULL;
727 NALLOC(pos, nvar);
728 NALLOC(group, nvar);
729 NALLOC(cnt, nvar);
730 NALLOC(rowgroup, P->NbConstraints);
732 M = Matrix_Alloc(P->NbConstraints, nvar);
733 for (i = 0; i < P->NbConstraints; ++i)
734 Vector_Copy(P->Constraint[i]+1, M->p[i], nvar);
735 left_hermite(M, &H, &Q, &U);
736 Matrix_Free(M);
737 Matrix_Free(Q);
739 for (i = 0; i < P->NbConstraints; ++i)
740 rowgroup[i] = -1;
741 for (i = 0, j = 0; i < H->NbColumns; ++i) {
742 for ( ; j < H->NbRows; ++j)
743 if (value_notzero_p(H->p[j][i]))
744 break;
745 pos[i] = j;
747 for (i = 0; i < nvar; ++i) {
748 group[i] = i;
749 cnt[i] = 1;
751 for (i = 0; i < H->NbColumns && cnt[0] < nvar; ++i) {
752 if (pos[i] == H->NbRows)
753 continue; /* A line direction */
754 if (rowgroup[pos[i]] == -1)
755 rowgroup[pos[i]] = i;
756 for (j = pos[i]+1; j < H->NbRows; ++j) {
757 if (value_zero_p(H->p[j][i]))
758 continue;
759 if (rowgroup[j] != -1)
760 continue;
761 rowgroup[j] = group[i];
762 for (k = i+1; k < H->NbColumns && j >= pos[k]; ++k) {
763 update_group(group, cnt, k);
764 update_group(group, cnt, i);
765 if (group[k] != group[i] && value_notzero_p(H->p[j][k])) {
766 assert(cnt[group[k]] != 0);
767 assert(cnt[group[i]] != 0);
768 if (group[i] < group[k]) {
769 cnt[group[i]] += cnt[group[k]];
770 cnt[group[k]] = 0;
771 group[group[k]] = group[i];
772 } else {
773 cnt[group[k]] += cnt[group[i]];
774 cnt[group[i]] = 0;
775 group[group[i]] = group[k];
781 for (i = 1; i < nvar; ++i)
782 update_group(group, cnt, i);
784 if (cnt[0] != nvar) {
785 /* Extract out pure context constraints separately */
786 Polyhedron **next = &F;
787 int tot_d = 0;
788 if (T)
789 *T = Matrix_Alloc(nvar, nvar);
790 for (i = nparam ? -1 : 0; i < nvar; ++i) {
791 int d;
793 if (i == -1) {
794 for (j = 0, k = 0; j < P->NbConstraints; ++j)
795 if (rowgroup[j] == -1) {
796 if (First_Non_Zero(P->Constraint[j]+1+nvar,
797 nparam) == -1)
798 rowgroup[j] = -2;
799 else
800 ++k;
802 if (k == 0)
803 continue;
804 d = 0;
805 } else {
806 if (cnt[i] == 0)
807 continue;
808 d = cnt[i];
809 for (j = 0, k = 0; j < P->NbConstraints; ++j)
810 if (rowgroup[j] >= 0 && group[rowgroup[j]] == i) {
811 rowgroup[j] = i;
812 ++k;
816 if (T)
817 for (j = 0; j < nvar; ++j) {
818 int l, m;
819 for (l = 0, m = 0; m < d; ++l) {
820 if (group[l] != i)
821 continue;
822 value_assign((*T)->p[j][tot_d+m++], U->p[j][l]);
826 M = Matrix_Alloc(k, d+nparam+2);
827 for (j = 0, k = 0; j < P->NbConstraints; ++j) {
828 int l, m;
829 if (rowgroup[j] != i)
830 continue;
831 value_assign(M->p[k][0], P->Constraint[j][0]);
832 for (l = 0, m = 0; m < d; ++l) {
833 if (group[l] != i)
834 continue;
835 value_assign(M->p[k][1+m++], H->p[j][l]);
837 Vector_Copy(P->Constraint[j]+1+nvar, M->p[k]+1+m, nparam+1);
838 ++k;
840 *next = Constraints2Polyhedron(M, NbMaxRays);
841 next = &(*next)->next;
842 Matrix_Free(M);
843 tot_d += d;
846 Matrix_Free(U);
847 Matrix_Free(H);
848 free(pos);
849 free(group);
850 free(cnt);
851 free(rowgroup);
852 return F;
855 /* Computes the intersection of the contexts of a list of factors */
856 Polyhedron *Factor_Context(Polyhedron *F, unsigned nparam, unsigned MaxRays)
858 Polyhedron *Q;
859 Polyhedron *C = NULL;
861 for (Q = F; Q; Q = Q->next) {
862 Polyhedron *QC = Q;
863 Polyhedron *next = Q->next;
864 Q->next = NULL;
866 if (Q->Dimension != nparam)
867 QC = Polyhedron_Project(Q, nparam);
869 if (!C)
870 C = Q == QC ? Polyhedron_Copy(QC) : QC;
871 else {
872 Polyhedron *C2 = C;
873 C = DomainIntersection(C, QC, MaxRays);
874 Polyhedron_Free(C2);
875 if (QC != Q)
876 Polyhedron_Free(QC);
878 Q->next = next;
880 return C;
884 * Project on final dim dimensions
886 Polyhedron* Polyhedron_Project(Polyhedron *P, int dim)
888 int i;
889 int remove = P->Dimension - dim;
890 Matrix *T;
891 Polyhedron *I;
893 if (P->Dimension == dim)
894 return Polyhedron_Copy(P);
896 T = Matrix_Alloc(dim+1, P->Dimension+1);
897 for (i = 0; i < dim+1; ++i)
898 value_set_si(T->p[i][i+remove], 1);
899 I = Polyhedron_Image(P, T, P->NbConstraints);
900 Matrix_Free(T);
901 return I;
904 /* Constructs a new constraint that ensures that
905 * the first constraint is (strictly) smaller than
906 * the second.
908 static void smaller_constraint(Value *a, Value *b, Value *c, int pos, int shift,
909 int len, int strict, Value *tmp)
911 value_oppose(*tmp, b[pos+1]);
912 value_set_si(c[0], 1);
913 Vector_Combine(a+1+shift, b+1+shift, c+1, *tmp, a[pos+1], len-shift-1);
914 if (strict)
915 value_decrement(c[len-shift-1], c[len-shift-1]);
916 ConstraintSimplify(c, c, len-shift, tmp);
920 /* For each pair of lower and upper bounds on the first variable,
921 * calls fn with the set of constraints on the remaining variables
922 * where these bounds are active, i.e., (stricly) larger/smaller than
923 * the other lower/upper bounds, the lower and upper bound and the
924 * call back data.
926 * If the first variable is equal to an affine combination of the
927 * other variables then fn is called with both lower and upper
928 * pointing to the corresponding equality.
930 * If there is no lower (or upper) bound, then NULL is passed
931 * as the corresponding bound.
933 void for_each_lower_upper_bound(Polyhedron *P,
934 for_each_lower_upper_bound_init init,
935 for_each_lower_upper_bound_fn fn,
936 void *cb_data)
938 unsigned dim = P->Dimension;
939 Matrix *M;
940 int *pos;
941 int i, p, n, z;
942 int k, l, k2, l2, q;
943 Value g;
945 if (value_zero_p(P->Constraint[0][0]) &&
946 value_notzero_p(P->Constraint[0][1])) {
947 M = Matrix_Alloc(P->NbConstraints-1, dim-1+2);
948 for (i = 1; i < P->NbConstraints; ++i) {
949 value_assign(M->p[i-1][0], P->Constraint[i][0]);
950 Vector_Copy(P->Constraint[i]+2, M->p[i-1]+1, dim);
952 if (init)
953 init(1, cb_data);
954 fn(M, P->Constraint[0], P->Constraint[0], cb_data);
955 Matrix_Free(M);
956 return;
959 value_init(g);
960 pos = ALLOCN(int, P->NbConstraints);
962 for (i = 0, z = 0; i < P->NbConstraints; ++i)
963 if (value_zero_p(P->Constraint[i][1]))
964 pos[P->NbConstraints-1 - z++] = i;
965 /* put those with positive coefficients first; number: p */
966 for (i = 0, p = 0, n = P->NbConstraints-z-1; i < P->NbConstraints; ++i)
967 if (value_pos_p(P->Constraint[i][1]))
968 pos[p++] = i;
969 else if (value_neg_p(P->Constraint[i][1]))
970 pos[n--] = i;
971 n = P->NbConstraints-z-p;
973 if (init)
974 init(p*n, cb_data);
976 M = Matrix_Alloc((p ? p-1 : 0) + (n ? n-1 : 0) + z + 1, dim-1+2);
977 for (i = 0; i < z; ++i) {
978 value_assign(M->p[i][0], P->Constraint[pos[P->NbConstraints-1 - i]][0]);
979 Vector_Copy(P->Constraint[pos[P->NbConstraints-1 - i]]+2,
980 M->p[i]+1, dim);
982 for (k = p ? 0 : -1; k < p; ++k) {
983 for (k2 = 0; k2 < p; ++k2) {
984 if (k2 == k)
985 continue;
986 q = 1 + z + k2 - (k2 > k);
987 smaller_constraint(
988 P->Constraint[pos[k]],
989 P->Constraint[pos[k2]],
990 M->p[q], 0, 1, dim+2, k2 > k, &g);
992 for (l = n ? p : p-1; l < p+n; ++l) {
993 Value *lower;
994 Value *upper;
995 for (l2 = p; l2 < p+n; ++l2) {
996 if (l2 == l)
997 continue;
998 q = 1 + z + l2-1 - (l2 > l);
999 smaller_constraint(
1000 P->Constraint[pos[l2]],
1001 P->Constraint[pos[l]],
1002 M->p[q], 0, 1, dim+2, l2 > l, &g);
1004 if (p && n)
1005 smaller_constraint(P->Constraint[pos[k]],
1006 P->Constraint[pos[l]],
1007 M->p[z], 0, 1, dim+2, 0, &g);
1008 lower = p ? P->Constraint[pos[k]] : NULL;
1009 upper = n ? P->Constraint[pos[l]] : NULL;
1010 fn(M, lower, upper, cb_data);
1013 Matrix_Free(M);
1015 free(pos);
1016 value_clear(g);
1019 struct section { Polyhedron * D; evalue E; };
1021 struct PLL_data {
1022 int nd;
1023 unsigned MaxRays;
1024 Polyhedron *C;
1025 evalue mone;
1026 struct section *s;
1029 static void PLL_init(unsigned n, void *cb_data)
1031 struct PLL_data *data = (struct PLL_data *)cb_data;
1033 data->s = ALLOCN(struct section, n);
1036 /* Computes ceil(-coef/abs(d)) */
1037 static evalue* bv_ceil3(Value *coef, int len, Value d, Polyhedron *P)
1039 Value t;
1040 evalue *EP;
1041 Vector *val = Vector_Alloc(len);
1043 value_init(t);
1044 Vector_Oppose(coef, val->p, len);
1045 value_absolute(t, d);
1047 EP = ceiling(val->p, t, len-1, P);
1049 value_clear(t);
1050 Vector_Free(val);
1052 return EP;
1055 static void PLL_cb(Matrix *M, Value *lower, Value *upper, void *cb_data)
1057 struct PLL_data *data = (struct PLL_data *)cb_data;
1058 unsigned dim = M->NbColumns-1;
1059 Matrix *M2;
1060 Polyhedron *T;
1061 evalue *L, *U;
1063 assert(lower);
1064 assert(upper);
1066 M2 = Matrix_Copy(M);
1067 T = Constraints2Polyhedron(M2, data->MaxRays);
1068 Matrix_Free(M2);
1069 data->s[data->nd].D = DomainIntersection(T, data->C, data->MaxRays);
1070 Domain_Free(T);
1072 POL_ENSURE_VERTICES(data->s[data->nd].D);
1073 if (emptyQ(data->s[data->nd].D)) {
1074 Polyhedron_Free(data->s[data->nd].D);
1075 return;
1077 L = bv_ceil3(lower+1+1, dim-1+1, lower[0+1], data->s[data->nd].D);
1078 U = bv_ceil3(upper+1+1, dim-1+1, upper[0+1], data->s[data->nd].D);
1079 eadd(L, U);
1080 eadd(&data->mone, U);
1081 emul(&data->mone, U);
1082 data->s[data->nd].E = *U;
1083 evalue_free(L);
1084 free(U);
1085 ++data->nd;
1088 static evalue *ParamLine_Length_mod(Polyhedron *P, Polyhedron *C, unsigned MaxRays)
1090 unsigned dim = P->Dimension;
1091 unsigned nvar = dim - C->Dimension;
1092 struct PLL_data data;
1093 evalue *F;
1094 int k;
1096 assert(nvar == 1);
1098 value_init(data.mone.d);
1099 evalue_set_si(&data.mone, -1, 1);
1101 data.nd = 0;
1102 data.MaxRays = MaxRays;
1103 data.C = C;
1104 for_each_lower_upper_bound(P, PLL_init, PLL_cb, &data);
1106 free_evalue_refs(&data.mone);
1108 if (data.nd == 0) {
1109 free(data.s);
1110 return evalue_zero();
1113 F = ALLOC(evalue);
1114 value_init(F->d);
1115 value_set_si(F->d, 0);
1116 F->x.p = new_enode(partition, 2*data.nd, dim-nvar);
1117 for (k = 0; k < data.nd; ++k) {
1118 EVALUE_SET_DOMAIN(F->x.p->arr[2*k], data.s[k].D);
1119 value_clear(F->x.p->arr[2*k+1].d);
1120 F->x.p->arr[2*k+1] = data.s[k].E;
1122 free(data.s);
1124 return F;
1127 evalue* ParamLine_Length(Polyhedron *P, Polyhedron *C,
1128 struct barvinok_options *options)
1130 evalue* tmp;
1131 tmp = ParamLine_Length_mod(P, C, options->MaxRays);
1132 if (options->lookup_table) {
1133 evalue_mod2table(tmp, C->Dimension);
1134 reduce_evalue(tmp);
1136 return tmp;
1139 Bool isIdentity(Matrix *M)
1141 unsigned i, j;
1142 if (M->NbRows != M->NbColumns)
1143 return False;
1145 for (i = 0;i < M->NbRows; i ++)
1146 for (j = 0; j < M->NbColumns; j ++)
1147 if (i == j) {
1148 if(value_notone_p(M->p[i][j]))
1149 return False;
1150 } else {
1151 if(value_notzero_p(M->p[i][j]))
1152 return False;
1154 return True;
1157 void Param_Polyhedron_Print(FILE* DST, Param_Polyhedron *PP,
1158 const char **param_names)
1160 Param_Domain *P;
1161 Param_Vertices *V;
1163 for(P=PP->D;P;P=P->next) {
1165 /* prints current val. dom. */
1166 fprintf(DST, "---------------------------------------\n");
1167 fprintf(DST, "Domain :\n");
1168 Print_Domain(DST, P->Domain, param_names);
1170 /* scan the vertices */
1171 fprintf(DST, "Vertices :\n");
1172 FORALL_PVertex_in_ParamPolyhedron(V,P,PP) {
1174 /* prints each vertex */
1175 Print_Vertex(DST, V->Vertex, param_names);
1176 fprintf(DST, "\n");
1178 END_FORALL_PVertex_in_ParamPolyhedron;
1182 void Enumeration_Print(FILE *Dst, Enumeration *en, const char **params)
1184 for (; en; en = en->next) {
1185 Print_Domain(Dst, en->ValidityDomain, params);
1186 print_evalue(Dst, &en->EP, params);
1190 void Enumeration_Free(Enumeration *en)
1192 Enumeration *ee;
1194 while( en )
1196 free_evalue_refs( &(en->EP) );
1197 Domain_Free( en->ValidityDomain );
1198 ee = en ->next;
1199 free( en );
1200 en = ee;
1204 void Enumeration_mod2table(Enumeration *en, unsigned nparam)
1206 for (; en; en = en->next) {
1207 evalue_mod2table(&en->EP, nparam);
1208 reduce_evalue(&en->EP);
1212 size_t Enumeration_size(Enumeration *en)
1214 size_t s = 0;
1216 for (; en; en = en->next) {
1217 s += domain_size(en->ValidityDomain);
1218 s += evalue_size(&en->EP);
1220 return s;
1223 /* Check whether every set in D2 is included in some set of D1 */
1224 int DomainIncludes(Polyhedron *D1, Polyhedron *D2)
1226 for ( ; D2; D2 = D2->next) {
1227 Polyhedron *P1;
1228 for (P1 = D1; P1; P1 = P1->next)
1229 if (PolyhedronIncludes(P1, D2))
1230 break;
1231 if (!P1)
1232 return 0;
1234 return 1;
1237 int line_minmax(Polyhedron *I, Value *min, Value *max)
1239 int i;
1241 if (I->NbEq >= 1) {
1242 value_oppose(I->Constraint[0][2], I->Constraint[0][2]);
1243 /* There should never be a remainder here */
1244 if (value_pos_p(I->Constraint[0][1]))
1245 mpz_fdiv_q(*min, I->Constraint[0][2], I->Constraint[0][1]);
1246 else
1247 mpz_fdiv_q(*min, I->Constraint[0][2], I->Constraint[0][1]);
1248 value_assign(*max, *min);
1249 } else for (i = 0; i < I->NbConstraints; ++i) {
1250 if (value_zero_p(I->Constraint[i][1])) {
1251 Polyhedron_Free(I);
1252 return 0;
1255 value_oppose(I->Constraint[i][2], I->Constraint[i][2]);
1256 if (value_pos_p(I->Constraint[i][1]))
1257 mpz_cdiv_q(*min, I->Constraint[i][2], I->Constraint[i][1]);
1258 else
1259 mpz_fdiv_q(*max, I->Constraint[i][2], I->Constraint[i][1]);
1261 Polyhedron_Free(I);
1262 return 1;
1265 int DomainContains(Polyhedron *P, Value *list_args, int len,
1266 unsigned MaxRays, int set)
1268 int i;
1269 Value m;
1271 if (P->Dimension == len)
1272 return in_domain(P, list_args);
1274 assert(set); // assume list_args is large enough
1275 assert((P->Dimension - len) % 2 == 0);
1276 value_init(m);
1277 for (i = 0; i < P->Dimension - len; i += 2) {
1278 int j, k;
1279 for (j = 0 ; j < P->NbEq; ++j)
1280 if (value_notzero_p(P->Constraint[j][1+len+i]))
1281 break;
1282 assert(j < P->NbEq);
1283 value_absolute(m, P->Constraint[j][1+len+i]);
1284 k = First_Non_Zero(P->Constraint[j]+1, len);
1285 assert(k != -1);
1286 assert(First_Non_Zero(P->Constraint[j]+1+k+1, len - k - 1) == -1);
1287 mpz_fdiv_q(list_args[len+i], list_args[k], m);
1288 mpz_fdiv_r(list_args[len+i+1], list_args[k], m);
1290 value_clear(m);
1292 return in_domain(P, list_args);
1295 Polyhedron *DomainConcat(Polyhedron *head, Polyhedron *tail)
1297 Polyhedron *S;
1298 if (!head)
1299 return tail;
1300 for (S = head; S->next; S = S->next)
1302 S->next = tail;
1303 return head;
1306 evalue *barvinok_lexsmaller_ev(Polyhedron *P, Polyhedron *D, unsigned dim,
1307 Polyhedron *C, unsigned MaxRays)
1309 evalue *ranking;
1310 Polyhedron *RC, *RD, *Q;
1311 unsigned nparam = dim + C->Dimension;
1312 unsigned exist;
1313 Polyhedron *CA;
1315 RC = LexSmaller(P, D, dim, C, MaxRays);
1316 RD = RC->next;
1317 RC->next = NULL;
1319 exist = RD->Dimension - nparam - dim;
1320 CA = align_context(RC, RD->Dimension, MaxRays);
1321 Q = DomainIntersection(RD, CA, MaxRays);
1322 Polyhedron_Free(CA);
1323 Domain_Free(RD);
1324 Polyhedron_Free(RC);
1325 RD = Q;
1327 for (Q = RD; Q; Q = Q->next) {
1328 evalue *t;
1329 Polyhedron *next = Q->next;
1330 Q->next = 0;
1332 t = barvinok_enumerate_e(Q, exist, nparam, MaxRays);
1334 if (Q == RD)
1335 ranking = t;
1336 else {
1337 eadd(t, ranking);
1338 evalue_free(t);
1341 Q->next = next;
1344 Domain_Free(RD);
1346 return ranking;
1349 Enumeration *barvinok_lexsmaller(Polyhedron *P, Polyhedron *D, unsigned dim,
1350 Polyhedron *C, unsigned MaxRays)
1352 evalue *EP = barvinok_lexsmaller_ev(P, D, dim, C, MaxRays);
1354 return partition2enumeration(EP);
1357 /* "align" matrix to have nrows by inserting
1358 * the necessary number of rows and an equal number of columns in front
1360 Matrix *align_matrix(Matrix *M, int nrows)
1362 int i;
1363 int newrows = nrows - M->NbRows;
1364 Matrix *M2 = Matrix_Alloc(nrows, newrows + M->NbColumns);
1365 for (i = 0; i < newrows; ++i)
1366 value_set_si(M2->p[i][i], 1);
1367 for (i = 0; i < M->NbRows; ++i)
1368 Vector_Copy(M->p[i], M2->p[newrows+i]+newrows, M->NbColumns);
1369 return M2;
1372 static void print_varlist(FILE *out, int n, char **names)
1374 int i;
1375 fprintf(out, "[");
1376 for (i = 0; i < n; ++i) {
1377 if (i)
1378 fprintf(out, ",");
1379 fprintf(out, "%s", names[i]);
1381 fprintf(out, "]");
1384 static void print_term(FILE *out, Value v, int pos, int dim, int nparam,
1385 char **iter_names, char **param_names, int *first)
1387 if (value_zero_p(v)) {
1388 if (first && *first && pos >= dim + nparam)
1389 fprintf(out, "0");
1390 return;
1393 if (first) {
1394 if (!*first && value_pos_p(v))
1395 fprintf(out, "+");
1396 *first = 0;
1398 if (pos < dim + nparam) {
1399 if (value_mone_p(v))
1400 fprintf(out, "-");
1401 else if (!value_one_p(v))
1402 value_print(out, VALUE_FMT, v);
1403 if (pos < dim)
1404 fprintf(out, "%s", iter_names[pos]);
1405 else
1406 fprintf(out, "%s", param_names[pos-dim]);
1407 } else
1408 value_print(out, VALUE_FMT, v);
1411 char **util_generate_names(int n, const char *prefix)
1413 int i;
1414 int len = (prefix ? strlen(prefix) : 0) + 10;
1415 char **names = ALLOCN(char*, n);
1416 if (!names) {
1417 fprintf(stderr, "ERROR: memory overflow.\n");
1418 exit(1);
1420 for (i = 0; i < n; ++i) {
1421 names[i] = ALLOCN(char, len);
1422 if (!names[i]) {
1423 fprintf(stderr, "ERROR: memory overflow.\n");
1424 exit(1);
1426 if (!prefix)
1427 snprintf(names[i], len, "%d", i);
1428 else
1429 snprintf(names[i], len, "%s%d", prefix, i);
1432 return names;
1435 void util_free_names(int n, char **names)
1437 int i;
1438 for (i = 0; i < n; ++i)
1439 free(names[i]);
1440 free(names);
1443 void Polyhedron_pprint(FILE *out, Polyhedron *P, int dim, int nparam,
1444 char **iter_names, char **param_names)
1446 int i, j;
1447 Value tmp;
1449 assert(dim + nparam == P->Dimension);
1451 value_init(tmp);
1453 fprintf(out, "{ ");
1454 if (nparam) {
1455 print_varlist(out, nparam, param_names);
1456 fprintf(out, " -> ");
1458 print_varlist(out, dim, iter_names);
1459 fprintf(out, " : ");
1461 if (emptyQ2(P))
1462 fprintf(out, "FALSE");
1463 else for (i = 0; i < P->NbConstraints; ++i) {
1464 int first = 1;
1465 int v = First_Non_Zero(P->Constraint[i]+1, P->Dimension);
1466 if (v == -1 && value_pos_p(P->Constraint[i][0]))
1467 continue;
1468 if (i)
1469 fprintf(out, " && ");
1470 if (v == -1 && value_notzero_p(P->Constraint[i][1+P->Dimension]))
1471 fprintf(out, "FALSE");
1472 else if (value_pos_p(P->Constraint[i][v+1])) {
1473 print_term(out, P->Constraint[i][v+1], v, dim, nparam,
1474 iter_names, param_names, NULL);
1475 if (value_zero_p(P->Constraint[i][0]))
1476 fprintf(out, " = ");
1477 else
1478 fprintf(out, " >= ");
1479 for (j = v+1; j <= dim+nparam; ++j) {
1480 value_oppose(tmp, P->Constraint[i][1+j]);
1481 print_term(out, tmp, j, dim, nparam,
1482 iter_names, param_names, &first);
1484 } else {
1485 value_oppose(tmp, P->Constraint[i][1+v]);
1486 print_term(out, tmp, v, dim, nparam,
1487 iter_names, param_names, NULL);
1488 fprintf(out, " <= ");
1489 for (j = v+1; j <= dim+nparam; ++j)
1490 print_term(out, P->Constraint[i][1+j], j, dim, nparam,
1491 iter_names, param_names, &first);
1495 fprintf(out, " }\n");
1497 value_clear(tmp);
1500 /* Construct a cone over P with P placed at x_d = 1, with
1501 * x_d the coordinate of an extra dimension
1503 * It's probably a mistake to depend so much on the internal
1504 * representation. We should probably simply compute the
1505 * vertices/facets first.
1507 Polyhedron *Cone_over_Polyhedron(Polyhedron *P)
1509 unsigned NbConstraints = 0;
1510 unsigned NbRays = 0;
1511 Polyhedron *C;
1512 int i;
1514 if (POL_HAS(P, POL_INEQUALITIES))
1515 NbConstraints = P->NbConstraints + 1;
1516 if (POL_HAS(P, POL_POINTS))
1517 NbRays = P->NbRays + 1;
1519 C = Polyhedron_Alloc(P->Dimension+1, NbConstraints, NbRays);
1520 if (POL_HAS(P, POL_INEQUALITIES)) {
1521 C->NbEq = P->NbEq;
1522 for (i = 0; i < P->NbConstraints; ++i)
1523 Vector_Copy(P->Constraint[i], C->Constraint[i], P->Dimension+2);
1524 /* n >= 0 */
1525 value_set_si(C->Constraint[P->NbConstraints][0], 1);
1526 value_set_si(C->Constraint[P->NbConstraints][1+P->Dimension], 1);
1528 if (POL_HAS(P, POL_POINTS)) {
1529 C->NbBid = P->NbBid;
1530 for (i = 0; i < P->NbRays; ++i)
1531 Vector_Copy(P->Ray[i], C->Ray[i], P->Dimension+2);
1532 /* vertex 0 */
1533 value_set_si(C->Ray[P->NbRays][0], 1);
1534 value_set_si(C->Ray[P->NbRays][1+C->Dimension], 1);
1536 POL_SET(C, POL_VALID);
1537 if (POL_HAS(P, POL_INEQUALITIES))
1538 POL_SET(C, POL_INEQUALITIES);
1539 if (POL_HAS(P, POL_POINTS))
1540 POL_SET(C, POL_POINTS);
1541 if (POL_HAS(P, POL_VERTICES))
1542 POL_SET(C, POL_VERTICES);
1543 return C;
1546 /* Returns a (dim+nparam+1)x((dim-n)+nparam+1) matrix
1547 * mapping the transformed subspace back to the original space.
1548 * n is the number of equalities involving the variables
1549 * (i.e., not purely the parameters).
1550 * The remaining n coordinates in the transformed space would
1551 * have constant (parametric) values and are therefore not
1552 * included in the variables of the new space.
1554 Matrix *compress_variables(Matrix *Equalities, unsigned nparam)
1556 unsigned dim = (Equalities->NbColumns-2) - nparam;
1557 Matrix *M, *H, *Q, *U, *C, *ratH, *invH, *Ul, *T1, *T2, *T;
1558 Value mone;
1559 int n, i, j;
1560 int ok;
1562 for (n = 0; n < Equalities->NbRows; ++n)
1563 if (First_Non_Zero(Equalities->p[n]+1, dim) == -1)
1564 break;
1565 if (n == 0)
1566 return Identity(dim+nparam+1);
1567 value_init(mone);
1568 value_set_si(mone, -1);
1569 M = Matrix_Alloc(n, dim);
1570 C = Matrix_Alloc(n+1, nparam+1);
1571 for (i = 0; i < n; ++i) {
1572 Vector_Copy(Equalities->p[i]+1, M->p[i], dim);
1573 Vector_Scale(Equalities->p[i]+1+dim, C->p[i], mone, nparam+1);
1575 value_set_si(C->p[n][nparam], 1);
1576 left_hermite(M, &H, &Q, &U);
1577 Matrix_Free(M);
1578 Matrix_Free(Q);
1579 value_clear(mone);
1581 ratH = Matrix_Alloc(n+1, n+1);
1582 invH = Matrix_Alloc(n+1, n+1);
1583 for (i = 0; i < n; ++i)
1584 Vector_Copy(H->p[i], ratH->p[i], n);
1585 value_set_si(ratH->p[n][n], 1);
1586 ok = Matrix_Inverse(ratH, invH);
1587 assert(ok);
1588 Matrix_Free(H);
1589 Matrix_Free(ratH);
1590 T1 = Matrix_Alloc(n+1, nparam+1);
1591 Matrix_Product(invH, C, T1);
1592 Matrix_Free(C);
1593 Matrix_Free(invH);
1594 if (value_notone_p(T1->p[n][nparam])) {
1595 for (i = 0; i < n; ++i) {
1596 if (!mpz_divisible_p(T1->p[i][nparam], T1->p[n][nparam])) {
1597 Matrix_Free(T1);
1598 Matrix_Free(U);
1599 return NULL;
1601 /* compress_params should have taken care of this */
1602 for (j = 0; j < nparam; ++j)
1603 assert(mpz_divisible_p(T1->p[i][j], T1->p[n][nparam]));
1604 Vector_AntiScale(T1->p[i], T1->p[i], T1->p[n][nparam], nparam+1);
1606 value_set_si(T1->p[n][nparam], 1);
1608 Ul = Matrix_Alloc(dim+1, n+1);
1609 for (i = 0; i < dim; ++i)
1610 Vector_Copy(U->p[i], Ul->p[i], n);
1611 value_set_si(Ul->p[dim][n], 1);
1612 T2 = Matrix_Alloc(dim+1, nparam+1);
1613 Matrix_Product(Ul, T1, T2);
1614 Matrix_Free(Ul);
1615 Matrix_Free(T1);
1617 T = Matrix_Alloc(dim+nparam+1, (dim-n)+nparam+1);
1618 for (i = 0; i < dim; ++i) {
1619 Vector_Copy(U->p[i]+n, T->p[i], dim-n);
1620 Vector_Copy(T2->p[i], T->p[i]+dim-n, nparam+1);
1622 for (i = 0; i < nparam+1; ++i)
1623 value_set_si(T->p[dim+i][(dim-n)+i], 1);
1624 assert(value_one_p(T2->p[dim][nparam]));
1625 Matrix_Free(U);
1626 Matrix_Free(T2);
1628 return T;
1631 /* Computes the left inverse of an affine embedding M and, if Eq is not NULL,
1632 * the equalities that define the affine subspace onto which M maps
1633 * its argument.
1635 Matrix *left_inverse(Matrix *M, Matrix **Eq)
1637 int i, ok;
1638 Matrix *L, *H, *Q, *U, *ratH, *invH, *Ut, *inv;
1639 Vector *t;
1641 if (M->NbColumns == 1) {
1642 inv = Matrix_Alloc(1, M->NbRows);
1643 value_set_si(inv->p[0][M->NbRows-1], 1);
1644 if (Eq) {
1645 *Eq = Matrix_Alloc(M->NbRows-1, 1+(M->NbRows-1)+1);
1646 for (i = 0; i < M->NbRows-1; ++i) {
1647 value_oppose((*Eq)->p[i][1+i], M->p[M->NbRows-1][0]);
1648 value_assign((*Eq)->p[i][1+(M->NbRows-1)], M->p[i][0]);
1651 return inv;
1653 if (Eq)
1654 *Eq = NULL;
1655 L = Matrix_Alloc(M->NbRows-1, M->NbColumns-1);
1656 for (i = 0; i < L->NbRows; ++i)
1657 Vector_Copy(M->p[i], L->p[i], L->NbColumns);
1658 right_hermite(L, &H, &U, &Q);
1659 Matrix_Free(L);
1660 Matrix_Free(Q);
1661 t = Vector_Alloc(U->NbColumns);
1662 for (i = 0; i < U->NbColumns; ++i)
1663 value_oppose(t->p[i], M->p[i][M->NbColumns-1]);
1664 if (Eq) {
1665 *Eq = Matrix_Alloc(H->NbRows - H->NbColumns, 2 + U->NbColumns);
1666 for (i = 0; i < H->NbRows - H->NbColumns; ++i) {
1667 Vector_Copy(U->p[H->NbColumns+i], (*Eq)->p[i]+1, U->NbColumns);
1668 Inner_Product(U->p[H->NbColumns+i], t->p, U->NbColumns,
1669 (*Eq)->p[i]+1+U->NbColumns);
1672 ratH = Matrix_Alloc(H->NbColumns+1, H->NbColumns+1);
1673 invH = Matrix_Alloc(H->NbColumns+1, H->NbColumns+1);
1674 for (i = 0; i < H->NbColumns; ++i)
1675 Vector_Copy(H->p[i], ratH->p[i], H->NbColumns);
1676 value_set_si(ratH->p[ratH->NbRows-1][ratH->NbColumns-1], 1);
1677 Matrix_Free(H);
1678 ok = Matrix_Inverse(ratH, invH);
1679 assert(ok);
1680 Matrix_Free(ratH);
1681 Ut = Matrix_Alloc(invH->NbRows, U->NbColumns+1);
1682 for (i = 0; i < Ut->NbRows-1; ++i) {
1683 Vector_Copy(U->p[i], Ut->p[i], U->NbColumns);
1684 Inner_Product(U->p[i], t->p, U->NbColumns, &Ut->p[i][Ut->NbColumns-1]);
1686 Matrix_Free(U);
1687 Vector_Free(t);
1688 value_set_si(Ut->p[Ut->NbRows-1][Ut->NbColumns-1], 1);
1689 inv = Matrix_Alloc(invH->NbRows, Ut->NbColumns);
1690 Matrix_Product(invH, Ut, inv);
1691 Matrix_Free(Ut);
1692 Matrix_Free(invH);
1693 return inv;
1696 /* Check whether all rays are revlex positive in the parameters
1698 int Polyhedron_has_revlex_positive_rays(Polyhedron *P, unsigned nparam)
1700 int r;
1701 for (r = 0; r < P->NbRays; ++r) {
1702 int i;
1703 if (value_notzero_p(P->Ray[r][P->Dimension+1]))
1704 continue;
1705 for (i = P->Dimension-1; i >= P->Dimension-nparam; --i) {
1706 if (value_neg_p(P->Ray[r][i+1]))
1707 return 0;
1708 if (value_pos_p(P->Ray[r][i+1]))
1709 break;
1711 /* A ray independent of the parameters */
1712 if (i < P->Dimension-nparam)
1713 return 0;
1715 return 1;
1718 static Polyhedron *Recession_Cone(Polyhedron *P, unsigned nparam, unsigned MaxRays)
1720 int i;
1721 unsigned nvar = P->Dimension - nparam;
1722 Matrix *M = Matrix_Alloc(P->NbConstraints, 1 + nvar + 1);
1723 Polyhedron *R;
1724 for (i = 0; i < P->NbConstraints; ++i)
1725 Vector_Copy(P->Constraint[i], M->p[i], 1+nvar);
1726 R = Constraints2Polyhedron(M, MaxRays);
1727 Matrix_Free(M);
1728 return R;
1731 int Polyhedron_is_unbounded(Polyhedron *P, unsigned nparam, unsigned MaxRays)
1733 int i;
1734 int is_unbounded;
1735 Polyhedron *R = Recession_Cone(P, nparam, MaxRays);
1736 POL_ENSURE_VERTICES(R);
1737 if (R->NbBid == 0)
1738 for (i = 0; i < R->NbRays; ++i)
1739 if (value_zero_p(R->Ray[i][1+R->Dimension]))
1740 break;
1741 is_unbounded = R->NbBid > 0 || i < R->NbRays;
1742 Polyhedron_Free(R);
1743 return is_unbounded;
1746 static void SwapColumns(Value **V, int n, int i, int j)
1748 int r;
1750 for (r = 0; r < n; ++r)
1751 value_swap(V[r][i], V[r][j]);
1754 void Polyhedron_ExchangeColumns(Polyhedron *P, int Column1, int Column2)
1756 SwapColumns(P->Constraint, P->NbConstraints, Column1, Column2);
1757 SwapColumns(P->Ray, P->NbRays, Column1, Column2);
1758 if (P->NbEq) {
1759 Matrix M;
1760 Polyhedron_Matrix_View(P, &M, P->NbConstraints);
1761 Gauss(&M, P->NbEq, P->Dimension+1);
1765 /* perform transposition inline; assumes M is a square matrix */
1766 void Matrix_Transposition(Matrix *M)
1768 int i, j;
1770 assert(M->NbRows == M->NbColumns);
1771 for (i = 0; i < M->NbRows; ++i)
1772 for (j = i+1; j < M->NbColumns; ++j)
1773 value_swap(M->p[i][j], M->p[j][i]);
1776 /* Matrix "view" of first rows rows */
1777 void Polyhedron_Matrix_View(Polyhedron *P, Matrix *M, unsigned rows)
1779 M->NbRows = rows;
1780 M->NbColumns = P->Dimension+2;
1781 M->p_Init = P->p_Init;
1782 M->p = P->Constraint;
1785 int Last_Non_Zero(Value *p, unsigned len)
1787 int i;
1789 for (i = len - 1; i >= 0; --i)
1790 if (value_notzero_p(p[i]))
1791 return i;
1793 return -1;