barvinok_count: stop counting if one of the factors is empty
[barvinok.git] / barvinok.cc
blobe5e6dcd487df9a9052d27584d273b362d9e9dcab
1 #include <assert.h>
2 #include <iostream>
3 #include <vector>
4 #include <deque>
5 #include <string>
6 #include <sstream>
7 #include <gmp.h>
8 #include <NTL/mat_ZZ.h>
9 #include <NTL/LLL.h>
10 #include <barvinok/util.h>
11 extern "C" {
12 #include <polylib/polylibgmp.h>
13 #include <barvinok/evalue.h>
14 #include "piputil.h"
16 #include "config.h"
17 #include <barvinok/barvinok.h>
18 #include <barvinok/genfun.h>
19 #include "conversion.h"
20 #include "decomposer.h"
21 #include "lattice_point.h"
22 #include "reduce_domain.h"
24 #ifdef NTL_STD_CXX
25 using namespace NTL;
26 #endif
27 using std::cerr;
28 using std::cout;
29 using std::endl;
30 using std::vector;
31 using std::deque;
32 using std::string;
33 using std::ostringstream;
35 #define ALLOC(t,p) p = (t*)malloc(sizeof(*p))
37 static void rays(mat_ZZ& r, Polyhedron *C)
39 unsigned dim = C->NbRays - 1; /* don't count zero vertex */
40 assert(C->NbRays - 1 == C->Dimension);
41 r.SetDims(dim, dim);
42 ZZ tmp;
44 int i, c;
45 for (i = 0, c = 0; i < dim; ++i)
46 if (value_zero_p(C->Ray[i][dim+1])) {
47 for (int j = 0; j < dim; ++j) {
48 value2zz(C->Ray[i][j+1], tmp);
49 r[j][c] = tmp;
51 ++c;
55 class dpoly {
56 public:
57 vec_ZZ coeff;
58 dpoly(int d, ZZ& degree, int offset = 0) {
59 coeff.SetLength(d+1);
61 int min = d + offset;
62 if (degree >= 0 && degree < ZZ(INIT_VAL, min))
63 min = to_int(degree);
65 ZZ c = ZZ(INIT_VAL, 1);
66 if (!offset)
67 coeff[0] = c;
68 for (int i = 1; i <= min; ++i) {
69 c *= (degree -i + 1);
70 c /= i;
71 coeff[i-offset] = c;
74 void operator *= (dpoly& f) {
75 assert(coeff.length() == f.coeff.length());
76 vec_ZZ old = coeff;
77 coeff = f.coeff[0] * coeff;
78 for (int i = 1; i < coeff.length(); ++i)
79 for (int j = 0; i+j < coeff.length(); ++j)
80 coeff[i+j] += f.coeff[i] * old[j];
82 void div(dpoly& d, mpq_t count, ZZ& sign) {
83 int len = coeff.length();
84 Value tmp;
85 value_init(tmp);
86 mpq_t* c = new mpq_t[coeff.length()];
87 mpq_t qtmp;
88 mpq_init(qtmp);
89 for (int i = 0; i < len; ++i) {
90 mpq_init(c[i]);
91 zz2value(coeff[i], tmp);
92 mpq_set_z(c[i], tmp);
94 for (int j = 1; j <= i; ++j) {
95 zz2value(d.coeff[j], tmp);
96 mpq_set_z(qtmp, tmp);
97 mpq_mul(qtmp, qtmp, c[i-j]);
98 mpq_sub(c[i], c[i], qtmp);
101 zz2value(d.coeff[0], tmp);
102 mpq_set_z(qtmp, tmp);
103 mpq_div(c[i], c[i], qtmp);
105 if (sign == -1)
106 mpq_sub(count, count, c[len-1]);
107 else
108 mpq_add(count, count, c[len-1]);
110 value_clear(tmp);
111 mpq_clear(qtmp);
112 for (int i = 0; i < len; ++i)
113 mpq_clear(c[i]);
114 delete [] c;
118 class dpoly_n {
119 public:
120 Matrix *coeff;
121 ~dpoly_n() {
122 Matrix_Free(coeff);
124 dpoly_n(int d, ZZ& degree_0, ZZ& degree_1, int offset = 0) {
125 Value d0, d1;
126 value_init(d0);
127 value_init(d1);
128 zz2value(degree_0, d0);
129 zz2value(degree_1, d1);
130 coeff = Matrix_Alloc(d+1, d+1+1);
131 value_set_si(coeff->p[0][0], 1);
132 value_set_si(coeff->p[0][d+1], 1);
133 for (int i = 1; i <= d; ++i) {
134 value_multiply(coeff->p[i][0], coeff->p[i-1][0], d0);
135 Vector_Combine(coeff->p[i-1], coeff->p[i-1]+1, coeff->p[i]+1,
136 d1, d0, i);
137 value_set_si(coeff->p[i][d+1], i);
138 value_multiply(coeff->p[i][d+1], coeff->p[i][d+1], coeff->p[i-1][d+1]);
139 value_decrement(d0, d0);
141 value_clear(d0);
142 value_clear(d1);
144 void div(dpoly& d, Vector *count, ZZ& sign) {
145 int len = coeff->NbRows;
146 Matrix * c = Matrix_Alloc(coeff->NbRows, coeff->NbColumns);
147 Value tmp;
148 value_init(tmp);
149 for (int i = 0; i < len; ++i) {
150 Vector_Copy(coeff->p[i], c->p[i], len+1);
151 for (int j = 1; j <= i; ++j) {
152 zz2value(d.coeff[j], tmp);
153 value_multiply(tmp, tmp, c->p[i][len]);
154 value_oppose(tmp, tmp);
155 Vector_Combine(c->p[i], c->p[i-j], c->p[i],
156 c->p[i-j][len], tmp, len);
157 value_multiply(c->p[i][len], c->p[i][len], c->p[i-j][len]);
159 zz2value(d.coeff[0], tmp);
160 value_multiply(c->p[i][len], c->p[i][len], tmp);
162 if (sign == -1) {
163 value_set_si(tmp, -1);
164 Vector_Scale(c->p[len-1], count->p, tmp, len);
165 value_assign(count->p[len], c->p[len-1][len]);
166 } else
167 Vector_Copy(c->p[len-1], count->p, len+1);
168 Vector_Normalize(count->p, len+1);
169 value_clear(tmp);
170 Matrix_Free(c);
174 struct dpoly_r_term {
175 int *powers;
176 ZZ coeff;
179 /* len: number of elements in c
180 * each element in c is the coefficient of a power of t
181 * in the MacLaurin expansion
183 struct dpoly_r {
184 vector< dpoly_r_term * > *c;
185 int len;
186 int dim;
187 ZZ denom;
189 void add_term(int i, int * powers, ZZ& coeff) {
190 if (coeff == 0)
191 return;
192 for (int k = 0; k < c[i].size(); ++k) {
193 if (memcmp(c[i][k]->powers, powers, dim * sizeof(int)) == 0) {
194 c[i][k]->coeff += coeff;
195 return;
198 dpoly_r_term *t = new dpoly_r_term;
199 t->powers = new int[dim];
200 memcpy(t->powers, powers, dim * sizeof(int));
201 t->coeff = coeff;
202 c[i].push_back(t);
204 dpoly_r(int len, int dim) {
205 denom = 1;
206 this->len = len;
207 this->dim = dim;
208 c = new vector< dpoly_r_term * > [len];
210 dpoly_r(dpoly& num, int dim) {
211 denom = 1;
212 len = num.coeff.length();
213 c = new vector< dpoly_r_term * > [len];
214 this->dim = dim;
215 int powers[dim];
216 memset(powers, 0, dim * sizeof(int));
218 for (int i = 0; i < len; ++i) {
219 ZZ coeff = num.coeff[i];
220 add_term(i, powers, coeff);
223 dpoly_r(dpoly& num, dpoly& den, int pos, int dim) {
224 denom = 1;
225 len = num.coeff.length();
226 c = new vector< dpoly_r_term * > [len];
227 this->dim = dim;
228 int powers[dim];
230 for (int i = 0; i < len; ++i) {
231 ZZ coeff = num.coeff[i];
232 memset(powers, 0, dim * sizeof(int));
233 powers[pos] = 1;
235 add_term(i, powers, coeff);
237 for (int j = 1; j <= i; ++j) {
238 for (int k = 0; k < c[i-j].size(); ++k) {
239 memcpy(powers, c[i-j][k]->powers, dim*sizeof(int));
240 powers[pos]++;
241 coeff = -den.coeff[j-1] * c[i-j][k]->coeff;
242 add_term(i, powers, coeff);
246 //dump();
248 dpoly_r(dpoly_r* num, dpoly& den, int pos, int dim) {
249 denom = num->denom;
250 len = num->len;
251 c = new vector< dpoly_r_term * > [len];
252 this->dim = dim;
253 int powers[dim];
254 ZZ coeff;
256 for (int i = 0 ; i < len; ++i) {
257 for (int k = 0; k < num->c[i].size(); ++k) {
258 memcpy(powers, num->c[i][k]->powers, dim*sizeof(int));
259 powers[pos]++;
260 add_term(i, powers, num->c[i][k]->coeff);
263 for (int j = 1; j <= i; ++j) {
264 for (int k = 0; k < c[i-j].size(); ++k) {
265 memcpy(powers, c[i-j][k]->powers, dim*sizeof(int));
266 powers[pos]++;
267 coeff = -den.coeff[j-1] * c[i-j][k]->coeff;
268 add_term(i, powers, coeff);
273 ~dpoly_r() {
274 for (int i = 0 ; i < len; ++i)
275 for (int k = 0; k < c[i].size(); ++k) {
276 delete [] c[i][k]->powers;
277 delete c[i][k];
279 delete [] c;
281 dpoly_r *div(dpoly& d) {
282 dpoly_r *rc = new dpoly_r(len, dim);
283 rc->denom = power(d.coeff[0], len);
284 ZZ inv_d = rc->denom / d.coeff[0];
285 ZZ coeff;
287 for (int i = 0; i < len; ++i) {
288 for (int k = 0; k < c[i].size(); ++k) {
289 coeff = c[i][k]->coeff * inv_d;
290 rc->add_term(i, c[i][k]->powers, coeff);
293 for (int j = 1; j <= i; ++j) {
294 for (int k = 0; k < rc->c[i-j].size(); ++k) {
295 coeff = - d.coeff[j] * rc->c[i-j][k]->coeff / d.coeff[0];
296 rc->add_term(i, rc->c[i-j][k]->powers, coeff);
300 return rc;
302 void dump(void) {
303 for (int i = 0; i < len; ++i) {
304 cerr << endl;
305 cerr << i << endl;
306 cerr << c[i].size() << endl;
307 for (int j = 0; j < c[i].size(); ++j) {
308 for (int k = 0; k < dim; ++k) {
309 cerr << c[i][j]->powers[k] << " ";
311 cerr << ": " << c[i][j]->coeff << "/" << denom << endl;
313 cerr << endl;
318 const int MAX_TRY=10;
320 * Searches for a vector that is not orthogonal to any
321 * of the rays in rays.
323 static void nonorthog(mat_ZZ& rays, vec_ZZ& lambda)
325 int dim = rays.NumCols();
326 bool found = false;
327 lambda.SetLength(dim);
328 if (dim == 0)
329 return;
331 for (int i = 2; !found && i <= 50*dim; i+=4) {
332 for (int j = 0; j < MAX_TRY; ++j) {
333 for (int k = 0; k < dim; ++k) {
334 int r = random_int(i)+2;
335 int v = (2*(r%2)-1) * (r >> 1);
336 lambda[k] = v;
338 int k = 0;
339 for (; k < rays.NumRows(); ++k)
340 if (lambda * rays[k] == 0)
341 break;
342 if (k == rays.NumRows()) {
343 found = true;
344 break;
348 assert(found);
351 static void randomvector(Polyhedron *P, vec_ZZ& lambda, int nvar)
353 Value tmp;
354 int max = 10 * 16;
355 unsigned int dim = P->Dimension;
356 value_init(tmp);
358 for (int i = 0; i < P->NbRays; ++i) {
359 for (int j = 1; j <= dim; ++j) {
360 value_absolute(tmp, P->Ray[i][j]);
361 int t = VALUE_TO_LONG(tmp) * 16;
362 if (t > max)
363 max = t;
366 for (int i = 0; i < P->NbConstraints; ++i) {
367 for (int j = 1; j <= dim; ++j) {
368 value_absolute(tmp, P->Constraint[i][j]);
369 int t = VALUE_TO_LONG(tmp) * 16;
370 if (t > max)
371 max = t;
374 value_clear(tmp);
376 lambda.SetLength(nvar);
377 for (int k = 0; k < nvar; ++k) {
378 int r = random_int(max*dim)+2;
379 int v = (2*(r%2)-1) * (max/2*dim + (r >> 1));
380 lambda[k] = v;
384 static void add_rays(mat_ZZ& rays, Polyhedron *i, int *r, int nvar = -1,
385 bool all = false)
387 unsigned dim = i->Dimension;
388 if (nvar == -1)
389 nvar = dim;
390 for (int k = 0; k < i->NbRays; ++k) {
391 if (!value_zero_p(i->Ray[k][dim+1]))
392 continue;
393 if (!all && nvar != dim && First_Non_Zero(i->Ray[k]+1, nvar) == -1)
394 continue;
395 values2zz(i->Ray[k]+1, rays[(*r)++], nvar);
399 static void mask_r(Matrix *f, int nr, Vector *lcm, int p, Vector *val, evalue *ev)
401 unsigned nparam = lcm->Size;
403 if (p == nparam) {
404 Vector * prod = Vector_Alloc(f->NbRows);
405 Matrix_Vector_Product(f, val->p, prod->p);
406 int isint = 1;
407 for (int i = 0; i < nr; ++i) {
408 value_modulus(prod->p[i], prod->p[i], f->p[i][nparam+1]);
409 isint &= value_zero_p(prod->p[i]);
411 value_set_si(ev->d, 1);
412 value_init(ev->x.n);
413 value_set_si(ev->x.n, isint);
414 Vector_Free(prod);
415 return;
418 Value tmp;
419 value_init(tmp);
420 if (value_one_p(lcm->p[p]))
421 mask_r(f, nr, lcm, p+1, val, ev);
422 else {
423 value_assign(tmp, lcm->p[p]);
424 value_set_si(ev->d, 0);
425 ev->x.p = new_enode(periodic, VALUE_TO_INT(tmp), p+1);
426 do {
427 value_decrement(tmp, tmp);
428 value_assign(val->p[p], tmp);
429 mask_r(f, nr, lcm, p+1, val, &ev->x.p->arr[VALUE_TO_INT(tmp)]);
430 } while (value_pos_p(tmp));
432 value_clear(tmp);
435 #ifdef USE_MODULO
436 static void mask(Matrix *f, evalue *factor)
438 int nr = f->NbRows, nc = f->NbColumns;
439 int n;
440 bool found = false;
441 for (n = 0; n < nr && value_notzero_p(f->p[n][nc-1]); ++n)
442 if (value_notone_p(f->p[n][nc-1]) &&
443 value_notmone_p(f->p[n][nc-1]))
444 found = true;
445 if (!found)
446 return;
448 evalue EP;
449 nr = n;
451 Value m;
452 value_init(m);
454 evalue EV;
455 value_init(EV.d);
456 value_init(EV.x.n);
457 value_set_si(EV.x.n, 1);
459 for (n = 0; n < nr; ++n) {
460 value_assign(m, f->p[n][nc-1]);
461 if (value_one_p(m) || value_mone_p(m))
462 continue;
464 int j = normal_mod(f->p[n], nc-1, &m);
465 if (j == nc-1) {
466 free_evalue_refs(factor);
467 value_init(factor->d);
468 evalue_set_si(factor, 0, 1);
469 break;
471 vec_ZZ row;
472 values2zz(f->p[n], row, nc-1);
473 ZZ g;
474 value2zz(m, g);
475 if (j < (nc-1)-1 && row[j] > g/2) {
476 for (int k = j; k < (nc-1); ++k)
477 if (row[k] != 0)
478 row[k] = g - row[k];
481 value_init(EP.d);
482 value_set_si(EP.d, 0);
483 EP.x.p = new_enode(relation, 2, 0);
484 value_clear(EP.x.p->arr[1].d);
485 EP.x.p->arr[1] = *factor;
486 evalue *ev = &EP.x.p->arr[0];
487 value_set_si(ev->d, 0);
488 ev->x.p = new_enode(fractional, 3, -1);
489 evalue_set_si(&ev->x.p->arr[1], 0, 1);
490 evalue_set_si(&ev->x.p->arr[2], 1, 1);
491 evalue *E = multi_monom(row);
492 value_assign(EV.d, m);
493 emul(&EV, E);
494 value_clear(ev->x.p->arr[0].d);
495 ev->x.p->arr[0] = *E;
496 delete E;
497 *factor = EP;
500 value_clear(m);
501 free_evalue_refs(&EV);
503 #else
507 static void mask(Matrix *f, evalue *factor)
509 int nr = f->NbRows, nc = f->NbColumns;
510 int n;
511 bool found = false;
512 for (n = 0; n < nr && value_notzero_p(f->p[n][nc-1]); ++n)
513 if (value_notone_p(f->p[n][nc-1]) &&
514 value_notmone_p(f->p[n][nc-1]))
515 found = true;
516 if (!found)
517 return;
519 Value tmp;
520 value_init(tmp);
521 nr = n;
522 unsigned np = nc - 2;
523 Vector *lcm = Vector_Alloc(np);
524 Vector *val = Vector_Alloc(nc);
525 Vector_Set(val->p, 0, nc);
526 value_set_si(val->p[np], 1);
527 Vector_Set(lcm->p, 1, np);
528 for (n = 0; n < nr; ++n) {
529 if (value_one_p(f->p[n][nc-1]) ||
530 value_mone_p(f->p[n][nc-1]))
531 continue;
532 for (int j = 0; j < np; ++j)
533 if (value_notzero_p(f->p[n][j])) {
534 Gcd(f->p[n][j], f->p[n][nc-1], &tmp);
535 value_division(tmp, f->p[n][nc-1], tmp);
536 value_lcm(tmp, lcm->p[j], &lcm->p[j]);
539 evalue EP;
540 value_init(EP.d);
541 mask_r(f, nr, lcm, 0, val, &EP);
542 value_clear(tmp);
543 Vector_Free(val);
544 Vector_Free(lcm);
545 emul(&EP,factor);
546 free_evalue_refs(&EP);
548 #endif
550 /* This structure encodes the power of the term in a rational generating function.
552 * Either E == NULL or constant = 0
553 * If E != NULL, then the power is E
554 * If E == NULL, then the power is coeff * param[pos] + constant
556 struct term_info {
557 evalue *E;
558 ZZ constant;
559 ZZ coeff;
560 int pos;
563 /* Returns the power of (t+1) in the term of a rational generating function,
564 * i.e., the scalar product of the actual lattice point and lambda.
565 * The lattice point is the unique lattice point in the fundamental parallelepiped
566 * of the unimodual cone i shifted to the parametric vertex V.
568 * PD is the parameter domain, which, if != NULL, may be used to simply the
569 * resulting expression.
571 * The result is returned in term.
573 void lattice_point(
574 Param_Vertices* V, Polyhedron *i, vec_ZZ& lambda, term_info* term,
575 Polyhedron *PD)
577 unsigned nparam = V->Vertex->NbColumns - 2;
578 unsigned dim = i->Dimension;
579 mat_ZZ vertex;
580 vertex.SetDims(V->Vertex->NbRows, nparam+1);
581 Value lcm, tmp;
582 value_init(lcm);
583 value_init(tmp);
584 value_set_si(lcm, 1);
585 for (int j = 0; j < V->Vertex->NbRows; ++j) {
586 value_lcm(lcm, V->Vertex->p[j][nparam+1], &lcm);
588 if (value_notone_p(lcm)) {
589 Matrix * mv = Matrix_Alloc(dim, nparam+1);
590 for (int j = 0 ; j < dim; ++j) {
591 value_division(tmp, lcm, V->Vertex->p[j][nparam+1]);
592 Vector_Scale(V->Vertex->p[j], mv->p[j], tmp, nparam+1);
595 term->E = lattice_point(i, lambda, mv, lcm, PD);
596 term->constant = 0;
598 Matrix_Free(mv);
599 value_clear(lcm);
600 value_clear(tmp);
601 return;
603 for (int i = 0; i < V->Vertex->NbRows; ++i) {
604 assert(value_one_p(V->Vertex->p[i][nparam+1])); // for now
605 values2zz(V->Vertex->p[i], vertex[i], nparam+1);
608 vec_ZZ num;
609 num = lambda * vertex;
611 int p = -1;
612 int nn = 0;
613 for (int j = 0; j < nparam; ++j)
614 if (num[j] != 0) {
615 ++nn;
616 p = j;
618 if (nn >= 2) {
619 term->E = multi_monom(num);
620 term->constant = 0;
621 } else {
622 term->E = NULL;
623 term->constant = num[nparam];
624 term->pos = p;
625 if (p != -1)
626 term->coeff = num[p];
629 value_clear(lcm);
630 value_clear(tmp);
633 static void normalize(ZZ& sign, ZZ& num, vec_ZZ& den)
635 unsigned dim = den.length();
637 int change = 0;
639 for (int j = 0; j < den.length(); ++j) {
640 if (den[j] > 0)
641 change ^= 1;
642 else {
643 den[j] = abs(den[j]);
644 num += den[j];
647 if (change)
648 sign = -sign;
651 /* input:
652 * f: the powers in the denominator for the remaining vars
653 * each row refers to a factor
654 * den_s: for each factor, the power of (s+1)
655 * sign
656 * num_s: powers in the numerator corresponding to the summed vars
657 * num_p: powers in the numerator corresponding to the remaining vars
658 * number of rays in cone: "dim" = "k"
659 * length of each ray: "dim" = "d"
660 * for now, it is assumed: k == d
661 * output:
662 * den_p: for each factor
663 * 0: independent of remaining vars
664 * 1: power corresponds to corresponding row in f
666 * all inputs are subject to change
668 static void normalize(ZZ& sign,
669 ZZ& num_s, vec_ZZ& num_p, vec_ZZ& den_s, vec_ZZ& den_p,
670 mat_ZZ& f)
672 unsigned dim = f.NumRows();
673 unsigned nparam = num_p.length();
674 unsigned nvar = dim - nparam;
676 int change = 0;
678 for (int j = 0; j < den_s.length(); ++j) {
679 if (den_s[j] == 0) {
680 den_p[j] = 1;
681 continue;
683 int k;
684 for (k = 0; k < nparam; ++k)
685 if (f[j][k] != 0)
686 break;
687 if (k < nparam) {
688 den_p[j] = 1;
689 if (den_s[j] > 0) {
690 f[j] = -f[j];
691 num_p += f[j];
693 } else
694 den_p[j] = 0;
695 if (den_s[j] > 0)
696 change ^= 1;
697 else {
698 den_s[j] = abs(den_s[j]);
699 num_s += den_s[j];
703 if (change)
704 sign = -sign;
707 struct counter : public polar_decomposer {
708 vec_ZZ lambda;
709 mat_ZZ rays;
710 vec_ZZ vertex;
711 vec_ZZ den;
712 ZZ sign;
713 ZZ num;
714 int j;
715 Polyhedron *P;
716 unsigned dim;
717 mpq_t count;
719 counter(Polyhedron *P) {
720 this->P = P;
721 dim = P->Dimension;
722 rays.SetDims(dim, dim);
723 den.SetLength(dim);
724 mpq_init(count);
727 void start(unsigned MaxRays);
729 ~counter() {
730 mpq_clear(count);
733 virtual void handle_polar(Polyhedron *P, int sign);
736 struct OrthogonalException {} Orthogonal;
738 void counter::handle_polar(Polyhedron *C, int s)
740 int r = 0;
741 assert(C->NbRays-1 == dim);
742 add_rays(rays, C, &r);
743 for (int k = 0; k < dim; ++k) {
744 if (lambda * rays[k] == 0)
745 throw Orthogonal;
748 sign = s;
750 lattice_point(P->Ray[j]+1, C, vertex);
751 num = vertex * lambda;
752 den = rays * lambda;
753 normalize(sign, num, den);
755 dpoly d(dim, num);
756 dpoly n(dim, den[0], 1);
757 for (int k = 1; k < dim; ++k) {
758 dpoly fact(dim, den[k], 1);
759 n *= fact;
761 d.div(n, count, sign);
764 void counter::start(unsigned MaxRays)
766 for (;;) {
767 try {
768 randomvector(P, lambda, dim);
769 for (j = 0; j < P->NbRays; ++j) {
770 Polyhedron *C = supporting_cone(P, j);
771 decompose(C, MaxRays);
773 break;
774 } catch (OrthogonalException &e) {
775 mpq_set_si(count, 0, 0);
780 /* base for non-parametric counting */
781 struct np_base : public polar_decomposer {
782 int current_vertex;
783 Polyhedron *P;
784 unsigned dim;
786 np_base(Polyhedron *P, unsigned dim) {
787 this->P = P;
788 this->dim = dim;
792 struct reducer : public np_base {
793 vec_ZZ vertex;
794 //vec_ZZ den;
795 ZZ sgn;
796 ZZ num;
797 ZZ one;
798 mpq_t tcount;
799 mpz_t tn;
800 mpz_t td;
801 int lower; // call base when only this many variables is left
803 reducer(Polyhedron *P) : np_base(P, P->Dimension) {
804 //den.SetLength(dim);
805 mpq_init(tcount);
806 mpz_init(tn);
807 mpz_init(td);
808 one = 1;
811 void start(unsigned MaxRays);
813 ~reducer() {
814 mpq_clear(tcount);
815 mpz_clear(tn);
816 mpz_clear(td);
819 virtual void handle_polar(Polyhedron *P, int sign);
820 void reduce(ZZ c, ZZ cd, vec_ZZ& num, mat_ZZ& den_f);
821 virtual void base(ZZ& c, ZZ& cd, vec_ZZ& num, mat_ZZ& den_f) = 0;
822 virtual void split(vec_ZZ& num, ZZ& num_s, vec_ZZ& num_p,
823 mat_ZZ& den_f, vec_ZZ& den_s, mat_ZZ& den_r) = 0;
826 void reducer::reduce(ZZ c, ZZ cd, vec_ZZ& num, mat_ZZ& den_f)
828 unsigned len = den_f.NumRows(); // number of factors in den
830 if (num.length() == lower) {
831 base(c, cd, num, den_f);
832 return;
834 assert(num.length() > 1);
836 vec_ZZ den_s;
837 mat_ZZ den_r;
838 ZZ num_s;
839 vec_ZZ num_p;
841 split(num, num_s, num_p, den_f, den_s, den_r);
843 vec_ZZ den_p;
844 den_p.SetLength(len);
846 normalize(c, num_s, num_p, den_s, den_p, den_r);
848 int only_param = 0; // k-r-s from text
849 int no_param = 0; // r from text
850 for (int k = 0; k < len; ++k) {
851 if (den_p[k] == 0)
852 ++no_param;
853 else if (den_s[k] == 0)
854 ++only_param;
856 if (no_param == 0) {
857 reduce(c, cd, num_p, den_r);
858 } else {
859 int k, l;
860 mat_ZZ pden;
861 pden.SetDims(only_param, den_r.NumCols());
863 for (k = 0, l = 0; k < len; ++k)
864 if (den_s[k] == 0)
865 pden[l++] = den_r[k];
867 for (k = 0; k < len; ++k)
868 if (den_p[k] == 0)
869 break;
871 dpoly n(no_param, num_s);
872 dpoly D(no_param, den_s[k], 1);
873 for ( ; ++k < len; )
874 if (den_p[k] == 0) {
875 dpoly fact(no_param, den_s[k], 1);
876 D *= fact;
879 if (no_param + only_param == len) {
880 mpq_set_si(tcount, 0, 1);
881 n.div(D, tcount, one);
883 ZZ qn, qd;
884 value2zz(mpq_numref(tcount), qn);
885 value2zz(mpq_denref(tcount), qd);
887 qn *= c;
888 qd *= cd;
890 if (qn != 0)
891 reduce(qn, qd, num_p, pden);
892 } else {
893 dpoly_r * r = 0;
895 for (k = 0; k < len; ++k) {
896 if (den_s[k] == 0 || den_p[k] == 0)
897 continue;
899 dpoly pd(no_param-1, den_s[k], 1);
901 int l;
902 for (l = 0; l < k; ++l)
903 if (den_r[l] == den_r[k])
904 break;
906 if (r == 0)
907 r = new dpoly_r(n, pd, l, len);
908 else {
909 dpoly_r *nr = new dpoly_r(r, pd, l, len);
910 delete r;
911 r = nr;
915 dpoly_r *rc = r->div(D);
917 rc->denom *= cd;
919 int common = pden.NumRows();
920 vector< dpoly_r_term * >& final = rc->c[rc->len-1];
921 int rows;
922 for (int j = 0; j < final.size(); ++j) {
923 if (final[j]->coeff == 0)
924 continue;
925 rows = common;
926 pden.SetDims(rows, pden.NumCols());
927 for (int k = 0; k < rc->dim; ++k) {
928 int n = final[j]->powers[k];
929 if (n == 0)
930 continue;
931 pden.SetDims(rows+n, pden.NumCols());
932 for (int l = 0; l < n; ++l)
933 pden[rows+l] = den_r[k];
934 rows += n;
936 final[j]->coeff *= c;
937 reduce(final[j]->coeff, rc->denom, num_p, pden);
940 delete rc;
941 delete r;
946 void reducer::handle_polar(Polyhedron *C, int s)
948 assert(C->NbRays-1 == dim);
950 sgn = s;
952 lattice_point(P->Ray[current_vertex]+1, C, vertex);
954 mat_ZZ den;
955 den.SetDims(dim, dim);
957 int r;
958 for (r = 0; r < dim; ++r)
959 values2zz(C->Ray[r]+1, den[r], dim);
961 reduce(sgn, one, vertex, den);
964 void reducer::start(unsigned MaxRays)
966 for (current_vertex = 0; current_vertex < P->NbRays; ++current_vertex) {
967 Polyhedron *C = supporting_cone(P, current_vertex);
968 decompose(C, MaxRays);
972 struct ireducer : public reducer {
973 ireducer(Polyhedron *P) : reducer(P) {}
975 virtual void split(vec_ZZ& num, ZZ& num_s, vec_ZZ& num_p,
976 mat_ZZ& den_f, vec_ZZ& den_s, mat_ZZ& den_r) {
977 unsigned len = den_f.NumRows(); // number of factors in den
978 unsigned d = num.length() - 1;
980 den_s.SetLength(len);
981 den_r.SetDims(len, d);
983 for (int r = 0; r < len; ++r) {
984 den_s[r] = den_f[r][0];
985 for (int k = 1; k <= d; ++k)
986 den_r[r][k-1] = den_f[r][k];
989 num_s = num[0];
990 num_p.SetLength(d);
991 for (int k = 1 ; k <= d; ++k)
992 num_p[k-1] = num[k];
996 // incremental counter
997 struct icounter : public ireducer {
998 mpq_t count;
1000 icounter(Polyhedron *P) : ireducer(P) {
1001 mpq_init(count);
1002 lower = 1;
1004 ~icounter() {
1005 mpq_clear(count);
1007 virtual void base(ZZ& c, ZZ& cd, vec_ZZ& num, mat_ZZ& den_f);
1010 void icounter::base(ZZ& c, ZZ& cd, vec_ZZ& num, mat_ZZ& den_f)
1012 int r;
1013 unsigned len = den_f.NumRows(); // number of factors in den
1014 vec_ZZ den_s;
1015 den_s.SetLength(len);
1016 ZZ num_s = num[0];
1017 for (r = 0; r < len; ++r)
1018 den_s[r] = den_f[r][0];
1019 normalize(c, num_s, den_s);
1021 dpoly n(len, num_s);
1022 dpoly D(len, den_s[0], 1);
1023 for (int k = 1; k < len; ++k) {
1024 dpoly fact(len, den_s[k], 1);
1025 D *= fact;
1027 mpq_set_si(tcount, 0, 1);
1028 n.div(D, tcount, one);
1029 zz2value(c, tn);
1030 zz2value(cd, td);
1031 mpz_mul(mpq_numref(tcount), mpq_numref(tcount), tn);
1032 mpz_mul(mpq_denref(tcount), mpq_denref(tcount), td);
1033 mpq_canonicalize(tcount);
1034 mpq_add(count, count, tcount);
1037 /* base for generating function counting */
1038 struct gf_base {
1039 np_base *base;
1040 gen_fun *gf;
1042 gf_base(np_base *npb, unsigned nparam) : base(npb) {
1043 gf = new gen_fun(Polyhedron_Project(base->P, nparam));
1045 void start(unsigned MaxRays);
1048 void gf_base::start(unsigned MaxRays)
1050 for (int i = 0; i < base->P->NbRays; ++i) {
1051 if (!value_pos_p(base->P->Ray[i][base->dim+1]))
1052 continue;
1054 Polyhedron *C = supporting_cone(base->P, i);
1055 base->current_vertex = i;
1056 base->decompose(C, MaxRays);
1060 struct partial_ireducer : public ireducer, public gf_base {
1061 partial_ireducer(Polyhedron *P, unsigned nparam) :
1062 ireducer(P), gf_base(this, nparam) {
1063 lower = nparam;
1065 ~partial_ireducer() {
1067 virtual void base(ZZ& c, ZZ& cd, vec_ZZ& num, mat_ZZ& den_f);
1068 /* we want to override the start method from reducer with the one from gf_base */
1069 void start(unsigned MaxRays) {
1070 gf_base::start(MaxRays);
1074 void partial_ireducer::base(ZZ& c, ZZ& cd, vec_ZZ& num, mat_ZZ& den_f)
1076 gf->add(c, cd, num, den_f);
1079 struct partial_reducer : public reducer, public gf_base {
1080 vec_ZZ lambda;
1081 vec_ZZ tmp;
1083 partial_reducer(Polyhedron *P, unsigned nparam) :
1084 reducer(P), gf_base(this, nparam) {
1085 lower = nparam;
1087 tmp.SetLength(dim - nparam);
1088 randomvector(P, lambda, dim - nparam);
1090 ~partial_reducer() {
1092 virtual void base(ZZ& c, ZZ& cd, vec_ZZ& num, mat_ZZ& den_f);
1093 /* we want to override the start method from reducer with the one from gf_base */
1094 void start(unsigned MaxRays) {
1095 gf_base::start(MaxRays);
1098 virtual void split(vec_ZZ& num, ZZ& num_s, vec_ZZ& num_p,
1099 mat_ZZ& den_f, vec_ZZ& den_s, mat_ZZ& den_r) {
1100 unsigned len = den_f.NumRows(); // number of factors in den
1101 unsigned nvar = tmp.length();
1103 den_s.SetLength(len);
1104 den_r.SetDims(len, lower);
1106 for (int r = 0; r < len; ++r) {
1107 for (int k = 0; k < nvar; ++k)
1108 tmp[k] = den_f[r][k];
1109 den_s[r] = tmp * lambda;
1111 for (int k = nvar; k < dim; ++k)
1112 den_r[r][k-nvar] = den_f[r][k];
1115 for (int k = 0; k < nvar; ++k)
1116 tmp[k] = num[k];
1117 num_s = tmp *lambda;
1118 num_p.SetLength(lower);
1119 for (int k = nvar ; k < dim; ++k)
1120 num_p[k-nvar] = num[k];
1124 void partial_reducer::base(ZZ& c, ZZ& cd, vec_ZZ& num, mat_ZZ& den_f)
1126 gf->add(c, cd, num, den_f);
1129 struct bfc_term_base {
1130 // the number of times a given factor appears in the denominator
1131 int *powers;
1132 mat_ZZ terms;
1134 bfc_term_base(int len) {
1135 powers = new int[len];
1138 virtual ~bfc_term_base() {
1139 delete [] powers;
1143 struct bfc_term : public bfc_term_base {
1144 vec_ZZ cn;
1145 vec_ZZ cd;
1147 bfc_term(int len) : bfc_term_base(len) {}
1150 struct bfe_term : public bfc_term_base {
1151 vector<evalue *> factors;
1153 bfe_term(int len) : bfc_term_base(len) {
1156 ~bfe_term() {
1157 for (int i = 0; i < factors.size(); ++i) {
1158 if (!factors[i])
1159 continue;
1160 free_evalue_refs(factors[i]);
1161 delete factors[i];
1166 typedef vector< bfc_term_base * > bfc_vec;
1168 struct bf_reducer;
1170 struct bf_base : public np_base {
1171 ZZ one;
1172 mpq_t tcount;
1173 mpz_t tn;
1174 mpz_t td;
1175 int lower; // call base when only this many variables is left
1177 bf_base(Polyhedron *P, unsigned dim) : np_base(P, dim) {
1178 mpq_init(tcount);
1179 mpz_init(tn);
1180 mpz_init(td);
1181 one = 1;
1184 ~bf_base() {
1185 mpq_clear(tcount);
1186 mpz_clear(tn);
1187 mpz_clear(td);
1190 void start(unsigned MaxRays);
1191 virtual void handle_polar(Polyhedron *P, int sign);
1192 int setup_factors(Polyhedron *P, mat_ZZ& factors, bfc_term_base* t, int s);
1194 bfc_term_base* find_bfc_term(bfc_vec& v, int *powers, int len);
1195 void add_term(bfc_term_base *t, vec_ZZ& num1, vec_ZZ& num);
1196 void add_term(bfc_term_base *t, vec_ZZ& num);
1198 void reduce(mat_ZZ& factors, bfc_vec& v);
1199 virtual void base(mat_ZZ& factors, bfc_vec& v) = 0;
1201 virtual bfc_term_base* new_bf_term(int len) = 0;
1202 virtual void set_factor(bfc_term_base *t, int k, int change) = 0;
1203 virtual void set_factor(bfc_term_base *t, int k, mpq_t &f, int change) = 0;
1204 virtual void set_factor(bfc_term_base *t, int k, ZZ& n, ZZ& d, int change) = 0;
1205 virtual void update_term(bfc_term_base *t, int i) = 0;
1206 virtual void insert_term(bfc_term_base *t, int i) = 0;
1207 virtual bool constant_vertex(int dim) = 0;
1208 virtual void cum(bf_reducer *bfr, bfc_term_base *t, int k,
1209 dpoly_r *r) = 0;
1212 static int lex_cmp(vec_ZZ& a, vec_ZZ& b)
1214 assert(a.length() == b.length());
1216 for (int j = 0; j < a.length(); ++j)
1217 if (a[j] != b[j])
1218 return a[j] < b[j] ? -1 : 1;
1219 return 0;
1222 void bf_base::add_term(bfc_term_base *t, vec_ZZ& num_orig, vec_ZZ& extra_num)
1224 vec_ZZ num;
1225 int d = num_orig.length();
1226 num.SetLength(d-1);
1227 for (int l = 0; l < d-1; ++l)
1228 num[l] = num_orig[l+1] + extra_num[l];
1230 add_term(t, num);
1233 void bf_base::add_term(bfc_term_base *t, vec_ZZ& num)
1235 int len = t->terms.NumRows();
1236 int i, r;
1237 for (i = 0; i < len; ++i) {
1238 r = lex_cmp(t->terms[i], num);
1239 if (r >= 0)
1240 break;
1242 if (i == len || r > 0) {
1243 t->terms.SetDims(len+1, num.length());
1244 insert_term(t, i);
1245 t->terms[i] = num;
1246 } else {
1247 // i < len && r == 0
1248 update_term(t, i);
1252 static void print_int_vector(int *v, int len, char *name)
1254 cerr << name << endl;
1255 for (int j = 0; j < len; ++j) {
1256 cerr << v[j] << " ";
1258 cerr << endl;
1261 static void print_bfc_terms(mat_ZZ& factors, bfc_vec& v)
1263 cerr << endl;
1264 cerr << "factors" << endl;
1265 cerr << factors << endl;
1266 for (int i = 0; i < v.size(); ++i) {
1267 cerr << "term: " << i << endl;
1268 print_int_vector(v[i]->powers, factors.NumRows(), "powers");
1269 cerr << "terms" << endl;
1270 cerr << v[i]->terms << endl;
1271 bfc_term* bfct = static_cast<bfc_term *>(v[i]);
1272 cerr << bfct->cn << endl;
1273 cerr << bfct->cd << endl;
1277 static void print_bfe_terms(mat_ZZ& factors, bfc_vec& v)
1279 cerr << endl;
1280 cerr << "factors" << endl;
1281 cerr << factors << endl;
1282 for (int i = 0; i < v.size(); ++i) {
1283 cerr << "term: " << i << endl;
1284 print_int_vector(v[i]->powers, factors.NumRows(), "powers");
1285 cerr << "terms" << endl;
1286 cerr << v[i]->terms << endl;
1287 bfe_term* bfet = static_cast<bfe_term *>(v[i]);
1288 for (int j = 0; j < v[i]->terms.NumRows(); ++j) {
1289 char * test[] = {"a", "b"};
1290 print_evalue(stderr, bfet->factors[j], test);
1291 fprintf(stderr, "\n");
1296 bfc_term_base* bf_base::find_bfc_term(bfc_vec& v, int *powers, int len)
1298 bfc_vec::iterator i;
1299 for (i = v.begin(); i != v.end(); ++i) {
1300 int j;
1301 for (j = 0; j < len; ++j)
1302 if ((*i)->powers[j] != powers[j])
1303 break;
1304 if (j == len)
1305 return (*i);
1306 if ((*i)->powers[j] > powers[j])
1307 break;
1310 bfc_term_base* t = new_bf_term(len);
1311 v.insert(i, t);
1312 memcpy(t->powers, powers, len * sizeof(int));
1314 return t;
1317 struct bf_reducer {
1318 mat_ZZ& factors;
1319 bfc_vec& v;
1320 bf_base *bf;
1322 unsigned nf;
1323 unsigned d;
1325 mat_ZZ nfactors;
1326 int *old2new;
1327 int *sign;
1328 unsigned int nnf;
1329 bfc_vec vn;
1331 vec_ZZ extra_num;
1332 int changes;
1333 int no_param; // r from text
1334 int only_param; // k-r-s from text
1335 int total_power; // k from text
1337 // created in compute_reduced_factors
1338 int *bpowers;
1339 // set in update_powers
1340 int *npowers;
1341 vec_ZZ l_extra_num;
1342 int l_changes;
1344 bf_reducer(mat_ZZ& factors, bfc_vec& v, bf_base *bf)
1345 : factors(factors), v(v), bf(bf) {
1346 nf = factors.NumRows();
1347 d = factors.NumCols();
1348 old2new = new int[nf];
1349 sign = new int[nf];
1351 extra_num.SetLength(d-1);
1353 ~bf_reducer() {
1354 delete [] old2new;
1355 delete [] sign;
1356 delete [] npowers;
1357 delete [] bpowers;
1360 void compute_reduced_factors();
1361 void compute_extra_num(int i);
1363 void reduce();
1365 void update_powers(int *powers, int len);
1368 void bf_reducer::compute_extra_num(int i)
1370 clear(extra_num);
1371 changes = 0;
1372 no_param = 0; // r from text
1373 only_param = 0; // k-r-s from text
1374 total_power = 0; // k from text
1376 for (int j = 0; j < nf; ++j) {
1377 if (v[i]->powers[j] == 0)
1378 continue;
1380 total_power += v[i]->powers[j];
1381 if (factors[j][0] == 0) {
1382 only_param += v[i]->powers[j];
1383 continue;
1386 if (old2new[j] == -1)
1387 no_param += v[i]->powers[j];
1388 else
1389 extra_num += -sign[j] * v[i]->powers[j] * nfactors[old2new[j]];
1390 changes += v[i]->powers[j];
1394 void bf_reducer::update_powers(int *powers, int len)
1396 for (int l = 0; l < nnf; ++l)
1397 npowers[l] = bpowers[l];
1399 l_extra_num = extra_num;
1400 l_changes = changes;
1402 for (int l = 0; l < len; ++l) {
1403 int n = powers[l];
1404 if (n == 0)
1405 continue;
1406 assert(old2new[l] != -1);
1408 npowers[old2new[l]] += n;
1409 // interpretation of sign has been inverted
1410 // since we inverted the power for specialization
1411 if (sign[l] == 1) {
1412 l_extra_num += n * nfactors[old2new[l]];
1413 l_changes += n;
1419 void bf_reducer::compute_reduced_factors()
1421 unsigned nf = factors.NumRows();
1422 unsigned d = factors.NumCols();
1423 nnf = 0;
1424 nfactors.SetDims(nnf, d-1);
1426 for (int i = 0; i < nf; ++i) {
1427 int j;
1428 int s = 1;
1429 for (j = 0; j < nnf; ++j) {
1430 int k;
1431 for (k = 1; k < d; ++k)
1432 if (factors[i][k] != 0 || nfactors[j][k-1] != 0)
1433 break;
1434 if (k < d && factors[i][k] == -nfactors[j][k-1])
1435 s = -1;
1436 for (; k < d; ++k)
1437 if (factors[i][k] != s * nfactors[j][k-1])
1438 break;
1439 if (k == d)
1440 break;
1442 old2new[i] = j;
1443 if (j == nnf) {
1444 int k;
1445 for (k = 1; k < d; ++k)
1446 if (factors[i][k] != 0)
1447 break;
1448 if (k < d) {
1449 if (factors[i][k] < 0)
1450 s = -1;
1451 nfactors.SetDims(++nnf, d-1);
1452 for (int k = 1; k < d; ++k)
1453 nfactors[j][k-1] = s * factors[i][k];
1454 } else
1455 old2new[i] = -1;
1457 sign[i] = s;
1459 npowers = new int[nnf];
1460 bpowers = new int[nnf];
1463 void bf_reducer::reduce()
1465 compute_reduced_factors();
1467 for (int i = 0; i < v.size(); ++i) {
1468 compute_extra_num(i);
1470 if (no_param == 0) {
1471 vec_ZZ extra_num;
1472 extra_num.SetLength(d-1);
1473 int changes = 0;
1474 int npowers[nnf];
1475 for (int k = 0; k < nnf; ++k)
1476 npowers[k] = 0;
1477 for (int k = 0; k < nf; ++k) {
1478 assert(old2new[k] != -1);
1479 npowers[old2new[k]] += v[i]->powers[k];
1480 if (sign[k] == -1) {
1481 extra_num += v[i]->powers[k] * nfactors[old2new[k]];
1482 changes += v[i]->powers[k];
1486 bfc_term_base * t = bf->find_bfc_term(vn, npowers, nnf);
1487 for (int k = 0; k < v[i]->terms.NumRows(); ++k) {
1488 bf->set_factor(v[i], k, changes % 2);
1489 bf->add_term(t, v[i]->terms[k], extra_num);
1491 } else {
1492 // powers of "constant" part
1493 for (int k = 0; k < nnf; ++k)
1494 bpowers[k] = 0;
1495 for (int k = 0; k < nf; ++k) {
1496 if (factors[k][0] != 0)
1497 continue;
1498 assert(old2new[k] != -1);
1499 bpowers[old2new[k]] += v[i]->powers[k];
1500 if (sign[k] == -1) {
1501 extra_num += v[i]->powers[k] * nfactors[old2new[k]];
1502 changes += v[i]->powers[k];
1506 int j;
1507 for (j = 0; j < nf; ++j)
1508 if (old2new[j] == -1 && v[i]->powers[j] > 0)
1509 break;
1511 dpoly D(no_param, factors[j][0], 1);
1512 for (int k = 1; k < v[i]->powers[j]; ++k) {
1513 dpoly fact(no_param, factors[j][0], 1);
1514 D *= fact;
1516 for ( ; ++j < nf; )
1517 if (old2new[j] == -1)
1518 for (int k = 0; k < v[i]->powers[j]; ++k) {
1519 dpoly fact(no_param, factors[j][0], 1);
1520 D *= fact;
1523 if (no_param + only_param == total_power &&
1524 bf->constant_vertex(d)) {
1525 bfc_term_base * t = NULL;
1526 vec_ZZ num;
1527 num.SetLength(d-1);
1528 ZZ cn;
1529 ZZ cd;
1530 for (int k = 0; k < v[i]->terms.NumRows(); ++k) {
1531 dpoly n(no_param, v[i]->terms[k][0]);
1532 mpq_set_si(bf->tcount, 0, 1);
1533 n.div(D, bf->tcount, bf->one);
1535 if (value_zero_p(mpq_numref(bf->tcount)))
1536 continue;
1538 if (!t)
1539 t = bf->find_bfc_term(vn, bpowers, nnf);
1540 bf->set_factor(v[i], k, bf->tcount, changes % 2);
1541 bf->add_term(t, v[i]->terms[k], extra_num);
1543 } else {
1544 for (int j = 0; j < v[i]->terms.NumRows(); ++j) {
1545 dpoly n(no_param, v[i]->terms[j][0]);
1547 dpoly_r * r = 0;
1548 if (no_param + only_param == total_power)
1549 r = new dpoly_r(n, nf);
1550 else
1551 for (int k = 0; k < nf; ++k) {
1552 if (v[i]->powers[k] == 0)
1553 continue;
1554 if (factors[k][0] == 0 || old2new[k] == -1)
1555 continue;
1557 dpoly pd(no_param-1, factors[k][0], 1);
1559 for (int l = 0; l < v[i]->powers[k]; ++l) {
1560 int q;
1561 for (q = 0; q < k; ++q)
1562 if (old2new[q] == old2new[k] &&
1563 sign[q] == sign[k])
1564 break;
1566 if (r == 0)
1567 r = new dpoly_r(n, pd, q, nf);
1568 else {
1569 dpoly_r *nr = new dpoly_r(r, pd, q, nf);
1570 delete r;
1571 r = nr;
1576 dpoly_r *rc = r->div(D);
1577 delete r;
1579 if (bf->constant_vertex(d)) {
1580 vector< dpoly_r_term * >& final = rc->c[rc->len-1];
1582 for (int k = 0; k < final.size(); ++k) {
1583 if (final[k]->coeff == 0)
1584 continue;
1586 update_powers(final[k]->powers, rc->dim);
1588 bfc_term_base * t = bf->find_bfc_term(vn, npowers, nnf);
1589 bf->set_factor(v[i], j, final[k]->coeff, rc->denom, l_changes % 2);
1590 bf->add_term(t, v[i]->terms[j], l_extra_num);
1592 } else
1593 bf->cum(this, v[i], j, rc);
1595 delete rc;
1599 delete v[i];
1604 void bf_base::reduce(mat_ZZ& factors, bfc_vec& v)
1606 assert(v.size() > 0);
1607 unsigned nf = factors.NumRows();
1608 unsigned d = factors.NumCols();
1610 if (d == lower)
1611 return base(factors, v);
1613 bf_reducer bfr(factors, v, this);
1615 bfr.reduce();
1617 if (bfr.vn.size() > 0)
1618 reduce(bfr.nfactors, bfr.vn);
1621 int bf_base::setup_factors(Polyhedron *C, mat_ZZ& factors,
1622 bfc_term_base* t, int s)
1624 factors.SetDims(dim, dim);
1626 int r;
1628 for (r = 0; r < dim; ++r)
1629 t->powers[r] = 1;
1631 for (r = 0; r < dim; ++r) {
1632 values2zz(C->Ray[r]+1, factors[r], dim);
1633 int k;
1634 for (k = 0; k < dim; ++k)
1635 if (factors[r][k] != 0)
1636 break;
1637 if (factors[r][k] < 0) {
1638 factors[r] = -factors[r];
1639 t->terms[0] += factors[r];
1640 s = -s;
1644 return s;
1647 void bf_base::handle_polar(Polyhedron *C, int s)
1649 bfc_term* t = new bfc_term(dim);
1650 vector< bfc_term_base * > v;
1651 v.push_back(t);
1653 assert(C->NbRays-1 == dim);
1655 t->cn.SetLength(1);
1656 t->cd.SetLength(1);
1658 t->terms.SetDims(1, dim);
1659 lattice_point(P->Ray[current_vertex]+1, C, t->terms[0]);
1661 // the elements of factors are always lexpositive
1662 mat_ZZ factors;
1663 s = setup_factors(C, factors, t, s);
1665 t->cn[0] = s;
1666 t->cd[0] = 1;
1668 reduce(factors, v);
1671 void bf_base::start(unsigned MaxRays)
1673 for (current_vertex = 0; current_vertex < P->NbRays; ++current_vertex) {
1674 Polyhedron *C = supporting_cone(P, current_vertex);
1675 decompose(C, MaxRays);
1679 struct bfcounter_base : public bf_base {
1680 ZZ cn;
1681 ZZ cd;
1683 bfcounter_base(Polyhedron *P) : bf_base(P, P->Dimension) {
1686 bfc_term_base* new_bf_term(int len) {
1687 bfc_term* t = new bfc_term(len);
1688 t->cn.SetLength(0);
1689 t->cd.SetLength(0);
1690 return t;
1693 virtual void set_factor(bfc_term_base *t, int k, int change) {
1694 bfc_term* bfct = static_cast<bfc_term *>(t);
1695 cn = bfct->cn[k];
1696 if (change)
1697 cn = -cn;
1698 cd = bfct->cd[k];
1701 virtual void set_factor(bfc_term_base *t, int k, mpq_t &f, int change) {
1702 bfc_term* bfct = static_cast<bfc_term *>(t);
1703 value2zz(mpq_numref(f), cn);
1704 value2zz(mpq_denref(f), cd);
1705 cn *= bfct->cn[k];
1706 if (change)
1707 cn = -cn;
1708 cd *= bfct->cd[k];
1711 virtual void set_factor(bfc_term_base *t, int k, ZZ& n, ZZ& d, int change) {
1712 bfc_term* bfct = static_cast<bfc_term *>(t);
1713 cn = bfct->cn[k] * n;
1714 if (change)
1715 cn = -cn;
1716 cd = bfct->cd[k] * d;
1719 virtual void insert_term(bfc_term_base *t, int i) {
1720 bfc_term* bfct = static_cast<bfc_term *>(t);
1721 int len = t->terms.NumRows()-1; // already increased by one
1723 bfct->cn.SetLength(len+1);
1724 bfct->cd.SetLength(len+1);
1725 for (int j = len; j > i; --j) {
1726 bfct->cn[j] = bfct->cn[j-1];
1727 bfct->cd[j] = bfct->cd[j-1];
1728 t->terms[j] = t->terms[j-1];
1730 bfct->cn[i] = cn;
1731 bfct->cd[i] = cd;
1734 virtual void update_term(bfc_term_base *t, int i) {
1735 bfc_term* bfct = static_cast<bfc_term *>(t);
1737 ZZ g = GCD(bfct->cd[i], cd);
1738 ZZ n = cn * bfct->cd[i]/g + bfct->cn[i] * cd/g;
1739 ZZ d = bfct->cd[i] * cd / g;
1740 bfct->cn[i] = n;
1741 bfct->cd[i] = d;
1744 virtual bool constant_vertex(int dim) { return true; }
1745 virtual void cum(bf_reducer *bfr, bfc_term_base *t, int k, dpoly_r *r) {
1746 assert(0);
1750 struct bfcounter : public bfcounter_base {
1751 mpq_t count;
1753 bfcounter(Polyhedron *P) : bfcounter_base(P) {
1754 mpq_init(count);
1755 lower = 1;
1757 ~bfcounter() {
1758 mpq_clear(count);
1760 virtual void base(mat_ZZ& factors, bfc_vec& v);
1763 void bfcounter::base(mat_ZZ& factors, bfc_vec& v)
1765 unsigned nf = factors.NumRows();
1767 for (int i = 0; i < v.size(); ++i) {
1768 bfc_term* bfct = static_cast<bfc_term *>(v[i]);
1769 int total_power = 0;
1770 // factor is always positive, so we always
1771 // change signs
1772 for (int k = 0; k < nf; ++k)
1773 total_power += v[i]->powers[k];
1775 int j;
1776 for (j = 0; j < nf; ++j)
1777 if (v[i]->powers[j] > 0)
1778 break;
1780 dpoly D(total_power, factors[j][0], 1);
1781 for (int k = 1; k < v[i]->powers[j]; ++k) {
1782 dpoly fact(total_power, factors[j][0], 1);
1783 D *= fact;
1785 for ( ; ++j < nf; )
1786 for (int k = 0; k < v[i]->powers[j]; ++k) {
1787 dpoly fact(total_power, factors[j][0], 1);
1788 D *= fact;
1791 for (int k = 0; k < v[i]->terms.NumRows(); ++k) {
1792 dpoly n(total_power, v[i]->terms[k][0]);
1793 mpq_set_si(tcount, 0, 1);
1794 n.div(D, tcount, one);
1795 if (total_power % 2)
1796 bfct->cn[k] = -bfct->cn[k];
1797 zz2value(bfct->cn[k], tn);
1798 zz2value(bfct->cd[k], td);
1800 mpz_mul(mpq_numref(tcount), mpq_numref(tcount), tn);
1801 mpz_mul(mpq_denref(tcount), mpq_denref(tcount), td);
1802 mpq_canonicalize(tcount);
1803 mpq_add(count, count, tcount);
1805 delete v[i];
1809 struct partial_bfcounter : public bfcounter_base, public gf_base {
1810 partial_bfcounter(Polyhedron *P, unsigned nparam) :
1811 bfcounter_base(P), gf_base(this, nparam) {
1812 lower = nparam;
1814 ~partial_bfcounter() {
1816 virtual void base(mat_ZZ& factors, bfc_vec& v);
1817 /* we want to override the start method from bf_base with the one from gf_base */
1818 void start(unsigned MaxRays) {
1819 gf_base::start(MaxRays);
1823 void partial_bfcounter::base(mat_ZZ& factors, bfc_vec& v)
1825 mat_ZZ den;
1826 unsigned nf = factors.NumRows();
1828 for (int i = 0; i < v.size(); ++i) {
1829 bfc_term* bfct = static_cast<bfc_term *>(v[i]);
1830 den.SetDims(0, lower);
1831 int total_power = 0;
1832 int p = 0;
1833 for (int j = 0; j < nf; ++j) {
1834 total_power += v[i]->powers[j];
1835 den.SetDims(total_power, lower);
1836 for (int k = 0; k < v[i]->powers[j]; ++k)
1837 den[p++] = factors[j];
1839 for (int j = 0; j < v[i]->terms.NumRows(); ++j)
1840 gf->add(bfct->cn[j], bfct->cd[j], v[i]->terms[j], den);
1841 delete v[i];
1846 typedef Polyhedron * Polyhedron_p;
1848 static void barvinok_count_f(Polyhedron *P, Value* result, unsigned NbMaxCons);
1850 void barvinok_count(Polyhedron *P, Value* result, unsigned NbMaxCons)
1852 unsigned dim;
1853 int allocated = 0;
1854 Polyhedron *Q;
1855 int r = 0;
1856 bool infinite = false;
1858 if (emptyQ2(P)) {
1859 value_set_si(*result, 0);
1860 return;
1862 if (P->NbBid == 0)
1863 for (; r < P->NbRays; ++r)
1864 if (value_zero_p(P->Ray[r][P->Dimension+1]))
1865 break;
1866 if (P->NbBid !=0 || r < P->NbRays) {
1867 value_set_si(*result, -1);
1868 return;
1870 if (P->NbEq != 0) {
1871 P = remove_equalities(P);
1872 if (emptyQ(P)) {
1873 Polyhedron_Free(P);
1874 value_set_si(*result, 0);
1875 return;
1877 allocated = 1;
1879 if (P->Dimension == 0) {
1880 /* Test whether the constraints are satisfied */
1881 POL_ENSURE_VERTICES(P);
1882 value_set_si(*result, !emptyQ(P));
1883 if (allocated)
1884 Polyhedron_Free(P);
1885 return;
1887 Q = Polyhedron_Factor(P, 0, NbMaxCons);
1888 if (Q) {
1889 if (allocated)
1890 Polyhedron_Free(P);
1891 P = Q;
1892 allocated = 1;
1895 barvinok_count_f(P, result, NbMaxCons);
1896 if (Q && P->next) {
1897 Value factor;
1898 value_init(factor);
1900 for (Q = P->next; Q; Q = Q->next) {
1901 barvinok_count_f(Q, &factor, NbMaxCons);
1902 if (value_neg_p(factor)) {
1903 infinite = true;
1904 continue;
1905 } else if (Q->next && value_zero_p(factor)) {
1906 value_set_si(*result, 0);
1907 break;
1909 value_multiply(*result, *result, factor);
1912 value_clear(factor);
1915 if (allocated)
1916 Domain_Free(P);
1917 if (infinite)
1918 value_set_si(*result, -1);
1921 static void barvinok_count_f(Polyhedron *P, Value* result, unsigned NbMaxCons)
1923 if (P->Dimension == 1)
1924 return Line_Length(P, result);
1926 int c = P->NbConstraints;
1927 POL_ENSURE_FACETS(P);
1928 if (c != P->NbConstraints || P->NbEq != 0)
1929 return barvinok_count(P, result, NbMaxCons);
1931 POL_ENSURE_VERTICES(P);
1933 #ifdef USE_INCREMENTAL_BF
1934 bfcounter cnt(P);
1935 #elif defined USE_INCREMENTAL_DF
1936 icounter cnt(P);
1937 #else
1938 counter cnt(P);
1939 #endif
1940 cnt.start(NbMaxCons);
1942 assert(value_one_p(&cnt.count[0]._mp_den));
1943 value_assign(*result, &cnt.count[0]._mp_num);
1946 static void uni_polynom(int param, Vector *c, evalue *EP)
1948 unsigned dim = c->Size-2;
1949 value_init(EP->d);
1950 value_set_si(EP->d,0);
1951 EP->x.p = new_enode(polynomial, dim+1, param+1);
1952 for (int j = 0; j <= dim; ++j)
1953 evalue_set(&EP->x.p->arr[j], c->p[j], c->p[dim+1]);
1956 static void multi_polynom(Vector *c, evalue* X, evalue *EP)
1958 unsigned dim = c->Size-2;
1959 evalue EC;
1961 value_init(EC.d);
1962 evalue_set(&EC, c->p[dim], c->p[dim+1]);
1964 value_init(EP->d);
1965 evalue_set(EP, c->p[dim], c->p[dim+1]);
1967 for (int i = dim-1; i >= 0; --i) {
1968 emul(X, EP);
1969 value_assign(EC.x.n, c->p[i]);
1970 eadd(&EC, EP);
1972 free_evalue_refs(&EC);
1975 Polyhedron *unfringe (Polyhedron *P, unsigned MaxRays)
1977 int len = P->Dimension+2;
1978 Polyhedron *T, *R = P;
1979 Value g;
1980 value_init(g);
1981 Vector *row = Vector_Alloc(len);
1982 value_set_si(row->p[0], 1);
1984 R = DomainConstraintSimplify(Polyhedron_Copy(P), MaxRays);
1986 Matrix *M = Matrix_Alloc(2, len-1);
1987 value_set_si(M->p[1][len-2], 1);
1988 for (int v = 0; v < P->Dimension; ++v) {
1989 value_set_si(M->p[0][v], 1);
1990 Polyhedron *I = Polyhedron_Image(P, M, 2+1);
1991 value_set_si(M->p[0][v], 0);
1992 for (int r = 0; r < I->NbConstraints; ++r) {
1993 if (value_zero_p(I->Constraint[r][0]))
1994 continue;
1995 if (value_zero_p(I->Constraint[r][1]))
1996 continue;
1997 if (value_one_p(I->Constraint[r][1]))
1998 continue;
1999 if (value_mone_p(I->Constraint[r][1]))
2000 continue;
2001 value_absolute(g, I->Constraint[r][1]);
2002 Vector_Set(row->p+1, 0, len-2);
2003 value_division(row->p[1+v], I->Constraint[r][1], g);
2004 mpz_fdiv_q(row->p[len-1], I->Constraint[r][2], g);
2005 T = R;
2006 R = AddConstraints(row->p, 1, R, MaxRays);
2007 if (T != P)
2008 Polyhedron_Free(T);
2010 Polyhedron_Free(I);
2012 Matrix_Free(M);
2013 Vector_Free(row);
2014 value_clear(g);
2015 return R;
2018 /* this procedure may have false negatives */
2019 static bool Polyhedron_is_infinite(Polyhedron *P, unsigned nparam)
2021 int r;
2022 for (r = 0; r < P->NbRays; ++r) {
2023 if (!value_zero_p(P->Ray[r][0]) &&
2024 !value_zero_p(P->Ray[r][P->Dimension+1]))
2025 continue;
2026 if (First_Non_Zero(P->Ray[r]+1+P->Dimension-nparam, nparam) == -1)
2027 return true;
2029 return false;
2032 /* Check whether all rays point in the positive directions
2033 * for the parameters
2035 static bool Polyhedron_has_positive_rays(Polyhedron *P, unsigned nparam)
2037 int r;
2038 for (r = 0; r < P->NbRays; ++r)
2039 if (value_zero_p(P->Ray[r][P->Dimension+1])) {
2040 int i;
2041 for (i = P->Dimension - nparam; i < P->Dimension; ++i)
2042 if (value_neg_p(P->Ray[r][i+1]))
2043 return false;
2045 return true;
2048 typedef evalue * evalue_p;
2050 struct enumerator : public polar_decomposer {
2051 vec_ZZ lambda;
2052 unsigned dim, nbV;
2053 evalue ** vE;
2054 int _i;
2055 mat_ZZ rays;
2056 vec_ZZ den;
2057 ZZ sign;
2058 Polyhedron *P;
2059 Param_Vertices *V;
2060 term_info num;
2061 Vector *c;
2062 mpq_t count;
2064 enumerator(Polyhedron *P, unsigned dim, unsigned nbV) {
2065 this->P = P;
2066 this->dim = dim;
2067 this->nbV = nbV;
2068 randomvector(P, lambda, dim);
2069 rays.SetDims(dim, dim);
2070 den.SetLength(dim);
2071 c = Vector_Alloc(dim+2);
2073 vE = new evalue_p[nbV];
2074 for (int j = 0; j < nbV; ++j)
2075 vE[j] = 0;
2077 mpq_init(count);
2080 void decompose_at(Param_Vertices *V, int _i, unsigned MaxRays) {
2081 Polyhedron *C = supporting_cone_p(P, V);
2082 this->_i = _i;
2083 this->V = V;
2085 vE[_i] = new evalue;
2086 value_init(vE[_i]->d);
2087 evalue_set_si(vE[_i], 0, 1);
2089 decompose(C, MaxRays);
2092 ~enumerator() {
2093 mpq_clear(count);
2094 Vector_Free(c);
2096 for (int j = 0; j < nbV; ++j)
2097 if (vE[j]) {
2098 free_evalue_refs(vE[j]);
2099 delete vE[j];
2101 delete [] vE;
2104 virtual void handle_polar(Polyhedron *P, int sign);
2107 void enumerator::handle_polar(Polyhedron *C, int s)
2109 int r = 0;
2110 assert(C->NbRays-1 == dim);
2111 add_rays(rays, C, &r);
2112 for (int k = 0; k < dim; ++k) {
2113 if (lambda * rays[k] == 0)
2114 throw Orthogonal;
2117 sign = s;
2119 lattice_point(V, C, lambda, &num, 0);
2120 den = rays * lambda;
2121 normalize(sign, num.constant, den);
2123 dpoly n(dim, den[0], 1);
2124 for (int k = 1; k < dim; ++k) {
2125 dpoly fact(dim, den[k], 1);
2126 n *= fact;
2128 if (num.E != NULL) {
2129 ZZ one(INIT_VAL, 1);
2130 dpoly_n d(dim, num.constant, one);
2131 d.div(n, c, sign);
2132 evalue EV;
2133 multi_polynom(c, num.E, &EV);
2134 eadd(&EV , vE[_i]);
2135 free_evalue_refs(&EV);
2136 free_evalue_refs(num.E);
2137 delete num.E;
2138 } else if (num.pos != -1) {
2139 dpoly_n d(dim, num.constant, num.coeff);
2140 d.div(n, c, sign);
2141 evalue EV;
2142 uni_polynom(num.pos, c, &EV);
2143 eadd(&EV , vE[_i]);
2144 free_evalue_refs(&EV);
2145 } else {
2146 mpq_set_si(count, 0, 1);
2147 dpoly d(dim, num.constant);
2148 d.div(n, count, sign);
2149 evalue EV;
2150 value_init(EV.d);
2151 evalue_set(&EV, &count[0]._mp_num, &count[0]._mp_den);
2152 eadd(&EV , vE[_i]);
2153 free_evalue_refs(&EV);
2157 struct enumerator_base {
2158 unsigned dim;
2159 evalue ** vE;
2160 evalue ** E_vertex;
2161 evalue mone;
2162 vertex_decomposer *vpd;
2164 enumerator_base(unsigned dim, vertex_decomposer *vpd)
2166 this->dim = dim;
2167 this->vpd = vpd;
2169 vE = new evalue_p[vpd->nbV];
2170 for (int j = 0; j < vpd->nbV; ++j)
2171 vE[j] = 0;
2173 E_vertex = new evalue_p[dim];
2175 value_init(mone.d);
2176 evalue_set_si(&mone, -1, 1);
2179 void decompose_at(Param_Vertices *V, int _i, unsigned MaxRays/*, Polyhedron *pVD*/) {
2180 //this->pVD = pVD;
2182 vE[_i] = new evalue;
2183 value_init(vE[_i]->d);
2184 evalue_set_si(vE[_i], 0, 1);
2186 vpd->decompose_at_vertex(V, _i, MaxRays);
2189 ~enumerator_base() {
2190 for (int j = 0; j < vpd->nbV; ++j)
2191 if (vE[j]) {
2192 free_evalue_refs(vE[j]);
2193 delete vE[j];
2195 delete [] vE;
2197 delete [] E_vertex;
2199 free_evalue_refs(&mone);
2202 evalue *E_num(int i, int d) {
2203 return E_vertex[i + (dim-d)];
2207 struct cumulator {
2208 evalue *factor;
2209 evalue *v;
2210 dpoly_r *r;
2212 cumulator(evalue *factor, evalue *v, dpoly_r *r) :
2213 factor(factor), v(v), r(r) {}
2215 void cumulate();
2217 virtual void add_term(int *powers, int len, evalue *f2) = 0;
2220 void cumulator::cumulate()
2222 evalue cum; // factor * 1 * E_num[0]/1 * (E_num[0]-1)/2 *...
2223 evalue f;
2224 evalue t; // E_num[0] - (m-1)
2225 #ifdef USE_MODULO
2226 evalue *cst;
2227 #else
2228 evalue mone;
2229 value_init(mone.d);
2230 evalue_set_si(&mone, -1, 1);
2231 #endif
2233 value_init(cum.d);
2234 evalue_copy(&cum, factor);
2235 value_init(f.d);
2236 value_init(f.x.n);
2237 value_set_si(f.d, 1);
2238 value_set_si(f.x.n, 1);
2239 value_init(t.d);
2240 evalue_copy(&t, v);
2242 #ifdef USE_MODULO
2243 for (cst = &t; value_zero_p(cst->d); ) {
2244 if (cst->x.p->type == fractional)
2245 cst = &cst->x.p->arr[1];
2246 else
2247 cst = &cst->x.p->arr[0];
2249 #endif
2251 for (int m = 0; m < r->len; ++m) {
2252 if (m > 0) {
2253 if (m > 1) {
2254 value_set_si(f.d, m);
2255 emul(&f, &cum);
2256 #ifdef USE_MODULO
2257 value_subtract(cst->x.n, cst->x.n, cst->d);
2258 #else
2259 eadd(&mone, &t);
2260 #endif
2262 emul(&t, &cum);
2264 vector< dpoly_r_term * >& current = r->c[r->len-1-m];
2265 for (int j = 0; j < current.size(); ++j) {
2266 if (current[j]->coeff == 0)
2267 continue;
2268 evalue *f2 = new evalue;
2269 value_init(f2->d);
2270 value_init(f2->x.n);
2271 zz2value(current[j]->coeff, f2->x.n);
2272 zz2value(r->denom, f2->d);
2273 emul(&cum, f2);
2275 add_term(current[j]->powers, r->dim, f2);
2278 free_evalue_refs(&f);
2279 free_evalue_refs(&t);
2280 free_evalue_refs(&cum);
2281 #ifndef USE_MODULO
2282 free_evalue_refs(&mone);
2283 #endif
2286 struct E_poly_term {
2287 int *powers;
2288 evalue *E;
2291 struct ie_cum : public cumulator {
2292 vector<E_poly_term *> terms;
2294 ie_cum(evalue *factor, evalue *v, dpoly_r *r) : cumulator(factor, v, r) {}
2296 virtual void add_term(int *powers, int len, evalue *f2);
2299 void ie_cum::add_term(int *powers, int len, evalue *f2)
2301 int k;
2302 for (k = 0; k < terms.size(); ++k) {
2303 if (memcmp(terms[k]->powers, powers, len * sizeof(int)) == 0) {
2304 eadd(f2, terms[k]->E);
2305 free_evalue_refs(f2);
2306 delete f2;
2307 break;
2310 if (k >= terms.size()) {
2311 E_poly_term *ET = new E_poly_term;
2312 ET->powers = new int[len];
2313 memcpy(ET->powers, powers, len * sizeof(int));
2314 ET->E = f2;
2315 terms.push_back(ET);
2319 struct ienumerator : public polar_decomposer, public vertex_decomposer,
2320 public enumerator_base {
2321 //Polyhedron *pVD;
2322 mat_ZZ den;
2323 vec_ZZ vertex;
2324 mpq_t tcount;
2326 ienumerator(Polyhedron *P, unsigned dim, unsigned nbV) :
2327 vertex_decomposer(P, nbV, this), enumerator_base(dim, this) {
2328 vertex.SetLength(dim);
2330 den.SetDims(dim, dim);
2331 mpq_init(tcount);
2334 ~ienumerator() {
2335 mpq_clear(tcount);
2338 virtual void handle_polar(Polyhedron *P, int sign);
2339 void reduce(evalue *factor, vec_ZZ& num, mat_ZZ& den_f);
2342 void ienumerator::reduce(
2343 evalue *factor, vec_ZZ& num, mat_ZZ& den_f)
2345 unsigned len = den_f.NumRows(); // number of factors in den
2346 unsigned dim = num.length();
2348 if (dim == 0) {
2349 eadd(factor, vE[vert]);
2350 return;
2353 vec_ZZ den_s;
2354 den_s.SetLength(len);
2355 mat_ZZ den_r;
2356 den_r.SetDims(len, dim-1);
2358 int r, k;
2360 for (r = 0; r < len; ++r) {
2361 den_s[r] = den_f[r][0];
2362 for (k = 0; k <= dim-1; ++k)
2363 if (k != 0)
2364 den_r[r][k-(k>0)] = den_f[r][k];
2367 ZZ num_s = num[0];
2368 vec_ZZ num_p;
2369 num_p.SetLength(dim-1);
2370 for (k = 0 ; k <= dim-1; ++k)
2371 if (k != 0)
2372 num_p[k-(k>0)] = num[k];
2374 vec_ZZ den_p;
2375 den_p.SetLength(len);
2377 ZZ one;
2378 one = 1;
2379 normalize(one, num_s, num_p, den_s, den_p, den_r);
2380 if (one != 1)
2381 emul(&mone, factor);
2383 int only_param = 0;
2384 int no_param = 0;
2385 for (int k = 0; k < len; ++k) {
2386 if (den_p[k] == 0)
2387 ++no_param;
2388 else if (den_s[k] == 0)
2389 ++only_param;
2391 if (no_param == 0) {
2392 reduce(factor, num_p, den_r);
2393 } else {
2394 int k, l;
2395 mat_ZZ pden;
2396 pden.SetDims(only_param, dim-1);
2398 for (k = 0, l = 0; k < len; ++k)
2399 if (den_s[k] == 0)
2400 pden[l++] = den_r[k];
2402 for (k = 0; k < len; ++k)
2403 if (den_p[k] == 0)
2404 break;
2406 dpoly n(no_param, num_s);
2407 dpoly D(no_param, den_s[k], 1);
2408 for ( ; ++k < len; )
2409 if (den_p[k] == 0) {
2410 dpoly fact(no_param, den_s[k], 1);
2411 D *= fact;
2414 dpoly_r * r = 0;
2415 // if no_param + only_param == len then all powers
2416 // below will be all zero
2417 if (no_param + only_param == len) {
2418 if (E_num(0, dim) != 0)
2419 r = new dpoly_r(n, len);
2420 else {
2421 mpq_set_si(tcount, 0, 1);
2422 one = 1;
2423 n.div(D, tcount, one);
2425 if (value_notzero_p(mpq_numref(tcount))) {
2426 evalue f;
2427 value_init(f.d);
2428 value_init(f.x.n);
2429 value_assign(f.x.n, mpq_numref(tcount));
2430 value_assign(f.d, mpq_denref(tcount));
2431 emul(&f, factor);
2432 reduce(factor, num_p, pden);
2433 free_evalue_refs(&f);
2435 return;
2437 } else {
2438 for (k = 0; k < len; ++k) {
2439 if (den_s[k] == 0 || den_p[k] == 0)
2440 continue;
2442 dpoly pd(no_param-1, den_s[k], 1);
2444 int l;
2445 for (l = 0; l < k; ++l)
2446 if (den_r[l] == den_r[k])
2447 break;
2449 if (r == 0)
2450 r = new dpoly_r(n, pd, l, len);
2451 else {
2452 dpoly_r *nr = new dpoly_r(r, pd, l, len);
2453 delete r;
2454 r = nr;
2458 dpoly_r *rc = r->div(D);
2459 delete r;
2460 r = rc;
2461 if (E_num(0, dim) == 0) {
2462 int common = pden.NumRows();
2463 vector< dpoly_r_term * >& final = r->c[r->len-1];
2464 int rows;
2465 evalue t;
2466 evalue f;
2467 value_init(f.d);
2468 value_init(f.x.n);
2469 zz2value(r->denom, f.d);
2470 for (int j = 0; j < final.size(); ++j) {
2471 if (final[j]->coeff == 0)
2472 continue;
2473 rows = common;
2474 for (int k = 0; k < r->dim; ++k) {
2475 int n = final[j]->powers[k];
2476 if (n == 0)
2477 continue;
2478 pden.SetDims(rows+n, pden.NumCols());
2479 for (int l = 0; l < n; ++l)
2480 pden[rows+l] = den_r[k];
2481 rows += n;
2483 value_init(t.d);
2484 evalue_copy(&t, factor);
2485 zz2value(final[j]->coeff, f.x.n);
2486 emul(&f, &t);
2487 reduce(&t, num_p, pden);
2488 free_evalue_refs(&t);
2490 free_evalue_refs(&f);
2491 } else {
2492 ie_cum cum(factor, E_num(0, dim), r);
2493 cum.cumulate();
2495 int common = pden.NumRows();
2496 int rows;
2497 for (int j = 0; j < cum.terms.size(); ++j) {
2498 rows = common;
2499 pden.SetDims(rows, pden.NumCols());
2500 for (int k = 0; k < r->dim; ++k) {
2501 int n = cum.terms[j]->powers[k];
2502 if (n == 0)
2503 continue;
2504 pden.SetDims(rows+n, pden.NumCols());
2505 for (int l = 0; l < n; ++l)
2506 pden[rows+l] = den_r[k];
2507 rows += n;
2509 reduce(cum.terms[j]->E, num_p, pden);
2510 free_evalue_refs(cum.terms[j]->E);
2511 delete cum.terms[j]->E;
2512 delete [] cum.terms[j]->powers;
2513 delete cum.terms[j];
2516 delete r;
2520 static int type_offset(enode *p)
2522 return p->type == fractional ? 1 :
2523 p->type == flooring ? 1 : 0;
2526 static int edegree(evalue *e)
2528 int d = 0;
2529 enode *p;
2531 if (value_notzero_p(e->d))
2532 return 0;
2534 p = e->x.p;
2535 int i = type_offset(p);
2536 if (p->size-i-1 > d)
2537 d = p->size - i - 1;
2538 for (; i < p->size; i++) {
2539 int d2 = edegree(&p->arr[i]);
2540 if (d2 > d)
2541 d = d2;
2543 return d;
2546 void ienumerator::handle_polar(Polyhedron *C, int s)
2548 assert(C->NbRays-1 == dim);
2550 lattice_point(V, C, vertex, E_vertex);
2552 int r;
2553 for (r = 0; r < dim; ++r)
2554 values2zz(C->Ray[r]+1, den[r], dim);
2556 evalue one;
2557 value_init(one.d);
2558 evalue_set_si(&one, s, 1);
2559 reduce(&one, vertex, den);
2560 free_evalue_refs(&one);
2562 for (int i = 0; i < dim; ++i)
2563 if (E_vertex[i]) {
2564 free_evalue_refs(E_vertex[i]);
2565 delete E_vertex[i];
2569 struct bfenumerator : public vertex_decomposer, public bf_base,
2570 public enumerator_base {
2571 evalue *factor;
2573 bfenumerator(Polyhedron *P, unsigned dim, unsigned nbV) :
2574 vertex_decomposer(P, nbV, this),
2575 bf_base(P, dim), enumerator_base(dim, this) {
2576 lower = 0;
2577 factor = NULL;
2580 ~bfenumerator() {
2583 virtual void handle_polar(Polyhedron *P, int sign);
2584 virtual void base(mat_ZZ& factors, bfc_vec& v);
2586 bfc_term_base* new_bf_term(int len) {
2587 bfe_term* t = new bfe_term(len);
2588 return t;
2591 virtual void set_factor(bfc_term_base *t, int k, int change) {
2592 bfe_term* bfet = static_cast<bfe_term *>(t);
2593 factor = bfet->factors[k];
2594 assert(factor != NULL);
2595 bfet->factors[k] = NULL;
2596 if (change)
2597 emul(&mone, factor);
2600 virtual void set_factor(bfc_term_base *t, int k, mpq_t &q, int change) {
2601 bfe_term* bfet = static_cast<bfe_term *>(t);
2602 factor = bfet->factors[k];
2603 assert(factor != NULL);
2604 bfet->factors[k] = NULL;
2606 evalue f;
2607 value_init(f.d);
2608 value_init(f.x.n);
2609 if (change)
2610 value_oppose(f.x.n, mpq_numref(q));
2611 else
2612 value_assign(f.x.n, mpq_numref(q));
2613 value_assign(f.d, mpq_denref(q));
2614 emul(&f, factor);
2617 virtual void set_factor(bfc_term_base *t, int k, ZZ& n, ZZ& d, int change) {
2618 bfe_term* bfet = static_cast<bfe_term *>(t);
2620 factor = new evalue;
2622 evalue f;
2623 value_init(f.d);
2624 value_init(f.x.n);
2625 zz2value(n, f.x.n);
2626 if (change)
2627 value_oppose(f.x.n, f.x.n);
2628 zz2value(d, f.d);
2630 value_init(factor->d);
2631 evalue_copy(factor, bfet->factors[k]);
2632 emul(&f, factor);
2635 void set_factor(evalue *f, int change) {
2636 if (change)
2637 emul(&mone, f);
2638 factor = f;
2641 virtual void insert_term(bfc_term_base *t, int i) {
2642 bfe_term* bfet = static_cast<bfe_term *>(t);
2643 int len = t->terms.NumRows()-1; // already increased by one
2645 bfet->factors.resize(len+1);
2646 for (int j = len; j > i; --j) {
2647 bfet->factors[j] = bfet->factors[j-1];
2648 t->terms[j] = t->terms[j-1];
2650 bfet->factors[i] = factor;
2651 factor = NULL;
2654 virtual void update_term(bfc_term_base *t, int i) {
2655 bfe_term* bfet = static_cast<bfe_term *>(t);
2657 eadd(factor, bfet->factors[i]);
2658 free_evalue_refs(factor);
2659 delete factor;
2662 virtual bool constant_vertex(int dim) { return E_num(0, dim) == 0; }
2664 virtual void cum(bf_reducer *bfr, bfc_term_base *t, int k, dpoly_r *r);
2667 struct bfe_cum : public cumulator {
2668 bfenumerator *bfe;
2669 bfc_term_base *told;
2670 int k;
2671 bf_reducer *bfr;
2673 bfe_cum(evalue *factor, evalue *v, dpoly_r *r, bf_reducer *bfr,
2674 bfc_term_base *t, int k, bfenumerator *e) :
2675 cumulator(factor, v, r), told(t), k(k),
2676 bfr(bfr), bfe(e) {
2679 virtual void add_term(int *powers, int len, evalue *f2);
2682 void bfe_cum::add_term(int *powers, int len, evalue *f2)
2684 bfr->update_powers(powers, len);
2686 bfc_term_base * t = bfe->find_bfc_term(bfr->vn, bfr->npowers, bfr->nnf);
2687 bfe->set_factor(f2, bfr->l_changes % 2);
2688 bfe->add_term(t, told->terms[k], bfr->l_extra_num);
2691 void bfenumerator::cum(bf_reducer *bfr, bfc_term_base *t, int k,
2692 dpoly_r *r)
2694 bfe_term* bfet = static_cast<bfe_term *>(t);
2695 bfe_cum cum(bfet->factors[k], E_num(0, bfr->d), r, bfr, t, k, this);
2696 cum.cumulate();
2699 void bfenumerator::base(mat_ZZ& factors, bfc_vec& v)
2701 for (int i = 0; i < v.size(); ++i) {
2702 assert(v[i]->terms.NumRows() == 1);
2703 evalue *factor = static_cast<bfe_term *>(v[i])->factors[0];
2704 eadd(factor, vE[vert]);
2705 delete v[i];
2709 void bfenumerator::handle_polar(Polyhedron *C, int s)
2711 assert(C->NbRays-1 == enumerator_base::dim);
2713 bfe_term* t = new bfe_term(enumerator_base::dim);
2714 vector< bfc_term_base * > v;
2715 v.push_back(t);
2717 t->factors.resize(1);
2719 t->terms.SetDims(1, enumerator_base::dim);
2720 lattice_point(V, C, t->terms[0], E_vertex);
2722 // the elements of factors are always lexpositive
2723 mat_ZZ factors;
2724 s = setup_factors(C, factors, t, s);
2726 t->factors[0] = new evalue;
2727 value_init(t->factors[0]->d);
2728 evalue_set_si(t->factors[0], s, 1);
2729 reduce(factors, v);
2731 for (int i = 0; i < enumerator_base::dim; ++i)
2732 if (E_vertex[i]) {
2733 free_evalue_refs(E_vertex[i]);
2734 delete E_vertex[i];
2738 #ifdef HAVE_CORRECT_VERTICES
2739 static inline Param_Polyhedron *Polyhedron2Param_SD(Polyhedron **Din,
2740 Polyhedron *Cin,int WS,Polyhedron **CEq,Matrix **CT)
2742 if (WS & POL_NO_DUAL)
2743 WS = 0;
2744 return Polyhedron2Param_SimplifiedDomain(Din, Cin, WS, CEq, CT);
2746 #else
2747 static Param_Polyhedron *Polyhedron2Param_SD(Polyhedron **Din,
2748 Polyhedron *Cin,int WS,Polyhedron **CEq,Matrix **CT)
2750 static char data[] = " 1 0 0 0 0 1 -18 "
2751 " 1 0 0 -20 0 19 1 "
2752 " 1 0 1 20 0 -20 16 "
2753 " 1 0 0 0 0 -1 19 "
2754 " 1 0 -1 0 0 0 4 "
2755 " 1 4 -20 0 0 -1 23 "
2756 " 1 -4 20 0 0 1 -22 "
2757 " 1 0 1 0 20 -20 16 "
2758 " 1 0 0 0 -20 19 1 ";
2759 static int checked = 0;
2760 if (!checked) {
2761 checked = 1;
2762 char *p = data;
2763 int n, v, i;
2764 Matrix *M = Matrix_Alloc(9, 7);
2765 for (i = 0; i < 9; ++i)
2766 for (int j = 0; j < 7; ++j) {
2767 sscanf(p, "%d%n", &v, &n);
2768 p += n;
2769 value_set_si(M->p[i][j], v);
2771 Polyhedron *P = Constraints2Polyhedron(M, 1024);
2772 Matrix_Free(M);
2773 Polyhedron *U = Universe_Polyhedron(1);
2774 Param_Polyhedron *PP = Polyhedron2Param_Domain(P, U, 1024);
2775 Polyhedron_Free(P);
2776 Polyhedron_Free(U);
2777 Param_Vertices *V;
2778 for (i = 0, V = PP->V; V; ++i, V = V->next)
2780 if (PP)
2781 Param_Polyhedron_Free(PP);
2782 if (i != 10) {
2783 fprintf(stderr, "WARNING: results may be incorrect\n");
2784 fprintf(stderr,
2785 "WARNING: use latest version of PolyLib to remove this warning\n");
2789 return Polyhedron2Param_SimplifiedDomain(Din, Cin, WS, CEq, CT);
2791 #endif
2793 static evalue* barvinok_enumerate_ev_f(Polyhedron *P, Polyhedron* C,
2794 unsigned MaxRays);
2796 /* Destroys C */
2797 static evalue* barvinok_enumerate_cst(Polyhedron *P, Polyhedron* C,
2798 unsigned MaxRays)
2800 evalue *eres;
2802 ALLOC(evalue, eres);
2803 value_init(eres->d);
2804 value_set_si(eres->d, 0);
2805 eres->x.p = new_enode(partition, 2, C->Dimension);
2806 EVALUE_SET_DOMAIN(eres->x.p->arr[0], DomainConstraintSimplify(C, MaxRays));
2807 value_set_si(eres->x.p->arr[1].d, 1);
2808 value_init(eres->x.p->arr[1].x.n);
2809 if (emptyQ(P))
2810 value_set_si(eres->x.p->arr[1].x.n, 0);
2811 else
2812 barvinok_count(P, &eres->x.p->arr[1].x.n, MaxRays);
2814 return eres;
2817 evalue* barvinok_enumerate_ev(Polyhedron *P, Polyhedron* C, unsigned MaxRays)
2819 //P = unfringe(P, MaxRays);
2820 Polyhedron *Corig = C;
2821 Polyhedron *CEq = NULL, *rVD, *CA;
2822 int r = 0;
2823 unsigned nparam = C->Dimension;
2824 evalue *eres;
2826 evalue factor;
2827 value_init(factor.d);
2828 evalue_set_si(&factor, 1, 1);
2830 CA = align_context(C, P->Dimension, MaxRays);
2831 P = DomainIntersection(P, CA, MaxRays);
2832 Polyhedron_Free(CA);
2834 /* for now */
2835 POL_ENSURE_FACETS(P);
2836 POL_ENSURE_VERTICES(P);
2837 POL_ENSURE_FACETS(C);
2838 POL_ENSURE_VERTICES(C);
2840 if (C->Dimension == 0 || emptyQ(P)) {
2841 constant:
2842 eres = barvinok_enumerate_cst(P, CEq ? CEq : Polyhedron_Copy(C),
2843 MaxRays);
2844 out:
2845 emul(&factor, eres);
2846 reduce_evalue(eres);
2847 free_evalue_refs(&factor);
2848 Domain_Free(P);
2849 if (C != Corig)
2850 Polyhedron_Free(C);
2852 return eres;
2854 if (Polyhedron_is_infinite(P, nparam))
2855 goto constant;
2857 if (P->NbEq != 0) {
2858 Matrix *f;
2859 P = remove_equalities_p(P, P->Dimension-nparam, &f);
2860 mask(f, &factor);
2861 Matrix_Free(f);
2863 if (P->Dimension == nparam) {
2864 CEq = P;
2865 P = Universe_Polyhedron(0);
2866 goto constant;
2869 Polyhedron *T = Polyhedron_Factor(P, nparam, MaxRays);
2870 if (T || (P->Dimension == nparam+1)) {
2871 Polyhedron *Q;
2872 Polyhedron *C2;
2873 for (Q = T ? T : P; Q; Q = Q->next) {
2874 Polyhedron *next = Q->next;
2875 Q->next = NULL;
2877 Polyhedron *QC = Q;
2878 if (Q->Dimension != C->Dimension)
2879 QC = Polyhedron_Project(Q, nparam);
2881 C2 = C;
2882 C = DomainIntersection(C, QC, MaxRays);
2883 if (C2 != Corig)
2884 Polyhedron_Free(C2);
2885 if (QC != Q)
2886 Polyhedron_Free(QC);
2888 Q->next = next;
2891 if (T) {
2892 Polyhedron_Free(P);
2893 P = T;
2894 if (T->Dimension == C->Dimension) {
2895 P = T->next;
2896 T->next = NULL;
2897 Polyhedron_Free(T);
2901 Polyhedron *next = P->next;
2902 P->next = NULL;
2903 eres = barvinok_enumerate_ev_f(P, C, MaxRays);
2904 P->next = next;
2906 if (P->next) {
2907 Polyhedron *Q;
2908 evalue *f;
2910 for (Q = P->next; Q; Q = Q->next) {
2911 Polyhedron *next = Q->next;
2912 Q->next = NULL;
2914 f = barvinok_enumerate_ev_f(Q, C, MaxRays);
2915 emul(f, eres);
2916 free_evalue_refs(f);
2917 free(f);
2919 Q->next = next;
2923 goto out;
2926 static evalue* barvinok_enumerate_ev_f(Polyhedron *P, Polyhedron* C,
2927 unsigned MaxRays)
2929 unsigned nparam = C->Dimension;
2931 if (P->Dimension - nparam == 1)
2932 return ParamLine_Length(P, C, MaxRays);
2934 Param_Polyhedron *PP = NULL;
2935 Polyhedron *CEq = NULL, *pVD;
2936 Matrix *CT = NULL;
2937 Param_Domain *D, *next;
2938 Param_Vertices *V;
2939 evalue *eres;
2940 Polyhedron *Porig = P;
2942 PP = Polyhedron2Param_SD(&P,C,MaxRays,&CEq,&CT);
2944 if (isIdentity(CT)) {
2945 Matrix_Free(CT);
2946 CT = NULL;
2947 } else {
2948 assert(CT->NbRows != CT->NbColumns);
2949 if (CT->NbRows == 1) { // no more parameters
2950 eres = barvinok_enumerate_cst(P, CEq, MaxRays);
2951 out:
2952 if (CT)
2953 Matrix_Free(CT);
2954 if (PP)
2955 Param_Polyhedron_Free(PP);
2956 if (P != Porig)
2957 Polyhedron_Free(P);
2959 return eres;
2961 nparam = CT->NbRows - 1;
2964 unsigned dim = P->Dimension - nparam;
2966 ALLOC(evalue, eres);
2967 value_init(eres->d);
2968 value_set_si(eres->d, 0);
2970 int nd;
2971 for (nd = 0, D=PP->D; D; ++nd, D=D->next);
2972 struct section { Polyhedron *D; evalue E; };
2973 section *s = new section[nd];
2974 Polyhedron **fVD = new Polyhedron_p[nd];
2976 try_again:
2977 #ifdef USE_INCREMENTAL_BF
2978 bfenumerator et(P, dim, PP->nbV);
2979 #elif defined USE_INCREMENTAL_DF
2980 ienumerator et(P, dim, PP->nbV);
2981 #else
2982 enumerator et(P, dim, PP->nbV);
2983 #endif
2985 for(nd = 0, D=PP->D; D; D=next) {
2986 next = D->next;
2988 Polyhedron *rVD = reduce_domain(D->Domain, CT, CEq,
2989 fVD, nd, MaxRays);
2990 if (!rVD)
2991 continue;
2993 pVD = CT ? DomainImage(rVD,CT,MaxRays) : rVD;
2995 value_init(s[nd].E.d);
2996 evalue_set_si(&s[nd].E, 0, 1);
2997 s[nd].D = rVD;
2999 FORALL_PVertex_in_ParamPolyhedron(V,D,PP) // _i is internal counter
3000 if (!et.vE[_i])
3001 try {
3002 et.decompose_at(V, _i, MaxRays);
3003 } catch (OrthogonalException &e) {
3004 if (rVD != pVD)
3005 Domain_Free(pVD);
3006 for (; nd >= 0; --nd) {
3007 free_evalue_refs(&s[nd].E);
3008 Domain_Free(s[nd].D);
3009 Domain_Free(fVD[nd]);
3011 goto try_again;
3013 eadd(et.vE[_i] , &s[nd].E);
3014 END_FORALL_PVertex_in_ParamPolyhedron;
3015 evalue_range_reduction_in_domain(&s[nd].E, pVD);
3017 if (CT)
3018 addeliminatedparams_evalue(&s[nd].E, CT);
3019 ++nd;
3020 if (rVD != pVD)
3021 Domain_Free(pVD);
3024 if (nd == 0)
3025 evalue_set_si(eres, 0, 1);
3026 else {
3027 eres->x.p = new_enode(partition, 2*nd, C->Dimension);
3028 for (int j = 0; j < nd; ++j) {
3029 EVALUE_SET_DOMAIN(eres->x.p->arr[2*j], s[j].D);
3030 value_clear(eres->x.p->arr[2*j+1].d);
3031 eres->x.p->arr[2*j+1] = s[j].E;
3032 Domain_Free(fVD[j]);
3035 delete [] s;
3036 delete [] fVD;
3038 if (CEq)
3039 Polyhedron_Free(CEq);
3040 goto out;
3043 Enumeration* barvinok_enumerate(Polyhedron *P, Polyhedron* C, unsigned MaxRays)
3045 evalue *EP = barvinok_enumerate_ev(P, C, MaxRays);
3047 return partition2enumeration(EP);
3050 static void SwapColumns(Value **V, int n, int i, int j)
3052 for (int r = 0; r < n; ++r)
3053 value_swap(V[r][i], V[r][j]);
3056 static void SwapColumns(Polyhedron *P, int i, int j)
3058 SwapColumns(P->Constraint, P->NbConstraints, i, j);
3059 SwapColumns(P->Ray, P->NbRays, i, j);
3062 /* Construct a constraint c from constraints l and u such that if
3063 * if constraint c holds then for each value of the other variables
3064 * there is at most one value of variable pos (position pos+1 in the constraints).
3066 * Given a lower and an upper bound
3067 * n_l v_i + <c_l,x> + c_l >= 0
3068 * -n_u v_i + <c_u,x> + c_u >= 0
3069 * the constructed constraint is
3071 * -(n_l<c_u,x> + n_u<c_l,x>) + (-n_l c_u - n_u c_l + n_l n_u - 1)
3073 * which is then simplified to remove the content of the non-constant coefficients
3075 * len is the total length of the constraints.
3076 * v is a temporary variable that can be used by this procedure
3078 static void negative_test_constraint(Value *l, Value *u, Value *c, int pos,
3079 int len, Value *v)
3081 value_oppose(*v, u[pos+1]);
3082 Vector_Combine(l+1, u+1, c+1, *v, l[pos+1], len-1);
3083 value_multiply(*v, *v, l[pos+1]);
3084 value_subtract(c[len-1], c[len-1], *v);
3085 value_set_si(*v, -1);
3086 Vector_Scale(c+1, c+1, *v, len-1);
3087 value_decrement(c[len-1], c[len-1]);
3088 ConstraintSimplify(c, c, len, v);
3091 static bool parallel_constraints(Value *l, Value *u, Value *c, int pos,
3092 int len)
3094 bool parallel;
3095 Value g1;
3096 Value g2;
3097 value_init(g1);
3098 value_init(g2);
3100 Vector_Gcd(&l[1+pos], len, &g1);
3101 Vector_Gcd(&u[1+pos], len, &g2);
3102 Vector_Combine(l+1+pos, u+1+pos, c+1, g2, g1, len);
3103 parallel = First_Non_Zero(c+1, len) == -1;
3105 value_clear(g1);
3106 value_clear(g2);
3108 return parallel;
3111 static void negative_test_constraint7(Value *l, Value *u, Value *c, int pos,
3112 int exist, int len, Value *v)
3114 Value g;
3115 value_init(g);
3117 Vector_Gcd(&u[1+pos], exist, v);
3118 Vector_Gcd(&l[1+pos], exist, &g);
3119 Vector_Combine(l+1, u+1, c+1, *v, g, len-1);
3120 value_multiply(*v, *v, g);
3121 value_subtract(c[len-1], c[len-1], *v);
3122 value_set_si(*v, -1);
3123 Vector_Scale(c+1, c+1, *v, len-1);
3124 value_decrement(c[len-1], c[len-1]);
3125 ConstraintSimplify(c, c, len, v);
3127 value_clear(g);
3130 /* Turns a x + b >= 0 into a x + b <= -1
3132 * len is the total length of the constraint.
3133 * v is a temporary variable that can be used by this procedure
3135 static void oppose_constraint(Value *c, int len, Value *v)
3137 value_set_si(*v, -1);
3138 Vector_Scale(c+1, c+1, *v, len-1);
3139 value_decrement(c[len-1], c[len-1]);
3142 /* Split polyhedron P into two polyhedra *pos and *neg, where
3143 * existential variable i has at most one solution for each
3144 * value of the other variables in *neg.
3146 * The splitting is performed using constraints l and u.
3148 * nvar: number of set variables
3149 * row: temporary vector that can be used by this procedure
3150 * f: temporary value that can be used by this procedure
3152 static bool SplitOnConstraint(Polyhedron *P, int i, int l, int u,
3153 int nvar, int MaxRays, Vector *row, Value& f,
3154 Polyhedron **pos, Polyhedron **neg)
3156 negative_test_constraint(P->Constraint[l], P->Constraint[u],
3157 row->p, nvar+i, P->Dimension+2, &f);
3158 *neg = AddConstraints(row->p, 1, P, MaxRays);
3160 /* We found an independent, but useless constraint
3161 * Maybe we should detect this earlier and not
3162 * mark the variable as INDEPENDENT
3164 if (emptyQ((*neg))) {
3165 Polyhedron_Free(*neg);
3166 return false;
3169 oppose_constraint(row->p, P->Dimension+2, &f);
3170 *pos = AddConstraints(row->p, 1, P, MaxRays);
3172 if (emptyQ((*pos))) {
3173 Polyhedron_Free(*neg);
3174 Polyhedron_Free(*pos);
3175 return false;
3178 return true;
3182 * unimodularly transform P such that constraint r is transformed
3183 * into a constraint that involves only a single (the first)
3184 * existential variable
3187 static Polyhedron *rotate_along(Polyhedron *P, int r, int nvar, int exist,
3188 unsigned MaxRays)
3190 Value g;
3191 value_init(g);
3193 Vector *row = Vector_Alloc(exist);
3194 Vector_Copy(P->Constraint[r]+1+nvar, row->p, exist);
3195 Vector_Gcd(row->p, exist, &g);
3196 if (value_notone_p(g))
3197 Vector_AntiScale(row->p, row->p, g, exist);
3198 value_clear(g);
3200 Matrix *M = unimodular_complete(row);
3201 Matrix *M2 = Matrix_Alloc(P->Dimension+1, P->Dimension+1);
3202 for (r = 0; r < nvar; ++r)
3203 value_set_si(M2->p[r][r], 1);
3204 for ( ; r < nvar+exist; ++r)
3205 Vector_Copy(M->p[r-nvar], M2->p[r]+nvar, exist);
3206 for ( ; r < P->Dimension+1; ++r)
3207 value_set_si(M2->p[r][r], 1);
3208 Polyhedron *T = Polyhedron_Image(P, M2, MaxRays);
3210 Matrix_Free(M2);
3211 Matrix_Free(M);
3212 Vector_Free(row);
3214 return T;
3217 /* Split polyhedron P into two polyhedra *pos and *neg, where
3218 * existential variable i has at most one solution for each
3219 * value of the other variables in *neg.
3221 * If independent is set, then the two constraints on which the
3222 * split will be performed need to be independent of the other
3223 * existential variables.
3225 * Return true if an appropriate split could be performed.
3227 * nvar: number of set variables
3228 * exist: number of existential variables
3229 * row: temporary vector that can be used by this procedure
3230 * f: temporary value that can be used by this procedure
3232 static bool SplitOnVar(Polyhedron *P, int i,
3233 int nvar, int exist, int MaxRays,
3234 Vector *row, Value& f, bool independent,
3235 Polyhedron **pos, Polyhedron **neg)
3237 int j;
3239 for (int l = P->NbEq; l < P->NbConstraints; ++l) {
3240 if (value_negz_p(P->Constraint[l][nvar+i+1]))
3241 continue;
3243 if (independent) {
3244 for (j = 0; j < exist; ++j)
3245 if (j != i && value_notzero_p(P->Constraint[l][nvar+j+1]))
3246 break;
3247 if (j < exist)
3248 continue;
3251 for (int u = P->NbEq; u < P->NbConstraints; ++u) {
3252 if (value_posz_p(P->Constraint[u][nvar+i+1]))
3253 continue;
3255 if (independent) {
3256 for (j = 0; j < exist; ++j)
3257 if (j != i && value_notzero_p(P->Constraint[u][nvar+j+1]))
3258 break;
3259 if (j < exist)
3260 continue;
3263 if (SplitOnConstraint(P, i, l, u, nvar, MaxRays, row, f, pos, neg)) {
3264 if (independent) {
3265 if (i != 0)
3266 SwapColumns(*neg, nvar+1, nvar+1+i);
3268 return true;
3273 return false;
3276 static bool double_bound_pair(Polyhedron *P, int nvar, int exist,
3277 int i, int l1, int l2,
3278 Polyhedron **pos, Polyhedron **neg)
3280 Value f;
3281 value_init(f);
3282 Vector *row = Vector_Alloc(P->Dimension+2);
3283 value_set_si(row->p[0], 1);
3284 value_oppose(f, P->Constraint[l1][nvar+i+1]);
3285 Vector_Combine(P->Constraint[l1]+1, P->Constraint[l2]+1,
3286 row->p+1,
3287 P->Constraint[l2][nvar+i+1], f,
3288 P->Dimension+1);
3289 ConstraintSimplify(row->p, row->p, P->Dimension+2, &f);
3290 *pos = AddConstraints(row->p, 1, P, 0);
3291 value_set_si(f, -1);
3292 Vector_Scale(row->p+1, row->p+1, f, P->Dimension+1);
3293 value_decrement(row->p[P->Dimension+1], row->p[P->Dimension+1]);
3294 *neg = AddConstraints(row->p, 1, P, 0);
3295 Vector_Free(row);
3296 value_clear(f);
3298 return !emptyQ((*pos)) && !emptyQ((*neg));
3301 static bool double_bound(Polyhedron *P, int nvar, int exist,
3302 Polyhedron **pos, Polyhedron **neg)
3304 for (int i = 0; i < exist; ++i) {
3305 int l1, l2;
3306 for (l1 = P->NbEq; l1 < P->NbConstraints; ++l1) {
3307 if (value_negz_p(P->Constraint[l1][nvar+i+1]))
3308 continue;
3309 for (l2 = l1 + 1; l2 < P->NbConstraints; ++l2) {
3310 if (value_negz_p(P->Constraint[l2][nvar+i+1]))
3311 continue;
3312 if (double_bound_pair(P, nvar, exist, i, l1, l2, pos, neg))
3313 return true;
3316 for (l1 = P->NbEq; l1 < P->NbConstraints; ++l1) {
3317 if (value_posz_p(P->Constraint[l1][nvar+i+1]))
3318 continue;
3319 if (l1 < P->NbConstraints)
3320 for (l2 = l1 + 1; l2 < P->NbConstraints; ++l2) {
3321 if (value_posz_p(P->Constraint[l2][nvar+i+1]))
3322 continue;
3323 if (double_bound_pair(P, nvar, exist, i, l1, l2, pos, neg))
3324 return true;
3327 return false;
3329 return false;
3332 enum constraint {
3333 ALL_POS = 1 << 0,
3334 ONE_NEG = 1 << 1,
3335 INDEPENDENT = 1 << 2,
3336 ROT_NEG = 1 << 3
3339 static evalue* enumerate_or(Polyhedron *D,
3340 unsigned exist, unsigned nparam, unsigned MaxRays)
3342 #ifdef DEBUG_ER
3343 fprintf(stderr, "\nER: Or\n");
3344 #endif /* DEBUG_ER */
3346 Polyhedron *N = D->next;
3347 D->next = 0;
3348 evalue *EP =
3349 barvinok_enumerate_e(D, exist, nparam, MaxRays);
3350 Polyhedron_Free(D);
3352 for (D = N; D; D = N) {
3353 N = D->next;
3354 D->next = 0;
3356 evalue *EN =
3357 barvinok_enumerate_e(D, exist, nparam, MaxRays);
3359 eor(EN, EP);
3360 free_evalue_refs(EN);
3361 free(EN);
3362 Polyhedron_Free(D);
3365 reduce_evalue(EP);
3367 return EP;
3370 static evalue* enumerate_sum(Polyhedron *P,
3371 unsigned exist, unsigned nparam, unsigned MaxRays)
3373 int nvar = P->Dimension - exist - nparam;
3374 int toswap = nvar < exist ? nvar : exist;
3375 for (int i = 0; i < toswap; ++i)
3376 SwapColumns(P, 1 + i, nvar+exist - i);
3377 nparam += nvar;
3379 #ifdef DEBUG_ER
3380 fprintf(stderr, "\nER: Sum\n");
3381 #endif /* DEBUG_ER */
3383 evalue *EP = barvinok_enumerate_e(P, exist, nparam, MaxRays);
3385 for (int i = 0; i < /* nvar */ nparam; ++i) {
3386 Matrix *C = Matrix_Alloc(1, 1 + nparam + 1);
3387 value_set_si(C->p[0][0], 1);
3388 evalue split;
3389 value_init(split.d);
3390 value_set_si(split.d, 0);
3391 split.x.p = new_enode(partition, 4, nparam);
3392 value_set_si(C->p[0][1+i], 1);
3393 Matrix *C2 = Matrix_Copy(C);
3394 EVALUE_SET_DOMAIN(split.x.p->arr[0],
3395 Constraints2Polyhedron(C2, MaxRays));
3396 Matrix_Free(C2);
3397 evalue_set_si(&split.x.p->arr[1], 1, 1);
3398 value_set_si(C->p[0][1+i], -1);
3399 value_set_si(C->p[0][1+nparam], -1);
3400 EVALUE_SET_DOMAIN(split.x.p->arr[2],
3401 Constraints2Polyhedron(C, MaxRays));
3402 evalue_set_si(&split.x.p->arr[3], 1, 1);
3403 emul(&split, EP);
3404 free_evalue_refs(&split);
3405 Matrix_Free(C);
3407 reduce_evalue(EP);
3408 evalue_range_reduction(EP);
3410 evalue_frac2floor(EP);
3412 evalue *sum = esum(EP, nvar);
3414 free_evalue_refs(EP);
3415 free(EP);
3416 EP = sum;
3418 evalue_range_reduction(EP);
3420 return EP;
3423 static evalue* split_sure(Polyhedron *P, Polyhedron *S,
3424 unsigned exist, unsigned nparam, unsigned MaxRays)
3426 int nvar = P->Dimension - exist - nparam;
3428 Matrix *M = Matrix_Alloc(exist, S->Dimension+2);
3429 for (int i = 0; i < exist; ++i)
3430 value_set_si(M->p[i][nvar+i+1], 1);
3431 Polyhedron *O = S;
3432 S = DomainAddRays(S, M, MaxRays);
3433 Polyhedron_Free(O);
3434 Polyhedron *F = DomainAddRays(P, M, MaxRays);
3435 Polyhedron *D = DomainDifference(F, S, MaxRays);
3436 O = D;
3437 D = Disjoint_Domain(D, 0, MaxRays);
3438 Polyhedron_Free(F);
3439 Domain_Free(O);
3440 Matrix_Free(M);
3442 M = Matrix_Alloc(P->Dimension+1-exist, P->Dimension+1);
3443 for (int j = 0; j < nvar; ++j)
3444 value_set_si(M->p[j][j], 1);
3445 for (int j = 0; j < nparam+1; ++j)
3446 value_set_si(M->p[nvar+j][nvar+exist+j], 1);
3447 Polyhedron *T = Polyhedron_Image(S, M, MaxRays);
3448 evalue *EP = barvinok_enumerate_e(T, 0, nparam, MaxRays);
3449 Polyhedron_Free(S);
3450 Polyhedron_Free(T);
3451 Matrix_Free(M);
3453 for (Polyhedron *Q = D; Q; Q = Q->next) {
3454 Polyhedron *N = Q->next;
3455 Q->next = 0;
3456 T = DomainIntersection(P, Q, MaxRays);
3457 evalue *E = barvinok_enumerate_e(T, exist, nparam, MaxRays);
3458 eadd(E, EP);
3459 free_evalue_refs(E);
3460 free(E);
3461 Polyhedron_Free(T);
3462 Q->next = N;
3464 Domain_Free(D);
3465 return EP;
3468 static evalue* enumerate_sure(Polyhedron *P,
3469 unsigned exist, unsigned nparam, unsigned MaxRays)
3471 int i;
3472 Polyhedron *S = P;
3473 int nvar = P->Dimension - exist - nparam;
3474 Value lcm;
3475 Value f;
3476 value_init(lcm);
3477 value_init(f);
3479 for (i = 0; i < exist; ++i) {
3480 Matrix *M = Matrix_Alloc(S->NbConstraints, S->Dimension+2);
3481 int c = 0;
3482 value_set_si(lcm, 1);
3483 for (int j = 0; j < S->NbConstraints; ++j) {
3484 if (value_negz_p(S->Constraint[j][1+nvar+i]))
3485 continue;
3486 if (value_one_p(S->Constraint[j][1+nvar+i]))
3487 continue;
3488 value_lcm(lcm, S->Constraint[j][1+nvar+i], &lcm);
3491 for (int j = 0; j < S->NbConstraints; ++j) {
3492 if (value_negz_p(S->Constraint[j][1+nvar+i]))
3493 continue;
3494 if (value_one_p(S->Constraint[j][1+nvar+i]))
3495 continue;
3496 value_division(f, lcm, S->Constraint[j][1+nvar+i]);
3497 Vector_Scale(S->Constraint[j], M->p[c], f, S->Dimension+2);
3498 value_subtract(M->p[c][S->Dimension+1],
3499 M->p[c][S->Dimension+1],
3500 lcm);
3501 value_increment(M->p[c][S->Dimension+1],
3502 M->p[c][S->Dimension+1]);
3503 ++c;
3505 Polyhedron *O = S;
3506 S = AddConstraints(M->p[0], c, S, MaxRays);
3507 if (O != P)
3508 Polyhedron_Free(O);
3509 Matrix_Free(M);
3510 if (emptyQ(S)) {
3511 Polyhedron_Free(S);
3512 value_clear(lcm);
3513 value_clear(f);
3514 return 0;
3517 value_clear(lcm);
3518 value_clear(f);
3520 #ifdef DEBUG_ER
3521 fprintf(stderr, "\nER: Sure\n");
3522 #endif /* DEBUG_ER */
3524 return split_sure(P, S, exist, nparam, MaxRays);
3527 static evalue* enumerate_sure2(Polyhedron *P,
3528 unsigned exist, unsigned nparam, unsigned MaxRays)
3530 int nvar = P->Dimension - exist - nparam;
3531 int r;
3532 for (r = 0; r < P->NbRays; ++r)
3533 if (value_one_p(P->Ray[r][0]) &&
3534 value_one_p(P->Ray[r][P->Dimension+1]))
3535 break;
3537 if (r >= P->NbRays)
3538 return 0;
3540 Matrix *M = Matrix_Alloc(nvar + 1 + nparam, P->Dimension+2);
3541 for (int i = 0; i < nvar; ++i)
3542 value_set_si(M->p[i][1+i], 1);
3543 for (int i = 0; i < nparam; ++i)
3544 value_set_si(M->p[i+nvar][1+nvar+exist+i], 1);
3545 Vector_Copy(P->Ray[r]+1+nvar, M->p[nvar+nparam]+1+nvar, exist);
3546 value_set_si(M->p[nvar+nparam][0], 1);
3547 value_set_si(M->p[nvar+nparam][P->Dimension+1], 1);
3548 Polyhedron * F = Rays2Polyhedron(M, MaxRays);
3549 Matrix_Free(M);
3551 Polyhedron *I = DomainIntersection(F, P, MaxRays);
3552 Polyhedron_Free(F);
3554 #ifdef DEBUG_ER
3555 fprintf(stderr, "\nER: Sure2\n");
3556 #endif /* DEBUG_ER */
3558 return split_sure(P, I, exist, nparam, MaxRays);
3561 static evalue* enumerate_cyclic(Polyhedron *P,
3562 unsigned exist, unsigned nparam,
3563 evalue * EP, int r, int p, unsigned MaxRays)
3565 int nvar = P->Dimension - exist - nparam;
3567 /* If EP in its fractional maps only contains references
3568 * to the remainder parameter with appropriate coefficients
3569 * then we could in principle avoid adding existentially
3570 * quantified variables to the validity domains.
3571 * We'd have to replace the remainder by m { p/m }
3572 * and multiply with an appropriate factor that is one
3573 * only in the appropriate range.
3574 * This last multiplication can be avoided if EP
3575 * has a single validity domain with no (further)
3576 * constraints on the remainder parameter
3579 Matrix *CT = Matrix_Alloc(nparam+1, nparam+3);
3580 Matrix *M = Matrix_Alloc(1, 1+nparam+3);
3581 for (int j = 0; j < nparam; ++j)
3582 if (j != p)
3583 value_set_si(CT->p[j][j], 1);
3584 value_set_si(CT->p[p][nparam+1], 1);
3585 value_set_si(CT->p[nparam][nparam+2], 1);
3586 value_set_si(M->p[0][1+p], -1);
3587 value_absolute(M->p[0][1+nparam], P->Ray[0][1+nvar+exist+p]);
3588 value_set_si(M->p[0][1+nparam+1], 1);
3589 Polyhedron *CEq = Constraints2Polyhedron(M, 1);
3590 Matrix_Free(M);
3591 addeliminatedparams_enum(EP, CT, CEq, MaxRays, nparam);
3592 Polyhedron_Free(CEq);
3593 Matrix_Free(CT);
3595 return EP;
3598 static void enumerate_vd_add_ray(evalue *EP, Matrix *Rays, unsigned MaxRays)
3600 if (value_notzero_p(EP->d))
3601 return;
3603 assert(EP->x.p->type == partition);
3604 assert(EP->x.p->pos == EVALUE_DOMAIN(EP->x.p->arr[0])->Dimension);
3605 for (int i = 0; i < EP->x.p->size/2; ++i) {
3606 Polyhedron *D = EVALUE_DOMAIN(EP->x.p->arr[2*i]);
3607 Polyhedron *N = DomainAddRays(D, Rays, MaxRays);
3608 EVALUE_SET_DOMAIN(EP->x.p->arr[2*i], N);
3609 Domain_Free(D);
3613 static evalue* enumerate_line(Polyhedron *P,
3614 unsigned exist, unsigned nparam, unsigned MaxRays)
3616 if (P->NbBid == 0)
3617 return 0;
3619 #ifdef DEBUG_ER
3620 fprintf(stderr, "\nER: Line\n");
3621 #endif /* DEBUG_ER */
3623 int nvar = P->Dimension - exist - nparam;
3624 int i, j;
3625 for (i = 0; i < nparam; ++i)
3626 if (value_notzero_p(P->Ray[0][1+nvar+exist+i]))
3627 break;
3628 assert(i < nparam);
3629 for (j = i+1; j < nparam; ++j)
3630 if (value_notzero_p(P->Ray[0][1+nvar+exist+i]))
3631 break;
3632 assert(j >= nparam); // for now
3634 Matrix *M = Matrix_Alloc(2, P->Dimension+2);
3635 value_set_si(M->p[0][0], 1);
3636 value_set_si(M->p[0][1+nvar+exist+i], 1);
3637 value_set_si(M->p[1][0], 1);
3638 value_set_si(M->p[1][1+nvar+exist+i], -1);
3639 value_absolute(M->p[1][1+P->Dimension], P->Ray[0][1+nvar+exist+i]);
3640 value_decrement(M->p[1][1+P->Dimension], M->p[1][1+P->Dimension]);
3641 Polyhedron *S = AddConstraints(M->p[0], 2, P, MaxRays);
3642 evalue *EP = barvinok_enumerate_e(S, exist, nparam, MaxRays);
3643 Polyhedron_Free(S);
3644 Matrix_Free(M);
3646 return enumerate_cyclic(P, exist, nparam, EP, 0, i, MaxRays);
3649 static int single_param_pos(Polyhedron*P, unsigned exist, unsigned nparam,
3650 int r)
3652 int nvar = P->Dimension - exist - nparam;
3653 if (First_Non_Zero(P->Ray[r]+1, nvar) != -1)
3654 return -1;
3655 int i = First_Non_Zero(P->Ray[r]+1+nvar+exist, nparam);
3656 if (i == -1)
3657 return -1;
3658 if (First_Non_Zero(P->Ray[r]+1+nvar+exist+1, nparam-i-1) != -1)
3659 return -1;
3660 return i;
3663 static evalue* enumerate_remove_ray(Polyhedron *P, int r,
3664 unsigned exist, unsigned nparam, unsigned MaxRays)
3666 #ifdef DEBUG_ER
3667 fprintf(stderr, "\nER: RedundantRay\n");
3668 #endif /* DEBUG_ER */
3670 Value one;
3671 value_init(one);
3672 value_set_si(one, 1);
3673 int len = P->NbRays-1;
3674 Matrix *M = Matrix_Alloc(2 * len, P->Dimension+2);
3675 Vector_Copy(P->Ray[0], M->p[0], r * (P->Dimension+2));
3676 Vector_Copy(P->Ray[r+1], M->p[r], (len-r) * (P->Dimension+2));
3677 for (int j = 0; j < P->NbRays; ++j) {
3678 if (j == r)
3679 continue;
3680 Vector_Combine(P->Ray[j], P->Ray[r], M->p[len+j-(j>r)],
3681 one, P->Ray[j][P->Dimension+1], P->Dimension+2);
3684 P = Rays2Polyhedron(M, MaxRays);
3685 Matrix_Free(M);
3686 evalue *EP = barvinok_enumerate_e(P, exist, nparam, MaxRays);
3687 Polyhedron_Free(P);
3688 value_clear(one);
3690 return EP;
3693 static evalue* enumerate_redundant_ray(Polyhedron *P,
3694 unsigned exist, unsigned nparam, unsigned MaxRays)
3696 assert(P->NbBid == 0);
3697 int nvar = P->Dimension - exist - nparam;
3698 Value m;
3699 value_init(m);
3701 for (int r = 0; r < P->NbRays; ++r) {
3702 if (value_notzero_p(P->Ray[r][P->Dimension+1]))
3703 continue;
3704 int i1 = single_param_pos(P, exist, nparam, r);
3705 if (i1 == -1)
3706 continue;
3707 for (int r2 = r+1; r2 < P->NbRays; ++r2) {
3708 if (value_notzero_p(P->Ray[r2][P->Dimension+1]))
3709 continue;
3710 int i2 = single_param_pos(P, exist, nparam, r2);
3711 if (i2 == -1)
3712 continue;
3713 if (i1 != i2)
3714 continue;
3716 value_division(m, P->Ray[r][1+nvar+exist+i1],
3717 P->Ray[r2][1+nvar+exist+i1]);
3718 value_multiply(m, m, P->Ray[r2][1+nvar+exist+i1]);
3719 /* r2 divides r => r redundant */
3720 if (value_eq(m, P->Ray[r][1+nvar+exist+i1])) {
3721 value_clear(m);
3722 return enumerate_remove_ray(P, r, exist, nparam, MaxRays);
3725 value_division(m, P->Ray[r2][1+nvar+exist+i1],
3726 P->Ray[r][1+nvar+exist+i1]);
3727 value_multiply(m, m, P->Ray[r][1+nvar+exist+i1]);
3728 /* r divides r2 => r2 redundant */
3729 if (value_eq(m, P->Ray[r2][1+nvar+exist+i1])) {
3730 value_clear(m);
3731 return enumerate_remove_ray(P, r2, exist, nparam, MaxRays);
3735 value_clear(m);
3736 return 0;
3739 static Polyhedron *upper_bound(Polyhedron *P,
3740 int pos, Value *max, Polyhedron **R)
3742 Value v;
3743 int r;
3744 value_init(v);
3746 *R = 0;
3747 Polyhedron *N;
3748 Polyhedron *B = 0;
3749 for (Polyhedron *Q = P; Q; Q = N) {
3750 N = Q->next;
3751 for (r = 0; r < P->NbRays; ++r) {
3752 if (value_zero_p(P->Ray[r][P->Dimension+1]) &&
3753 value_pos_p(P->Ray[r][1+pos]))
3754 break;
3756 if (r < P->NbRays) {
3757 Q->next = *R;
3758 *R = Q;
3759 continue;
3760 } else {
3761 Q->next = B;
3762 B = Q;
3764 for (r = 0; r < P->NbRays; ++r) {
3765 if (value_zero_p(P->Ray[r][P->Dimension+1]))
3766 continue;
3767 mpz_fdiv_q(v, P->Ray[r][1+pos], P->Ray[r][1+P->Dimension]);
3768 if ((!Q->next && r == 0) || value_gt(v, *max))
3769 value_assign(*max, v);
3772 value_clear(v);
3773 return B;
3776 static evalue* enumerate_ray(Polyhedron *P,
3777 unsigned exist, unsigned nparam, unsigned MaxRays)
3779 assert(P->NbBid == 0);
3780 int nvar = P->Dimension - exist - nparam;
3782 int r;
3783 for (r = 0; r < P->NbRays; ++r)
3784 if (value_zero_p(P->Ray[r][P->Dimension+1]))
3785 break;
3786 if (r >= P->NbRays)
3787 return 0;
3789 int r2;
3790 for (r2 = r+1; r2 < P->NbRays; ++r2)
3791 if (value_zero_p(P->Ray[r2][P->Dimension+1]))
3792 break;
3793 if (r2 < P->NbRays) {
3794 if (nvar > 0)
3795 return enumerate_sum(P, exist, nparam, MaxRays);
3798 #ifdef DEBUG_ER
3799 fprintf(stderr, "\nER: Ray\n");
3800 #endif /* DEBUG_ER */
3802 Value m;
3803 Value one;
3804 value_init(m);
3805 value_init(one);
3806 value_set_si(one, 1);
3807 int i = single_param_pos(P, exist, nparam, r);
3808 assert(i != -1); // for now;
3810 Matrix *M = Matrix_Alloc(P->NbRays, P->Dimension+2);
3811 for (int j = 0; j < P->NbRays; ++j) {
3812 Vector_Combine(P->Ray[j], P->Ray[r], M->p[j],
3813 one, P->Ray[j][P->Dimension+1], P->Dimension+2);
3815 Polyhedron *S = Rays2Polyhedron(M, MaxRays);
3816 Matrix_Free(M);
3817 Polyhedron *D = DomainDifference(P, S, MaxRays);
3818 Polyhedron_Free(S);
3819 // Polyhedron_Print(stderr, P_VALUE_FMT, D);
3820 assert(value_pos_p(P->Ray[r][1+nvar+exist+i])); // for now
3821 Polyhedron *R;
3822 D = upper_bound(D, nvar+exist+i, &m, &R);
3823 assert(D);
3824 Domain_Free(D);
3826 M = Matrix_Alloc(2, P->Dimension+2);
3827 value_set_si(M->p[0][0], 1);
3828 value_set_si(M->p[1][0], 1);
3829 value_set_si(M->p[0][1+nvar+exist+i], -1);
3830 value_set_si(M->p[1][1+nvar+exist+i], 1);
3831 value_assign(M->p[0][1+P->Dimension], m);
3832 value_oppose(M->p[1][1+P->Dimension], m);
3833 value_addto(M->p[1][1+P->Dimension], M->p[1][1+P->Dimension],
3834 P->Ray[r][1+nvar+exist+i]);
3835 value_decrement(M->p[1][1+P->Dimension], M->p[1][1+P->Dimension]);
3836 // Matrix_Print(stderr, P_VALUE_FMT, M);
3837 D = AddConstraints(M->p[0], 2, P, MaxRays);
3838 // Polyhedron_Print(stderr, P_VALUE_FMT, D);
3839 value_subtract(M->p[0][1+P->Dimension], M->p[0][1+P->Dimension],
3840 P->Ray[r][1+nvar+exist+i]);
3841 // Matrix_Print(stderr, P_VALUE_FMT, M);
3842 S = AddConstraints(M->p[0], 1, P, MaxRays);
3843 // Polyhedron_Print(stderr, P_VALUE_FMT, S);
3844 Matrix_Free(M);
3846 evalue *EP = barvinok_enumerate_e(D, exist, nparam, MaxRays);
3847 Polyhedron_Free(D);
3848 value_clear(one);
3849 value_clear(m);
3851 if (value_notone_p(P->Ray[r][1+nvar+exist+i]))
3852 EP = enumerate_cyclic(P, exist, nparam, EP, r, i, MaxRays);
3853 else {
3854 M = Matrix_Alloc(1, nparam+2);
3855 value_set_si(M->p[0][0], 1);
3856 value_set_si(M->p[0][1+i], 1);
3857 enumerate_vd_add_ray(EP, M, MaxRays);
3858 Matrix_Free(M);
3861 if (!emptyQ(S)) {
3862 evalue *E = barvinok_enumerate_e(S, exist, nparam, MaxRays);
3863 eadd(E, EP);
3864 free_evalue_refs(E);
3865 free(E);
3867 Polyhedron_Free(S);
3869 if (R) {
3870 assert(nvar == 0);
3871 evalue *ER = enumerate_or(R, exist, nparam, MaxRays);
3872 eor(ER, EP);
3873 free_evalue_refs(ER);
3874 free(ER);
3877 return EP;
3880 static evalue* enumerate_vd(Polyhedron **PA,
3881 unsigned exist, unsigned nparam, unsigned MaxRays)
3883 Polyhedron *P = *PA;
3884 int nvar = P->Dimension - exist - nparam;
3885 Param_Polyhedron *PP = NULL;
3886 Polyhedron *C = Universe_Polyhedron(nparam);
3887 Polyhedron *CEq;
3888 Matrix *CT;
3889 Polyhedron *PR = P;
3890 PP = Polyhedron2Param_SimplifiedDomain(&PR,C,MaxRays,&CEq,&CT);
3891 Polyhedron_Free(C);
3893 int nd;
3894 Param_Domain *D, *last;
3895 Value c;
3896 value_init(c);
3897 for (nd = 0, D=PP->D; D; D=D->next, ++nd)
3900 Polyhedron **VD = new Polyhedron_p[nd];
3901 Polyhedron **fVD = new Polyhedron_p[nd];
3902 for(nd = 0, D=PP->D; D; D=D->next) {
3903 Polyhedron *rVD = reduce_domain(D->Domain, CT, CEq,
3904 fVD, nd, MaxRays);
3905 if (!rVD)
3906 continue;
3908 VD[nd++] = rVD;
3909 last = D;
3912 evalue *EP = 0;
3914 if (nd == 0)
3915 EP = evalue_zero();
3917 /* This doesn't seem to have any effect */
3918 if (nd == 1) {
3919 Polyhedron *CA = align_context(VD[0], P->Dimension, MaxRays);
3920 Polyhedron *O = P;
3921 P = DomainIntersection(P, CA, MaxRays);
3922 if (O != *PA)
3923 Polyhedron_Free(O);
3924 Polyhedron_Free(CA);
3925 if (emptyQ(P))
3926 EP = evalue_zero();
3929 if (!EP && CT->NbColumns != CT->NbRows) {
3930 Polyhedron *CEqr = DomainImage(CEq, CT, MaxRays);
3931 Polyhedron *CA = align_context(CEqr, PR->Dimension, MaxRays);
3932 Polyhedron *I = DomainIntersection(PR, CA, MaxRays);
3933 Polyhedron_Free(CEqr);
3934 Polyhedron_Free(CA);
3935 #ifdef DEBUG_ER
3936 fprintf(stderr, "\nER: Eliminate\n");
3937 #endif /* DEBUG_ER */
3938 nparam -= CT->NbColumns - CT->NbRows;
3939 EP = barvinok_enumerate_e(I, exist, nparam, MaxRays);
3940 nparam += CT->NbColumns - CT->NbRows;
3941 addeliminatedparams_enum(EP, CT, CEq, MaxRays, nparam);
3942 Polyhedron_Free(I);
3944 if (PR != *PA)
3945 Polyhedron_Free(PR);
3946 PR = 0;
3948 if (!EP && nd > 1) {
3949 #ifdef DEBUG_ER
3950 fprintf(stderr, "\nER: VD\n");
3951 #endif /* DEBUG_ER */
3952 for (int i = 0; i < nd; ++i) {
3953 Polyhedron *CA = align_context(VD[i], P->Dimension, MaxRays);
3954 Polyhedron *I = DomainIntersection(P, CA, MaxRays);
3956 if (i == 0)
3957 EP = barvinok_enumerate_e(I, exist, nparam, MaxRays);
3958 else {
3959 evalue *E = barvinok_enumerate_e(I, exist, nparam, MaxRays);
3960 eadd(E, EP);
3961 free_evalue_refs(E);
3962 free(E);
3964 Polyhedron_Free(I);
3965 Polyhedron_Free(CA);
3969 for (int i = 0; i < nd; ++i) {
3970 Polyhedron_Free(VD[i]);
3971 Polyhedron_Free(fVD[i]);
3973 delete [] VD;
3974 delete [] fVD;
3975 value_clear(c);
3977 if (!EP && nvar == 0) {
3978 Value f;
3979 value_init(f);
3980 Param_Vertices *V, *V2;
3981 Matrix* M = Matrix_Alloc(1, P->Dimension+2);
3983 FORALL_PVertex_in_ParamPolyhedron(V, last, PP) {
3984 bool found = false;
3985 FORALL_PVertex_in_ParamPolyhedron(V2, last, PP) {
3986 if (V == V2) {
3987 found = true;
3988 continue;
3990 if (!found)
3991 continue;
3992 for (int i = 0; i < exist; ++i) {
3993 value_oppose(f, V->Vertex->p[i][nparam+1]);
3994 Vector_Combine(V->Vertex->p[i],
3995 V2->Vertex->p[i],
3996 M->p[0] + 1 + nvar + exist,
3997 V2->Vertex->p[i][nparam+1],
3999 nparam+1);
4000 int j;
4001 for (j = 0; j < nparam; ++j)
4002 if (value_notzero_p(M->p[0][1+nvar+exist+j]))
4003 break;
4004 if (j >= nparam)
4005 continue;
4006 ConstraintSimplify(M->p[0], M->p[0],
4007 P->Dimension+2, &f);
4008 value_set_si(M->p[0][0], 0);
4009 Polyhedron *para = AddConstraints(M->p[0], 1, P,
4010 MaxRays);
4011 if (emptyQ(para)) {
4012 Polyhedron_Free(para);
4013 continue;
4015 Polyhedron *pos, *neg;
4016 value_set_si(M->p[0][0], 1);
4017 value_decrement(M->p[0][P->Dimension+1],
4018 M->p[0][P->Dimension+1]);
4019 neg = AddConstraints(M->p[0], 1, P, MaxRays);
4020 value_set_si(f, -1);
4021 Vector_Scale(M->p[0]+1, M->p[0]+1, f,
4022 P->Dimension+1);
4023 value_decrement(M->p[0][P->Dimension+1],
4024 M->p[0][P->Dimension+1]);
4025 value_decrement(M->p[0][P->Dimension+1],
4026 M->p[0][P->Dimension+1]);
4027 pos = AddConstraints(M->p[0], 1, P, MaxRays);
4028 if (emptyQ(neg) && emptyQ(pos)) {
4029 Polyhedron_Free(para);
4030 Polyhedron_Free(pos);
4031 Polyhedron_Free(neg);
4032 continue;
4034 #ifdef DEBUG_ER
4035 fprintf(stderr, "\nER: Order\n");
4036 #endif /* DEBUG_ER */
4037 EP = barvinok_enumerate_e(para, exist, nparam, MaxRays);
4038 evalue *E;
4039 if (!emptyQ(pos)) {
4040 E = barvinok_enumerate_e(pos, exist, nparam, MaxRays);
4041 eadd(E, EP);
4042 free_evalue_refs(E);
4043 free(E);
4045 if (!emptyQ(neg)) {
4046 E = barvinok_enumerate_e(neg, exist, nparam, MaxRays);
4047 eadd(E, EP);
4048 free_evalue_refs(E);
4049 free(E);
4051 Polyhedron_Free(para);
4052 Polyhedron_Free(pos);
4053 Polyhedron_Free(neg);
4054 break;
4056 if (EP)
4057 break;
4058 } END_FORALL_PVertex_in_ParamPolyhedron;
4059 if (EP)
4060 break;
4061 } END_FORALL_PVertex_in_ParamPolyhedron;
4063 if (!EP) {
4064 /* Search for vertex coordinate to split on */
4065 /* First look for one independent of the parameters */
4066 FORALL_PVertex_in_ParamPolyhedron(V, last, PP) {
4067 for (int i = 0; i < exist; ++i) {
4068 int j;
4069 for (j = 0; j < nparam; ++j)
4070 if (value_notzero_p(V->Vertex->p[i][j]))
4071 break;
4072 if (j < nparam)
4073 continue;
4074 value_set_si(M->p[0][0], 1);
4075 Vector_Set(M->p[0]+1, 0, nvar+exist);
4076 Vector_Copy(V->Vertex->p[i],
4077 M->p[0] + 1 + nvar + exist, nparam+1);
4078 value_oppose(M->p[0][1+nvar+i],
4079 V->Vertex->p[i][nparam+1]);
4081 Polyhedron *pos, *neg;
4082 value_set_si(M->p[0][0], 1);
4083 value_decrement(M->p[0][P->Dimension+1],
4084 M->p[0][P->Dimension+1]);
4085 neg = AddConstraints(M->p[0], 1, P, MaxRays);
4086 value_set_si(f, -1);
4087 Vector_Scale(M->p[0]+1, M->p[0]+1, f,
4088 P->Dimension+1);
4089 value_decrement(M->p[0][P->Dimension+1],
4090 M->p[0][P->Dimension+1]);
4091 value_decrement(M->p[0][P->Dimension+1],
4092 M->p[0][P->Dimension+1]);
4093 pos = AddConstraints(M->p[0], 1, P, MaxRays);
4094 if (emptyQ(neg) || emptyQ(pos)) {
4095 Polyhedron_Free(pos);
4096 Polyhedron_Free(neg);
4097 continue;
4099 Polyhedron_Free(pos);
4100 value_increment(M->p[0][P->Dimension+1],
4101 M->p[0][P->Dimension+1]);
4102 pos = AddConstraints(M->p[0], 1, P, MaxRays);
4103 #ifdef DEBUG_ER
4104 fprintf(stderr, "\nER: Vertex\n");
4105 #endif /* DEBUG_ER */
4106 pos->next = neg;
4107 EP = enumerate_or(pos, exist, nparam, MaxRays);
4108 break;
4110 if (EP)
4111 break;
4112 } END_FORALL_PVertex_in_ParamPolyhedron;
4115 if (!EP) {
4116 /* Search for vertex coordinate to split on */
4117 /* Now look for one that depends on the parameters */
4118 FORALL_PVertex_in_ParamPolyhedron(V, last, PP) {
4119 for (int i = 0; i < exist; ++i) {
4120 value_set_si(M->p[0][0], 1);
4121 Vector_Set(M->p[0]+1, 0, nvar+exist);
4122 Vector_Copy(V->Vertex->p[i],
4123 M->p[0] + 1 + nvar + exist, nparam+1);
4124 value_oppose(M->p[0][1+nvar+i],
4125 V->Vertex->p[i][nparam+1]);
4127 Polyhedron *pos, *neg;
4128 value_set_si(M->p[0][0], 1);
4129 value_decrement(M->p[0][P->Dimension+1],
4130 M->p[0][P->Dimension+1]);
4131 neg = AddConstraints(M->p[0], 1, P, MaxRays);
4132 value_set_si(f, -1);
4133 Vector_Scale(M->p[0]+1, M->p[0]+1, f,
4134 P->Dimension+1);
4135 value_decrement(M->p[0][P->Dimension+1],
4136 M->p[0][P->Dimension+1]);
4137 value_decrement(M->p[0][P->Dimension+1],
4138 M->p[0][P->Dimension+1]);
4139 pos = AddConstraints(M->p[0], 1, P, MaxRays);
4140 if (emptyQ(neg) || emptyQ(pos)) {
4141 Polyhedron_Free(pos);
4142 Polyhedron_Free(neg);
4143 continue;
4145 Polyhedron_Free(pos);
4146 value_increment(M->p[0][P->Dimension+1],
4147 M->p[0][P->Dimension+1]);
4148 pos = AddConstraints(M->p[0], 1, P, MaxRays);
4149 #ifdef DEBUG_ER
4150 fprintf(stderr, "\nER: ParamVertex\n");
4151 #endif /* DEBUG_ER */
4152 pos->next = neg;
4153 EP = enumerate_or(pos, exist, nparam, MaxRays);
4154 break;
4156 if (EP)
4157 break;
4158 } END_FORALL_PVertex_in_ParamPolyhedron;
4161 Matrix_Free(M);
4162 value_clear(f);
4165 if (CEq)
4166 Polyhedron_Free(CEq);
4167 if (CT)
4168 Matrix_Free(CT);
4169 if (PP)
4170 Param_Polyhedron_Free(PP);
4171 *PA = P;
4173 return EP;
4176 #ifndef HAVE_PIPLIB
4177 evalue *barvinok_enumerate_pip(Polyhedron *P,
4178 unsigned exist, unsigned nparam, unsigned MaxRays)
4180 return 0;
4182 #else
4183 evalue *barvinok_enumerate_pip(Polyhedron *P,
4184 unsigned exist, unsigned nparam, unsigned MaxRays)
4186 int nvar = P->Dimension - exist - nparam;
4187 evalue *EP = evalue_zero();
4188 Polyhedron *Q, *N;
4190 #ifdef DEBUG_ER
4191 fprintf(stderr, "\nER: PIP\n");
4192 #endif /* DEBUG_ER */
4194 Polyhedron *D = pip_projectout(P, nvar, exist, nparam);
4195 for (Q = D; Q; Q = N) {
4196 N = Q->next;
4197 Q->next = 0;
4198 evalue *E;
4199 exist = Q->Dimension - nvar - nparam;
4200 E = barvinok_enumerate_e(Q, exist, nparam, MaxRays);
4201 Polyhedron_Free(Q);
4202 eadd(E, EP);
4203 free_evalue_refs(E);
4204 free(E);
4207 return EP;
4209 #endif
4212 static bool is_single(Value *row, int pos, int len)
4214 return First_Non_Zero(row, pos) == -1 &&
4215 First_Non_Zero(row+pos+1, len-pos-1) == -1;
4218 static evalue* barvinok_enumerate_e_r(Polyhedron *P,
4219 unsigned exist, unsigned nparam, unsigned MaxRays);
4221 #ifdef DEBUG_ER
4222 static int er_level = 0;
4224 evalue* barvinok_enumerate_e(Polyhedron *P,
4225 unsigned exist, unsigned nparam, unsigned MaxRays)
4227 fprintf(stderr, "\nER: level %i\n", er_level);
4229 Polyhedron_PrintConstraints(stderr, P_VALUE_FMT, P);
4230 fprintf(stderr, "\nE %d\nP %d\n", exist, nparam);
4231 ++er_level;
4232 P = DomainConstraintSimplify(Polyhedron_Copy(P), MaxRays);
4233 evalue *EP = barvinok_enumerate_e_r(P, exist, nparam, MaxRays);
4234 Polyhedron_Free(P);
4235 --er_level;
4236 return EP;
4238 #else
4239 evalue* barvinok_enumerate_e(Polyhedron *P,
4240 unsigned exist, unsigned nparam, unsigned MaxRays)
4242 P = DomainConstraintSimplify(Polyhedron_Copy(P), MaxRays);
4243 evalue *EP = barvinok_enumerate_e_r(P, exist, nparam, MaxRays);
4244 Polyhedron_Free(P);
4245 return EP;
4247 #endif
4249 static evalue* barvinok_enumerate_e_r(Polyhedron *P,
4250 unsigned exist, unsigned nparam, unsigned MaxRays)
4252 if (exist == 0) {
4253 Polyhedron *U = Universe_Polyhedron(nparam);
4254 evalue *EP = barvinok_enumerate_ev(P, U, MaxRays);
4255 //char *param_name[] = {"P", "Q", "R", "S", "T" };
4256 //print_evalue(stdout, EP, param_name);
4257 Polyhedron_Free(U);
4258 return EP;
4261 int nvar = P->Dimension - exist - nparam;
4262 int len = P->Dimension + 2;
4264 /* for now */
4265 POL_ENSURE_FACETS(P);
4266 POL_ENSURE_VERTICES(P);
4268 if (emptyQ(P))
4269 return evalue_zero();
4271 if (nvar == 0 && nparam == 0) {
4272 evalue *EP = evalue_zero();
4273 barvinok_count(P, &EP->x.n, MaxRays);
4274 if (value_pos_p(EP->x.n))
4275 value_set_si(EP->x.n, 1);
4276 return EP;
4279 int r;
4280 for (r = 0; r < P->NbRays; ++r)
4281 if (value_zero_p(P->Ray[r][0]) ||
4282 value_zero_p(P->Ray[r][P->Dimension+1])) {
4283 int i;
4284 for (i = 0; i < nvar; ++i)
4285 if (value_notzero_p(P->Ray[r][i+1]))
4286 break;
4287 if (i >= nvar)
4288 continue;
4289 for (i = nvar + exist; i < nvar + exist + nparam; ++i)
4290 if (value_notzero_p(P->Ray[r][i+1]))
4291 break;
4292 if (i >= nvar + exist + nparam)
4293 break;
4295 if (r < P->NbRays) {
4296 evalue *EP = evalue_zero();
4297 value_set_si(EP->x.n, -1);
4298 return EP;
4301 int first;
4302 for (r = 0; r < P->NbEq; ++r)
4303 if ((first = First_Non_Zero(P->Constraint[r]+1+nvar, exist)) != -1)
4304 break;
4305 if (r < P->NbEq) {
4306 if (First_Non_Zero(P->Constraint[r]+1+nvar+first+1,
4307 exist-first-1) != -1) {
4308 Polyhedron *T = rotate_along(P, r, nvar, exist, MaxRays);
4309 #ifdef DEBUG_ER
4310 fprintf(stderr, "\nER: Equality\n");
4311 #endif /* DEBUG_ER */
4312 evalue *EP = barvinok_enumerate_e(T, exist-1, nparam, MaxRays);
4313 Polyhedron_Free(T);
4314 return EP;
4315 } else {
4316 #ifdef DEBUG_ER
4317 fprintf(stderr, "\nER: Fixed\n");
4318 #endif /* DEBUG_ER */
4319 if (first == 0)
4320 return barvinok_enumerate_e(P, exist-1, nparam, MaxRays);
4321 else {
4322 Polyhedron *T = Polyhedron_Copy(P);
4323 SwapColumns(T, nvar+1, nvar+1+first);
4324 evalue *EP = barvinok_enumerate_e(T, exist-1, nparam, MaxRays);
4325 Polyhedron_Free(T);
4326 return EP;
4331 Vector *row = Vector_Alloc(len);
4332 value_set_si(row->p[0], 1);
4334 Value f;
4335 value_init(f);
4337 enum constraint* info = new constraint[exist];
4338 for (int i = 0; i < exist; ++i) {
4339 info[i] = ALL_POS;
4340 for (int l = P->NbEq; l < P->NbConstraints; ++l) {
4341 if (value_negz_p(P->Constraint[l][nvar+i+1]))
4342 continue;
4343 bool l_parallel = is_single(P->Constraint[l]+nvar+1, i, exist);
4344 for (int u = P->NbEq; u < P->NbConstraints; ++u) {
4345 if (value_posz_p(P->Constraint[u][nvar+i+1]))
4346 continue;
4347 bool lu_parallel = l_parallel ||
4348 is_single(P->Constraint[u]+nvar+1, i, exist);
4349 value_oppose(f, P->Constraint[u][nvar+i+1]);
4350 Vector_Combine(P->Constraint[l]+1, P->Constraint[u]+1, row->p+1,
4351 f, P->Constraint[l][nvar+i+1], len-1);
4352 if (!(info[i] & INDEPENDENT)) {
4353 int j;
4354 for (j = 0; j < exist; ++j)
4355 if (j != i && value_notzero_p(row->p[nvar+j+1]))
4356 break;
4357 if (j == exist) {
4358 //printf("independent: i: %d, l: %d, u: %d\n", i, l, u);
4359 info[i] = (constraint)(info[i] | INDEPENDENT);
4362 if (info[i] & ALL_POS) {
4363 value_addto(row->p[len-1], row->p[len-1],
4364 P->Constraint[l][nvar+i+1]);
4365 value_addto(row->p[len-1], row->p[len-1], f);
4366 value_multiply(f, f, P->Constraint[l][nvar+i+1]);
4367 value_subtract(row->p[len-1], row->p[len-1], f);
4368 value_decrement(row->p[len-1], row->p[len-1]);
4369 ConstraintSimplify(row->p, row->p, len, &f);
4370 value_set_si(f, -1);
4371 Vector_Scale(row->p+1, row->p+1, f, len-1);
4372 value_decrement(row->p[len-1], row->p[len-1]);
4373 Polyhedron *T = AddConstraints(row->p, 1, P, MaxRays);
4374 if (!emptyQ(T)) {
4375 //printf("not all_pos: i: %d, l: %d, u: %d\n", i, l, u);
4376 info[i] = (constraint)(info[i] ^ ALL_POS);
4378 //puts("pos remainder");
4379 //Polyhedron_Print(stdout, P_VALUE_FMT, T);
4380 Polyhedron_Free(T);
4382 if (!(info[i] & ONE_NEG)) {
4383 if (lu_parallel) {
4384 negative_test_constraint(P->Constraint[l],
4385 P->Constraint[u],
4386 row->p, nvar+i, len, &f);
4387 oppose_constraint(row->p, len, &f);
4388 Polyhedron *T = AddConstraints(row->p, 1, P, MaxRays);
4389 if (emptyQ(T)) {
4390 //printf("one_neg i: %d, l: %d, u: %d\n", i, l, u);
4391 info[i] = (constraint)(info[i] | ONE_NEG);
4393 //puts("neg remainder");
4394 //Polyhedron_Print(stdout, P_VALUE_FMT, T);
4395 Polyhedron_Free(T);
4396 } else if (!(info[i] & ROT_NEG)) {
4397 if (parallel_constraints(P->Constraint[l],
4398 P->Constraint[u],
4399 row->p, nvar, exist)) {
4400 negative_test_constraint7(P->Constraint[l],
4401 P->Constraint[u],
4402 row->p, nvar, exist,
4403 len, &f);
4404 oppose_constraint(row->p, len, &f);
4405 Polyhedron *T = AddConstraints(row->p, 1, P, MaxRays);
4406 if (emptyQ(T)) {
4407 // printf("rot_neg i: %d, l: %d, u: %d\n", i, l, u);
4408 info[i] = (constraint)(info[i] | ROT_NEG);
4409 r = l;
4411 //puts("neg remainder");
4412 //Polyhedron_Print(stdout, P_VALUE_FMT, T);
4413 Polyhedron_Free(T);
4417 if (!(info[i] & ALL_POS) && (info[i] & (ONE_NEG | ROT_NEG)))
4418 goto next;
4421 if (info[i] & ALL_POS)
4422 break;
4423 next:
4428 for (int i = 0; i < exist; ++i)
4429 printf("%i: %i\n", i, info[i]);
4431 for (int i = 0; i < exist; ++i)
4432 if (info[i] & ALL_POS) {
4433 #ifdef DEBUG_ER
4434 fprintf(stderr, "\nER: Positive\n");
4435 #endif /* DEBUG_ER */
4436 // Eliminate
4437 // Maybe we should chew off some of the fat here
4438 Matrix *M = Matrix_Alloc(P->Dimension, P->Dimension+1);
4439 for (int j = 0; j < P->Dimension; ++j)
4440 value_set_si(M->p[j][j + (j >= i+nvar)], 1);
4441 Polyhedron *T = Polyhedron_Image(P, M, MaxRays);
4442 Matrix_Free(M);
4443 evalue *EP = barvinok_enumerate_e(T, exist-1, nparam, MaxRays);
4444 Polyhedron_Free(T);
4445 value_clear(f);
4446 Vector_Free(row);
4447 delete [] info;
4448 return EP;
4450 for (int i = 0; i < exist; ++i)
4451 if (info[i] & ONE_NEG) {
4452 #ifdef DEBUG_ER
4453 fprintf(stderr, "\nER: Negative\n");
4454 #endif /* DEBUG_ER */
4455 Vector_Free(row);
4456 value_clear(f);
4457 delete [] info;
4458 if (i == 0)
4459 return barvinok_enumerate_e(P, exist-1, nparam, MaxRays);
4460 else {
4461 Polyhedron *T = Polyhedron_Copy(P);
4462 SwapColumns(T, nvar+1, nvar+1+i);
4463 evalue *EP = barvinok_enumerate_e(T, exist-1, nparam, MaxRays);
4464 Polyhedron_Free(T);
4465 return EP;
4468 for (int i = 0; i < exist; ++i)
4469 if (info[i] & ROT_NEG) {
4470 #ifdef DEBUG_ER
4471 fprintf(stderr, "\nER: Rotate\n");
4472 #endif /* DEBUG_ER */
4473 Vector_Free(row);
4474 value_clear(f);
4475 delete [] info;
4476 Polyhedron *T = rotate_along(P, r, nvar, exist, MaxRays);
4477 evalue *EP = barvinok_enumerate_e(T, exist-1, nparam, MaxRays);
4478 Polyhedron_Free(T);
4479 return EP;
4481 for (int i = 0; i < exist; ++i)
4482 if (info[i] & INDEPENDENT) {
4483 Polyhedron *pos, *neg;
4485 /* Find constraint again and split off negative part */
4487 if (SplitOnVar(P, i, nvar, exist, MaxRays,
4488 row, f, true, &pos, &neg)) {
4489 #ifdef DEBUG_ER
4490 fprintf(stderr, "\nER: Split\n");
4491 #endif /* DEBUG_ER */
4493 evalue *EP =
4494 barvinok_enumerate_e(neg, exist-1, nparam, MaxRays);
4495 evalue *E =
4496 barvinok_enumerate_e(pos, exist, nparam, MaxRays);
4497 eadd(E, EP);
4498 free_evalue_refs(E);
4499 free(E);
4500 Polyhedron_Free(neg);
4501 Polyhedron_Free(pos);
4502 value_clear(f);
4503 Vector_Free(row);
4504 delete [] info;
4505 return EP;
4508 delete [] info;
4510 Polyhedron *O = P;
4511 Polyhedron *F;
4513 evalue *EP;
4515 EP = enumerate_line(P, exist, nparam, MaxRays);
4516 if (EP)
4517 goto out;
4519 EP = barvinok_enumerate_pip(P, exist, nparam, MaxRays);
4520 if (EP)
4521 goto out;
4523 EP = enumerate_redundant_ray(P, exist, nparam, MaxRays);
4524 if (EP)
4525 goto out;
4527 EP = enumerate_sure(P, exist, nparam, MaxRays);
4528 if (EP)
4529 goto out;
4531 EP = enumerate_ray(P, exist, nparam, MaxRays);
4532 if (EP)
4533 goto out;
4535 EP = enumerate_sure2(P, exist, nparam, MaxRays);
4536 if (EP)
4537 goto out;
4539 F = unfringe(P, MaxRays);
4540 if (!PolyhedronIncludes(F, P)) {
4541 #ifdef DEBUG_ER
4542 fprintf(stderr, "\nER: Fringed\n");
4543 #endif /* DEBUG_ER */
4544 EP = barvinok_enumerate_e(F, exist, nparam, MaxRays);
4545 Polyhedron_Free(F);
4546 goto out;
4548 Polyhedron_Free(F);
4550 if (nparam)
4551 EP = enumerate_vd(&P, exist, nparam, MaxRays);
4552 if (EP)
4553 goto out2;
4555 if (nvar != 0) {
4556 EP = enumerate_sum(P, exist, nparam, MaxRays);
4557 goto out2;
4560 assert(nvar == 0);
4562 int i;
4563 Polyhedron *pos, *neg;
4564 for (i = 0; i < exist; ++i)
4565 if (SplitOnVar(P, i, nvar, exist, MaxRays,
4566 row, f, false, &pos, &neg))
4567 break;
4569 assert (i < exist);
4571 pos->next = neg;
4572 EP = enumerate_or(pos, exist, nparam, MaxRays);
4574 out2:
4575 if (O != P)
4576 Polyhedron_Free(P);
4578 out:
4579 value_clear(f);
4580 Vector_Free(row);
4581 return EP;
4584 static void split_param_compression(Matrix *CP, mat_ZZ& map, vec_ZZ& offset)
4586 Matrix *T = Transpose(CP);
4587 matrix2zz(T, map, T->NbRows-1, T->NbColumns-1);
4588 values2zz(T->p[T->NbRows-1], offset, T->NbColumns-1);
4589 Matrix_Free(T);
4593 * remove equalities that require a "compression" of the parameters
4595 #ifndef HAVE_COMPRESS_PARMS
4596 static Polyhedron *remove_more_equalities(Polyhedron *P, unsigned nparam,
4597 Matrix **CP, unsigned MaxRays)
4599 return P;
4601 #else
4602 static Polyhedron *remove_more_equalities(Polyhedron *P, unsigned nparam,
4603 Matrix **CP, unsigned MaxRays)
4605 Matrix *M, *T;
4606 Polyhedron *Q;
4608 /* compress_parms doesn't like equalities that only involve parameters */
4609 for (int i = 0; i < P->NbEq; ++i)
4610 assert(First_Non_Zero(P->Constraint[i]+1, P->Dimension-nparam) != -1);
4612 M = Matrix_Alloc(P->NbEq, P->Dimension+2);
4613 Vector_Copy(P->Constraint[0], M->p[0], P->NbEq * (P->Dimension+2));
4614 *CP = compress_parms(M, nparam);
4615 T = align_matrix(*CP, P->Dimension+1);
4616 Q = Polyhedron_Preimage(P, T, MaxRays);
4617 Polyhedron_Free(P);
4618 P = Q;
4619 P = remove_equalities_p(P, P->Dimension-nparam, NULL);
4620 Matrix_Free(T);
4621 Matrix_Free(M);
4622 return P;
4624 #endif
4626 gen_fun * barvinok_series(Polyhedron *P, Polyhedron* C, unsigned MaxRays)
4628 Matrix *CP = NULL;
4629 Polyhedron *CA;
4630 unsigned nparam = C->Dimension;
4632 CA = align_context(C, P->Dimension, MaxRays);
4633 P = DomainIntersection(P, CA, MaxRays);
4634 Polyhedron_Free(CA);
4636 if (emptyQ2(P)) {
4637 Polyhedron_Free(P);
4638 return new gen_fun;
4641 assert(!Polyhedron_is_infinite(P, nparam));
4642 assert(P->NbBid == 0);
4643 assert(Polyhedron_has_positive_rays(P, nparam));
4644 if (P->NbEq != 0)
4645 P = remove_equalities_p(P, P->Dimension-nparam, NULL);
4646 if (P->NbEq != 0)
4647 P = remove_more_equalities(P, nparam, &CP, MaxRays);
4648 assert(P->NbEq == 0);
4650 #ifdef USE_INCREMENTAL_BF
4651 partial_bfcounter red(P, nparam);
4652 #elif defined USE_INCREMENTAL_DF
4653 partial_ireducer red(P, nparam);
4654 #else
4655 partial_reducer red(P, nparam);
4656 #endif
4657 red.start(MaxRays);
4658 Polyhedron_Free(P);
4659 if (CP) {
4660 mat_ZZ map;
4661 vec_ZZ offset;
4662 split_param_compression(CP, map, offset);
4663 red.gf->substitute(CP, map, offset);
4664 Matrix_Free(CP);
4666 return red.gf;
4669 static Polyhedron *skew_into_positive_orthant(Polyhedron *D, unsigned nparam,
4670 unsigned MaxRays)
4672 Matrix *M = NULL;
4673 Value tmp;
4674 value_init(tmp);
4675 for (Polyhedron *P = D; P; P = P->next) {
4676 POL_ENSURE_VERTICES(P);
4677 assert(!Polyhedron_is_infinite(P, nparam));
4678 assert(P->NbBid == 0);
4679 assert(Polyhedron_has_positive_rays(P, nparam));
4681 for (int r = 0; r < P->NbRays; ++r) {
4682 if (value_notzero_p(P->Ray[r][P->Dimension+1]))
4683 continue;
4684 for (int i = 0; i < nparam; ++i) {
4685 int j;
4686 if (value_posz_p(P->Ray[r][i+1]))
4687 continue;
4688 if (!M) {
4689 M = Matrix_Alloc(D->Dimension+1, D->Dimension+1);
4690 for (int i = 0; i < D->Dimension+1; ++i)
4691 value_set_si(M->p[i][i], 1);
4692 } else {
4693 Inner_Product(P->Ray[r]+1, M->p[i], D->Dimension+1, &tmp);
4694 if (value_posz_p(tmp))
4695 continue;
4697 for (j = P->Dimension - nparam; j < P->Dimension; ++j)
4698 if (value_pos_p(P->Ray[r][j+1]))
4699 break;
4700 assert(j < P->Dimension);
4701 value_pdivision(tmp, P->Ray[r][j+1], P->Ray[r][i+1]);
4702 value_subtract(M->p[i][j], M->p[i][j], tmp);
4706 value_clear(tmp);
4707 if (M) {
4708 D = DomainImage(D, M, MaxRays);
4709 Matrix_Free(M);
4711 return D;
4714 evalue* barvinok_enumerate_union(Polyhedron *D, Polyhedron* C, unsigned MaxRays)
4716 evalue *EP;
4717 Polyhedron *conv, *D2;
4718 gen_fun *gf = NULL;
4719 unsigned nparam = C->Dimension;
4720 ZZ one, mone;
4721 one = 1;
4722 mone = -1;
4723 D2 = skew_into_positive_orthant(D, nparam, MaxRays);
4724 for (Polyhedron *P = D2; P; P = P->next) {
4725 assert(P->Dimension == D2->Dimension);
4726 POL_ENSURE_VERTICES(P);
4727 /* it doesn't matter which reducer we use, since we don't actually
4728 * reduce anything here
4730 partial_reducer red(P, P->Dimension);
4731 red.start(MaxRays);
4732 if (!gf)
4733 gf = red.gf;
4734 else {
4735 gf->add_union(red.gf, MaxRays);
4736 delete red.gf;
4739 /* we actually only need the convex union of the parameter space
4740 * but the reducer classes currently expect a polyhedron in
4741 * the combined space
4743 conv = DomainConvex(D2, MaxRays);
4744 #ifdef USE_INCREMENTAL_DF
4745 partial_ireducer red(conv, nparam);
4746 #else
4747 partial_reducer red(conv, nparam);
4748 #endif
4749 for (int i = 0; i < gf->term.size(); ++i) {
4750 for (int j = 0; j < gf->term[i]->n.power.NumRows(); ++j) {
4751 red.reduce(gf->term[i]->n.coeff[j][0], gf->term[i]->n.coeff[j][1],
4752 gf->term[i]->n.power[j], gf->term[i]->d.power);
4755 delete gf;
4756 if (D != D2)
4757 Domain_Free(D2);
4758 Polyhedron_Free(conv);
4759 EP = *red.gf;
4760 delete red.gf;
4761 return EP;