8 #include <NTL/mat_ZZ.h>
12 #include <polylib/polylibgmp.h>
13 #include "ev_operations.h"
26 using std::ostringstream
;
28 #define ALLOC(p) (((long *) (p))[0])
29 #define SIZE(p) (((long *) (p))[1])
30 #define DATA(p) ((mp_limb_t *) (((long *) (p)) + 2))
32 static void value2zz(Value v
, ZZ
& z
)
34 int sa
= v
[0]._mp_size
;
35 int abs_sa
= sa
< 0 ? -sa
: sa
;
37 _ntl_gsetlength(&z
.rep
, abs_sa
);
38 mp_limb_t
* adata
= DATA(z
.rep
);
39 for (int i
= 0; i
< abs_sa
; ++i
)
40 adata
[i
] = v
[0]._mp_d
[i
];
44 static void zz2value(ZZ
& z
, Value
& v
)
52 int abs_sa
= sa
< 0 ? -sa
: sa
;
54 mp_limb_t
* adata
= DATA(z
.rep
);
55 mpz_realloc2(v
, __GMP_BITS_PER_MP_LIMB
* abs_sa
);
56 for (int i
= 0; i
< abs_sa
; ++i
)
57 v
[0]._mp_d
[i
] = adata
[i
];
62 #define ALLOC(p) p = (typeof(p))malloc(sizeof(*p))
65 * We just ignore the last column and row
66 * If the final element is not equal to one
67 * then the result will actually be a multiple of the input
69 static void matrix2zz(Matrix
*M
, mat_ZZ
& m
, unsigned nr
, unsigned nc
)
73 for (int i
= 0; i
< nr
; ++i
) {
74 // assert(value_one_p(M->p[i][M->NbColumns - 1]));
75 for (int j
= 0; j
< nc
; ++j
) {
76 value2zz(M
->p
[i
][j
], m
[i
][j
]);
81 static void values2zz(Value
*p
, vec_ZZ
& v
, int len
)
85 for (int i
= 0; i
< len
; ++i
) {
92 static void zz2values(vec_ZZ
& v
, Value
*p
)
94 for (int i
= 0; i
< v
.length(); ++i
)
98 static void rays(mat_ZZ
& r
, Polyhedron
*C
)
100 unsigned dim
= C
->NbRays
- 1; /* don't count zero vertex */
101 assert(C
->NbRays
- 1 == C
->Dimension
);
106 for (i
= 0, c
= 0; i
< dim
; ++i
)
107 if (value_zero_p(C
->Ray
[i
][dim
+1])) {
108 for (int j
= 0; j
< dim
; ++j
) {
109 value2zz(C
->Ray
[i
][j
+1], tmp
);
116 static Matrix
* rays(Polyhedron
*C
)
118 unsigned dim
= C
->NbRays
- 1; /* don't count zero vertex */
119 assert(C
->NbRays
- 1 == C
->Dimension
);
121 Matrix
*M
= Matrix_Alloc(dim
+1, dim
+1);
125 for (i
= 0, c
= 0; i
<= dim
&& c
< dim
; ++i
)
126 if (value_zero_p(C
->Ray
[i
][dim
+1])) {
127 Vector_Copy(C
->Ray
[i
] + 1, M
->p
[c
], dim
);
128 value_set_si(M
->p
[c
++][dim
], 0);
131 value_set_si(M
->p
[dim
][dim
], 1);
136 static Matrix
* rays2(Polyhedron
*C
)
138 unsigned dim
= C
->NbRays
- 1; /* don't count zero vertex */
139 assert(C
->NbRays
- 1 == C
->Dimension
);
141 Matrix
*M
= Matrix_Alloc(dim
, dim
);
145 for (i
= 0, c
= 0; i
<= dim
&& c
< dim
; ++i
)
146 if (value_zero_p(C
->Ray
[i
][dim
+1]))
147 Vector_Copy(C
->Ray
[i
] + 1, M
->p
[c
++], dim
);
154 * Returns the largest absolute value in the vector
156 static ZZ
max(vec_ZZ
& v
)
159 for (int i
= 1; i
< v
.length(); ++i
)
169 Rays
= Matrix_Copy(M
);
172 cone(Polyhedron
*C
) {
173 Cone
= Polyhedron_Copy(C
);
179 matrix2zz(Rays
, A
, Rays
->NbRows
- 1, Rays
->NbColumns
- 1);
180 det
= determinant(A
);
187 Vector
* short_vector(vec_ZZ
& lambda
) {
188 Matrix
*M
= Matrix_Copy(Rays
);
189 Matrix
*inv
= Matrix_Alloc(M
->NbRows
, M
->NbColumns
);
190 int ok
= Matrix_Inverse(M
, inv
);
197 matrix2zz(inv
, B
, inv
->NbRows
- 1, inv
->NbColumns
- 1);
198 long r
= LLL(det2
, B
, U
);
202 for (int i
= 1; i
< B
.NumRows(); ++i
) {
214 Vector
*z
= Vector_Alloc(U
[index
].length()+1);
216 zz2values(U
[index
], z
->p
);
217 value_set_si(z
->p
[U
[index
].length()], 0);
221 Polyhedron
*C
= poly();
223 for (i
= 0; i
< C
->NbConstraints
; ++i
) {
224 Inner_Product(z
->p
, C
->Constraint
[i
]+1, z
->Size
-1, &tmp
);
225 if (value_pos_p(tmp
))
228 if (i
== C
->NbConstraints
) {
229 value_set_si(tmp
, -1);
230 Vector_Scale(z
->p
, z
->p
, tmp
, z
->Size
-1);
237 Polyhedron_Free(Cone
);
243 Matrix
*M
= Matrix_Alloc(Rays
->NbRows
+1, Rays
->NbColumns
+1);
244 for (int i
= 0; i
< Rays
->NbRows
; ++i
) {
245 Vector_Copy(Rays
->p
[i
], M
->p
[i
]+1, Rays
->NbColumns
);
246 value_set_si(M
->p
[i
][0], 1);
248 Vector_Set(M
->p
[Rays
->NbRows
]+1, 0, Rays
->NbColumns
-1);
249 value_set_si(M
->p
[Rays
->NbRows
][0], 1);
250 value_set_si(M
->p
[Rays
->NbRows
][Rays
->NbColumns
], 1);
251 Cone
= Rays2Polyhedron(M
, M
->NbRows
+1);
252 assert(Cone
->NbConstraints
== Cone
->NbRays
);
266 dpoly(int d
, ZZ
& degree
, int offset
= 0) {
267 coeff
.SetLength(d
+1);
269 int min
= d
+ offset
;
270 if (degree
< ZZ(INIT_VAL
, min
))
271 min
= to_int(degree
);
273 ZZ c
= ZZ(INIT_VAL
, 1);
276 for (int i
= 1; i
<= min
; ++i
) {
277 c
*= (degree
-i
+ 1);
282 void operator *= (dpoly
& f
) {
283 assert(coeff
.length() == f
.coeff
.length());
285 coeff
= f
.coeff
[0] * coeff
;
286 for (int i
= 1; i
< coeff
.length(); ++i
)
287 for (int j
= 0; i
+j
< coeff
.length(); ++j
)
288 coeff
[i
+j
] += f
.coeff
[i
] * old
[j
];
290 void div(dpoly
& d
, mpq_t count
, ZZ
& sign
) {
291 int len
= coeff
.length();
294 mpq_t
* c
= new mpq_t
[coeff
.length()];
297 for (int i
= 0; i
< len
; ++i
) {
299 zz2value(coeff
[i
], tmp
);
300 mpq_set_z(c
[i
], tmp
);
302 for (int j
= 1; j
<= i
; ++j
) {
303 zz2value(d
.coeff
[j
], tmp
);
304 mpq_set_z(qtmp
, tmp
);
305 mpq_mul(qtmp
, qtmp
, c
[i
-j
]);
306 mpq_sub(c
[i
], c
[i
], qtmp
);
309 zz2value(d
.coeff
[0], tmp
);
310 mpq_set_z(qtmp
, tmp
);
311 mpq_div(c
[i
], c
[i
], qtmp
);
314 mpq_sub(count
, count
, c
[len
-1]);
316 mpq_add(count
, count
, c
[len
-1]);
320 for (int i
= 0; i
< len
; ++i
)
332 dpoly_n(int d
, ZZ
& degree_0
, ZZ
& degree_1
, int offset
= 0) {
336 zz2value(degree_0
, d0
);
337 zz2value(degree_1
, d1
);
338 coeff
= Matrix_Alloc(d
+1, d
+1+1);
339 value_set_si(coeff
->p
[0][0], 1);
340 value_set_si(coeff
->p
[0][d
+1], 1);
341 for (int i
= 1; i
<= d
; ++i
) {
342 value_multiply(coeff
->p
[i
][0], coeff
->p
[i
-1][0], d0
);
343 Vector_Combine(coeff
->p
[i
-1], coeff
->p
[i
-1]+1, coeff
->p
[i
]+1,
345 value_set_si(coeff
->p
[i
][d
+1], i
);
346 value_multiply(coeff
->p
[i
][d
+1], coeff
->p
[i
][d
+1], coeff
->p
[i
-1][d
+1]);
347 value_decrement(d0
, d0
);
352 void div(dpoly
& d
, Vector
*count
, ZZ
& sign
) {
353 int len
= coeff
->NbRows
;
354 Matrix
* c
= Matrix_Alloc(coeff
->NbRows
, coeff
->NbColumns
);
357 for (int i
= 0; i
< len
; ++i
) {
358 Vector_Copy(coeff
->p
[i
], c
->p
[i
], len
+1);
359 for (int j
= 1; j
<= i
; ++j
) {
360 zz2value(d
.coeff
[j
], tmp
);
361 value_multiply(tmp
, tmp
, c
->p
[i
][len
]);
362 value_oppose(tmp
, tmp
);
363 Vector_Combine(c
->p
[i
], c
->p
[i
-j
], c
->p
[i
],
364 c
->p
[i
-j
][len
], tmp
, len
);
365 value_multiply(c
->p
[i
][len
], c
->p
[i
][len
], c
->p
[i
-j
][len
]);
367 zz2value(d
.coeff
[0], tmp
);
368 value_multiply(c
->p
[i
][len
], c
->p
[i
][len
], tmp
);
371 value_set_si(tmp
, -1);
372 Vector_Scale(c
->p
[len
-1], count
->p
, tmp
, len
);
373 value_assign(count
->p
[len
], c
->p
[len
-1][len
]);
375 Vector_Copy(c
->p
[len
-1], count
->p
, len
+1);
376 Vector_Normalize(count
->p
, len
+1);
383 * Barvinok's Decomposition of a simplicial cone
385 * Returns two lists of polyhedra
387 void barvinok_decompose(Polyhedron
*C
, Polyhedron
**ppos
, Polyhedron
**pneg
)
389 Polyhedron
*pos
= *ppos
, *neg
= *pneg
;
390 vector
<cone
*> nonuni
;
391 cone
* c
= new cone(C
);
398 Polyhedron
*p
= Polyhedron_Copy(c
->Cone
);
404 while (!nonuni
.empty()) {
407 Vector
* v
= c
->short_vector(lambda
);
408 for (int i
= 0; i
< c
->Rays
->NbRows
- 1; ++i
) {
411 Matrix
* M
= Matrix_Copy(c
->Rays
);
412 Vector_Copy(v
->p
, M
->p
[i
], v
->Size
);
413 cone
* pc
= new cone(M
);
414 assert (pc
->det
!= 0);
415 if (abs(pc
->det
) > 1) {
416 assert(abs(pc
->det
) < abs(c
->det
));
417 nonuni
.push_back(pc
);
419 Polyhedron
*p
= pc
->poly();
421 if (sign(pc
->det
) == s
) {
440 * Returns a single list of npos "positive" cones followed by nneg
442 * The input cone is freed
444 void decompose(Polyhedron
*cone
, Polyhedron
**parts
, int *npos
, int *nneg
, unsigned MaxRays
)
446 Polyhedron_Polarize(cone
);
447 if (cone
->NbRays
- 1 != cone
->Dimension
) {
448 Polyhedron
*tmp
= cone
;
449 cone
= triangularize_cone(cone
, MaxRays
);
450 Polyhedron_Free(tmp
);
452 Polyhedron
*polpos
= NULL
, *polneg
= NULL
;
453 *npos
= 0; *nneg
= 0;
454 for (Polyhedron
*Polar
= cone
; Polar
; Polar
= Polar
->next
)
455 barvinok_decompose(Polar
, &polpos
, &polneg
);
458 for (Polyhedron
*i
= polpos
; i
; i
= i
->next
) {
459 Polyhedron_Polarize(i
);
463 for (Polyhedron
*i
= polneg
; i
; i
= i
->next
) {
464 Polyhedron_Polarize(i
);
475 const int MAX_TRY
=10;
477 * Searches for a vector that is not othogonal to any
478 * of the rays in rays.
480 static void nonorthog(mat_ZZ
& rays
, vec_ZZ
& lambda
)
482 int dim
= rays
.NumCols();
484 lambda
.SetLength(dim
);
485 for (int i
= 2; !found
&& i
<= 50*dim
; i
+=4) {
486 for (int j
= 0; j
< MAX_TRY
; ++j
) {
487 for (int k
= 0; k
< dim
; ++k
) {
488 int r
= random_int(i
)+2;
489 int v
= (2*(r
%2)-1) * (r
>> 1);
493 for (; k
< rays
.NumRows(); ++k
)
494 if (lambda
* rays
[k
] == 0)
496 if (k
== rays
.NumRows()) {
505 static void add_rays(mat_ZZ
& rays
, Polyhedron
*i
, int *r
)
507 unsigned dim
= i
->Dimension
;
508 for (int k
= 0; k
< i
->NbRays
; ++k
) {
509 if (!value_zero_p(i
->Ray
[k
][dim
+1]))
511 values2zz(i
->Ray
[k
]+1, rays
[(*r
)++], dim
);
515 void lattice_point(Value
* values
, Polyhedron
*i
, vec_ZZ
& lambda
, ZZ
& num
)
518 unsigned dim
= i
->Dimension
;
519 if(!value_one_p(values
[dim
])) {
520 Matrix
* Rays
= rays(i
);
521 Matrix
*inv
= Matrix_Alloc(Rays
->NbRows
, Rays
->NbColumns
);
522 int ok
= Matrix_Inverse(Rays
, inv
);
526 Vector
*lambda
= Vector_Alloc(dim
+1);
527 Vector_Matrix_Product(values
, inv
, lambda
->p
);
529 for (int j
= 0; j
< dim
; ++j
)
530 mpz_cdiv_q(lambda
->p
[j
], lambda
->p
[j
], lambda
->p
[dim
]);
531 value_set_si(lambda
->p
[dim
], 1);
532 Vector
*A
= Vector_Alloc(dim
+1);
533 Vector_Matrix_Product(lambda
->p
, Rays
, A
->p
);
536 values2zz(A
->p
, vertex
, dim
);
539 values2zz(values
, vertex
, dim
);
541 num
= vertex
* lambda
;
544 static evalue
*term(int param
, ZZ
& c
, Value
*den
= NULL
)
546 evalue
*EP
= new evalue();
548 value_set_si(EP
->d
,0);
549 EP
->x
.p
= new_enode(polynomial
, 2, param
+ 1);
550 evalue_set_si(&EP
->x
.p
->arr
[0], 0, 1);
551 value_init(EP
->x
.p
->arr
[1].x
.n
);
553 value_set_si(EP
->x
.p
->arr
[1].d
, 1);
555 value_assign(EP
->x
.p
->arr
[1].d
, *den
);
556 zz2value(c
, EP
->x
.p
->arr
[1].x
.n
);
560 static void vertex_period(
561 Polyhedron
*i
, vec_ZZ
& lambda
, Matrix
*T
,
562 Value lcm
, int p
, Vector
*val
,
563 evalue
*E
, evalue
* ev
,
566 unsigned nparam
= T
->NbRows
- 1;
567 unsigned dim
= i
->Dimension
;
573 Vector
* values
= Vector_Alloc(dim
+ 1);
574 Vector_Matrix_Product(val
->p
, T
, values
->p
);
575 value_assign(values
->p
[dim
], lcm
);
576 lattice_point(values
->p
, i
, lambda
, num
);
581 zz2value(num
, ev
->x
.n
);
582 value_assign(ev
->d
, lcm
);
589 values2zz(T
->p
[p
], vertex
, dim
);
590 nump
= vertex
* lambda
;
591 if (First_Non_Zero(val
->p
, p
) == -1) {
592 value_assign(tmp
, lcm
);
593 evalue
*ET
= term(p
, nump
, &tmp
);
595 free_evalue_refs(ET
);
599 value_assign(tmp
, lcm
);
600 if (First_Non_Zero(T
->p
[p
], dim
) != -1)
601 Vector_Gcd(T
->p
[p
], dim
, &tmp
);
603 if (value_lt(tmp
, lcm
)) {
606 value_division(tmp
, lcm
, tmp
);
607 value_set_si(ev
->d
, 0);
608 ev
->x
.p
= new_enode(periodic
, VALUE_TO_INT(tmp
), p
+1);
609 value2zz(tmp
, count
);
611 value_decrement(tmp
, tmp
);
613 ZZ new_offset
= offset
- count
* nump
;
614 value_assign(val
->p
[p
], tmp
);
615 vertex_period(i
, lambda
, T
, lcm
, p
+1, val
, E
,
616 &ev
->x
.p
->arr
[VALUE_TO_INT(tmp
)], new_offset
);
617 } while (value_pos_p(tmp
));
619 vertex_period(i
, lambda
, T
, lcm
, p
+1, val
, E
, ev
, offset
);
623 static void mask_r(Matrix
*f
, int nr
, Vector
*lcm
, int p
, Vector
*val
, evalue
*ev
)
625 unsigned nparam
= lcm
->Size
;
628 Vector
* prod
= Vector_Alloc(f
->NbRows
);
629 Matrix_Vector_Product(f
, val
->p
, prod
->p
);
631 for (int i
= 0; i
< nr
; ++i
) {
632 value_modulus(prod
->p
[i
], prod
->p
[i
], f
->p
[i
][nparam
+1]);
633 isint
&= value_zero_p(prod
->p
[i
]);
635 value_set_si(ev
->d
, 1);
637 value_set_si(ev
->x
.n
, isint
);
644 if (value_one_p(lcm
->p
[p
]))
645 mask_r(f
, nr
, lcm
, p
+1, val
, ev
);
647 value_assign(tmp
, lcm
->p
[p
]);
648 value_set_si(ev
->d
, 0);
649 ev
->x
.p
= new_enode(periodic
, VALUE_TO_INT(tmp
), p
+1);
651 value_decrement(tmp
, tmp
);
652 value_assign(val
->p
[p
], tmp
);
653 mask_r(f
, nr
, lcm
, p
+1, val
, &ev
->x
.p
->arr
[VALUE_TO_INT(tmp
)]);
654 } while (value_pos_p(tmp
));
659 static evalue
*multi_monom(vec_ZZ
& p
)
661 evalue
*X
= new evalue();
664 unsigned nparam
= p
.length()-1;
665 zz2value(p
[nparam
], X
->x
.n
);
666 value_set_si(X
->d
, 1);
667 for (int i
= 0; i
< nparam
; ++i
) {
670 evalue
*T
= term(i
, p
[i
]);
679 * Check whether mapping polyhedron P on the affine combination
680 * num yields a range that has a fixed quotient on integer
682 * If zero is true, then we are only interested in the quotient
683 * for the cases where the remainder is zero.
684 * Returns NULL if false and a newly allocated value if true.
686 static Value
*fixed_quotient(Polyhedron
*P
, vec_ZZ
& num
, Value d
, bool zero
)
689 int len
= num
.length();
690 Matrix
*T
= Matrix_Alloc(2, len
);
691 zz2values(num
, T
->p
[0]);
692 value_set_si(T
->p
[1][len
-1], 1);
693 Polyhedron
*I
= Polyhedron_Image(P
, T
, P
->NbConstraints
);
697 for (i
= 0; i
< I
->NbRays
; ++i
)
698 if (value_zero_p(I
->Ray
[i
][2])) {
706 int bounded
= line_minmax(I
, &min
, &max
);
710 mpz_cdiv_q(min
, min
, d
);
712 mpz_fdiv_q(min
, min
, d
);
713 mpz_fdiv_q(max
, max
, d
);
715 if (value_eq(min
, max
)) {
718 value_assign(*ret
, min
);
726 * Normalize linear expression coef modulo m
727 * Removes common factor and reduces coefficients
728 * Returns index of first non-zero coefficient or len
730 static int normal_mod(Value
*coef
, int len
, Value
*m
)
735 Vector_Gcd(coef
, len
, &gcd
);
737 Vector_AntiScale(coef
, coef
, gcd
, len
);
739 value_division(*m
, *m
, gcd
);
746 for (j
= 0; j
< len
; ++j
)
747 mpz_fdiv_r(coef
[j
], coef
[j
], *m
);
748 for (j
= 0; j
< len
; ++j
)
749 if (value_notzero_p(coef
[j
]))
756 static void mask(Matrix
*f
, evalue
*factor
)
758 int nr
= f
->NbRows
, nc
= f
->NbColumns
;
761 for (n
= 0; n
< nr
&& value_notzero_p(f
->p
[n
][nc
-1]); ++n
)
762 if (value_notone_p(f
->p
[n
][nc
-1]) &&
763 value_notmone_p(f
->p
[n
][nc
-1]))
777 value_set_si(EV
.x
.n
, 1);
779 for (n
= 0; n
< nr
; ++n
) {
780 value_assign(m
, f
->p
[n
][nc
-1]);
781 if (value_one_p(m
) || value_mone_p(m
))
784 int j
= normal_mod(f
->p
[n
], nc
-1, &m
);
786 free_evalue_refs(factor
);
787 value_init(factor
->d
);
788 evalue_set_si(factor
, 0, 1);
792 values2zz(f
->p
[n
], row
, nc
-1);
795 if (j
< (nc
-1)-1 && row
[j
] > g
/2) {
796 for (int k
= j
; k
< (nc
-1); ++k
)
802 value_set_si(EP
.d
, 0);
803 EP
.x
.p
= new_enode(relation
, 2, 0);
804 value_clear(EP
.x
.p
->arr
[1].d
);
805 EP
.x
.p
->arr
[1] = *factor
;
806 evalue
*ev
= &EP
.x
.p
->arr
[0];
807 value_set_si(ev
->d
, 0);
808 ev
->x
.p
= new_enode(fractional
, 3, -1);
809 evalue_set_si(&ev
->x
.p
->arr
[1], 0, 1);
810 evalue_set_si(&ev
->x
.p
->arr
[2], 1, 1);
811 evalue
*E
= multi_monom(row
);
812 value_assign(EV
.d
, m
);
814 value_clear(ev
->x
.p
->arr
[0].d
);
815 ev
->x
.p
->arr
[0] = *E
;
821 free_evalue_refs(&EV
);
827 static void mask(Matrix
*f
, evalue
*factor
)
829 int nr
= f
->NbRows
, nc
= f
->NbColumns
;
832 for (n
= 0; n
< nr
&& value_notzero_p(f
->p
[n
][nc
-1]); ++n
)
833 if (value_notone_p(f
->p
[n
][nc
-1]) &&
834 value_notmone_p(f
->p
[n
][nc
-1]))
842 unsigned np
= nc
- 2;
843 Vector
*lcm
= Vector_Alloc(np
);
844 Vector
*val
= Vector_Alloc(nc
);
845 Vector_Set(val
->p
, 0, nc
);
846 value_set_si(val
->p
[np
], 1);
847 Vector_Set(lcm
->p
, 1, np
);
848 for (n
= 0; n
< nr
; ++n
) {
849 if (value_one_p(f
->p
[n
][nc
-1]) ||
850 value_mone_p(f
->p
[n
][nc
-1]))
852 for (int j
= 0; j
< np
; ++j
)
853 if (value_notzero_p(f
->p
[n
][j
])) {
854 Gcd(f
->p
[n
][j
], f
->p
[n
][nc
-1], &tmp
);
855 value_division(tmp
, f
->p
[n
][nc
-1], tmp
);
856 value_lcm(tmp
, lcm
->p
[j
], &lcm
->p
[j
]);
861 mask_r(f
, nr
, lcm
, 0, val
, &EP
);
866 free_evalue_refs(&EP
);
877 static bool mod_needed(Polyhedron
*PD
, vec_ZZ
& num
, Value d
, evalue
*E
)
879 Value
*q
= fixed_quotient(PD
, num
, d
, false);
884 value_oppose(*q
, *q
);
887 value_set_si(EV
.d
, 1);
889 value_multiply(EV
.x
.n
, *q
, d
);
891 free_evalue_refs(&EV
);
897 static void ceil_mod(Value
*coef
, int len
, Value d
, ZZ
& f
, evalue
*EP
, Polyhedron
*PD
)
903 Vector_Scale(coef
, coef
, m
, len
);
906 int j
= normal_mod(coef
, len
, &m
);
914 values2zz(coef
, num
, len
);
921 evalue_set_si(&tmp
, 0, 1);
925 while (j
< len
-1 && (num
[j
] == g
/2 || num
[j
] == 0))
927 if ((j
< len
-1 && num
[j
] > g
/2) || (j
== len
-1 && num
[j
] >= (g
+1)/2)) {
928 for (int k
= j
; k
< len
-1; ++k
)
931 num
[len
-1] = g
- 1 - num
[len
-1];
932 value_assign(tmp
.d
, m
);
934 zz2value(t
, tmp
.x
.n
);
940 ZZ t
= num
[len
-1] * f
;
941 zz2value(t
, tmp
.x
.n
);
942 value_assign(tmp
.d
, m
);
945 evalue
*E
= multi_monom(num
);
949 if (PD
&& !mod_needed(PD
, num
, m
, E
)) {
952 value_assign(EV
.d
, m
);
957 value_set_si(EV
.x
.n
, 1);
958 value_assign(EV
.d
, m
);
961 value_set_si(EV
.d
, 0);
962 EV
.x
.p
= new_enode(fractional
, 3, -1);
963 evalue_copy(&EV
.x
.p
->arr
[0], E
);
964 evalue_set_si(&EV
.x
.p
->arr
[1], 0, 1);
965 value_init(EV
.x
.p
->arr
[2].x
.n
);
966 zz2value(f
, EV
.x
.p
->arr
[2].x
.n
);
967 value_set_si(EV
.x
.p
->arr
[2].d
, 1);
972 free_evalue_refs(&EV
);
977 free_evalue_refs(&tmp
);
983 evalue
* bv_ceil3(Value
*coef
, int len
, Value d
, Polyhedron
*P
)
985 Vector
*val
= Vector_Alloc(len
);
990 Vector_Scale(coef
, val
->p
, t
, len
);
991 value_absolute(t
, d
);
994 values2zz(val
->p
, num
, len
);
995 evalue
*EP
= multi_monom(num
);
1000 value_set_si(tmp
.x
.n
, 1);
1001 value_assign(tmp
.d
, t
);
1007 ceil_mod(val
->p
, len
, t
, one
, EP
, P
);
1010 /* copy EP to malloc'ed evalue */
1016 free_evalue_refs(&tmp
);
1023 evalue
* lattice_point(
1024 Polyhedron
*i
, vec_ZZ
& lambda
, Matrix
*W
, Value lcm
, Polyhedron
*PD
)
1026 unsigned nparam
= W
->NbColumns
- 1;
1028 Matrix
* Rays
= rays2(i
);
1029 Matrix
*T
= Transpose(Rays
);
1030 Matrix
*T2
= Matrix_Copy(T
);
1031 Matrix
*inv
= Matrix_Alloc(T2
->NbRows
, T2
->NbColumns
);
1032 int ok
= Matrix_Inverse(T2
, inv
);
1037 matrix2zz(W
, vertex
, W
->NbRows
, W
->NbColumns
);
1040 num
= lambda
* vertex
;
1042 evalue
*EP
= multi_monom(num
);
1046 value_init(tmp
.x
.n
);
1047 value_set_si(tmp
.x
.n
, 1);
1048 value_assign(tmp
.d
, lcm
);
1052 Matrix
*L
= Matrix_Alloc(inv
->NbRows
, W
->NbColumns
);
1053 Matrix_Product(inv
, W
, L
);
1056 matrix2zz(T
, RT
, T
->NbRows
, T
->NbColumns
);
1059 vec_ZZ p
= lambda
* RT
;
1061 for (int i
= 0; i
< L
->NbRows
; ++i
) {
1062 ceil_mod(L
->p
[i
], nparam
+1, lcm
, p
[i
], EP
, PD
);
1068 free_evalue_refs(&tmp
);
1072 evalue
* lattice_point(
1073 Polyhedron
*i
, vec_ZZ
& lambda
, Matrix
*W
, Value lcm
, Polyhedron
*PD
)
1075 Matrix
*T
= Transpose(W
);
1076 unsigned nparam
= T
->NbRows
- 1;
1078 evalue
*EP
= new evalue();
1080 evalue_set_si(EP
, 0, 1);
1083 Vector
*val
= Vector_Alloc(nparam
+1);
1084 value_set_si(val
->p
[nparam
], 1);
1085 ZZ
offset(INIT_VAL
, 0);
1087 vertex_period(i
, lambda
, T
, lcm
, 0, val
, EP
, &ev
, offset
);
1090 free_evalue_refs(&ev
);
1101 Param_Vertices
* V
, Polyhedron
*i
, vec_ZZ
& lambda
, term_info
* term
,
1104 unsigned nparam
= V
->Vertex
->NbColumns
- 2;
1105 unsigned dim
= i
->Dimension
;
1107 vertex
.SetDims(V
->Vertex
->NbRows
, nparam
+1);
1111 value_set_si(lcm
, 1);
1112 for (int j
= 0; j
< V
->Vertex
->NbRows
; ++j
) {
1113 value_lcm(lcm
, V
->Vertex
->p
[j
][nparam
+1], &lcm
);
1115 if (value_notone_p(lcm
)) {
1116 Matrix
* mv
= Matrix_Alloc(dim
, nparam
+1);
1117 for (int j
= 0 ; j
< dim
; ++j
) {
1118 value_division(tmp
, lcm
, V
->Vertex
->p
[j
][nparam
+1]);
1119 Vector_Scale(V
->Vertex
->p
[j
], mv
->p
[j
], tmp
, nparam
+1);
1122 term
->E
= lattice_point(i
, lambda
, mv
, lcm
, PD
);
1130 for (int i
= 0; i
< V
->Vertex
->NbRows
; ++i
) {
1131 assert(value_one_p(V
->Vertex
->p
[i
][nparam
+1])); // for now
1132 values2zz(V
->Vertex
->p
[i
], vertex
[i
], nparam
+1);
1136 num
= lambda
* vertex
;
1140 for (int j
= 0; j
< nparam
; ++j
)
1146 term
->E
= multi_monom(num
);
1150 term
->constant
= num
[nparam
];
1153 term
->coeff
= num
[p
];
1160 void normalize(Polyhedron
*i
, vec_ZZ
& lambda
, ZZ
& sign
, ZZ
& num
, vec_ZZ
& den
)
1162 unsigned dim
= i
->Dimension
;
1166 rays
.SetDims(dim
, dim
);
1167 add_rays(rays
, i
, &r
);
1168 den
= rays
* lambda
;
1171 for (int j
= 0; j
< den
.length(); ++j
) {
1175 den
[j
] = abs(den
[j
]);
1183 void barvinok_count(Polyhedron
*P
, Value
* result
, unsigned NbMaxCons
)
1185 Polyhedron
** vcone
;
1188 sign
.SetLength(ncone
);
1196 value_set_si(*result
, 0);
1200 for (; r
< P
->NbRays
; ++r
)
1201 if (value_zero_p(P
->Ray
[r
][P
->Dimension
+1]))
1203 if (P
->NbBid
!=0 || r
< P
->NbRays
) {
1204 value_set_si(*result
, -1);
1208 P
= remove_equalities(P
);
1211 value_set_si(*result
, 0);
1217 value_set_si(factor
, 1);
1218 Q
= Polyhedron_Reduce(P
, &factor
);
1225 if (P
->Dimension
== 0) {
1226 value_assign(*result
, factor
);
1229 value_clear(factor
);
1234 vcone
= new (Polyhedron
*)[P
->NbRays
];
1236 for (int j
= 0; j
< P
->NbRays
; ++j
) {
1238 Polyhedron
*C
= supporting_cone(P
, j
);
1239 decompose(C
, &vcone
[j
], &npos
, &nneg
, NbMaxCons
);
1240 ncone
+= npos
+ nneg
;
1241 sign
.SetLength(ncone
);
1242 for (int k
= 0; k
< npos
; ++k
)
1243 sign
[ncone
-nneg
-k
-1] = 1;
1244 for (int k
= 0; k
< nneg
; ++k
)
1245 sign
[ncone
-k
-1] = -1;
1249 rays
.SetDims(ncone
* dim
, dim
);
1251 for (int j
= 0; j
< P
->NbRays
; ++j
) {
1252 for (Polyhedron
*i
= vcone
[j
]; i
; i
= i
->next
) {
1253 assert(i
->NbRays
-1 == dim
);
1254 add_rays(rays
, i
, &r
);
1258 nonorthog(rays
, lambda
);
1262 num
.SetLength(ncone
);
1263 den
.SetDims(ncone
,dim
);
1266 for (int j
= 0; j
< P
->NbRays
; ++j
) {
1267 for (Polyhedron
*i
= vcone
[j
]; i
; i
= i
->next
) {
1268 lattice_point(P
->Ray
[j
]+1, i
, lambda
, num
[f
]);
1269 normalize(i
, lambda
, sign
[f
], num
[f
], den
[f
]);
1274 for (int j
= 1; j
< num
.length(); ++j
)
1277 for (int j
= 0; j
< num
.length(); ++j
)
1283 for (int j
= 0; j
< P
->NbRays
; ++j
) {
1284 for (Polyhedron
*i
= vcone
[j
]; i
; i
= i
->next
) {
1285 dpoly
d(dim
, num
[f
]);
1286 dpoly
n(dim
, den
[f
][0], 1);
1287 for (int k
= 1; k
< dim
; ++k
) {
1288 dpoly
fact(dim
, den
[f
][k
], 1);
1291 d
.div(n
, count
, sign
[f
]);
1295 assert(value_one_p(&count
[0]._mp_den
));
1296 value_multiply(*result
, &count
[0]._mp_num
, factor
);
1299 for (int j
= 0; j
< P
->NbRays
; ++j
)
1300 Domain_Free(vcone
[j
]);
1306 value_clear(factor
);
1309 static void uni_polynom(int param
, Vector
*c
, evalue
*EP
)
1311 unsigned dim
= c
->Size
-2;
1313 value_set_si(EP
->d
,0);
1314 EP
->x
.p
= new_enode(polynomial
, dim
+1, param
+1);
1315 for (int j
= 0; j
<= dim
; ++j
)
1316 evalue_set(&EP
->x
.p
->arr
[j
], c
->p
[j
], c
->p
[dim
+1]);
1319 static void multi_polynom(Vector
*c
, evalue
* X
, evalue
*EP
)
1321 unsigned dim
= c
->Size
-2;
1325 evalue_set(&EC
, c
->p
[dim
], c
->p
[dim
+1]);
1328 evalue_set(EP
, c
->p
[dim
], c
->p
[dim
+1]);
1330 for (int i
= dim
-1; i
>= 0; --i
) {
1332 value_assign(EC
.x
.n
, c
->p
[i
]);
1335 free_evalue_refs(&EC
);
1338 Polyhedron
*unfringe (Polyhedron
*P
, unsigned MaxRays
)
1340 int len
= P
->Dimension
+2;
1341 Polyhedron
*T
, *R
= P
;
1344 Vector
*row
= Vector_Alloc(len
);
1345 value_set_si(row
->p
[0], 1);
1347 R
= DomainConstraintSimplify(Polyhedron_Copy(P
), MaxRays
);
1349 Matrix
*M
= Matrix_Alloc(2, len
-1);
1350 value_set_si(M
->p
[1][len
-2], 1);
1351 for (int v
= 0; v
< P
->Dimension
; ++v
) {
1352 value_set_si(M
->p
[0][v
], 1);
1353 Polyhedron
*I
= Polyhedron_Image(P
, M
, 2+1);
1354 value_set_si(M
->p
[0][v
], 0);
1355 for (int r
= 0; r
< I
->NbConstraints
; ++r
) {
1356 if (value_zero_p(I
->Constraint
[r
][0]))
1358 if (value_zero_p(I
->Constraint
[r
][1]))
1360 if (value_one_p(I
->Constraint
[r
][1]))
1362 if (value_mone_p(I
->Constraint
[r
][1]))
1364 value_absolute(g
, I
->Constraint
[r
][1]);
1365 Vector_Set(row
->p
+1, 0, len
-2);
1366 value_division(row
->p
[1+v
], I
->Constraint
[r
][1], g
);
1367 mpz_fdiv_q(row
->p
[len
-1], I
->Constraint
[r
][2], g
);
1369 R
= AddConstraints(row
->p
, 1, R
, MaxRays
);
1381 static Polyhedron
*reduce_domain(Polyhedron
*D
, Matrix
*CT
, Polyhedron
*CEq
,
1382 Polyhedron
**fVD
, int nd
, unsigned MaxRays
)
1387 Dt
= CT
? DomainPreimage(D
, CT
, MaxRays
) : D
;
1388 Polyhedron
*rVD
= DomainIntersection(Dt
, CEq
, MaxRays
);
1390 /* if rVD is empty or too small in geometric dimension */
1391 if(!rVD
|| emptyQ(rVD
) ||
1392 (rVD
->Dimension
-rVD
->NbEq
< Dt
->Dimension
-Dt
->NbEq
-CEq
->NbEq
)) {
1397 return 0; /* empty validity domain */
1403 fVD
[nd
] = Domain_Copy(rVD
);
1404 for (int i
= 0 ; i
< nd
; ++i
) {
1405 Polyhedron
*I
= DomainIntersection(fVD
[nd
], fVD
[i
], MaxRays
);
1410 Polyhedron
*F
= DomainSimplify(I
, fVD
[nd
], MaxRays
);
1412 Polyhedron
*T
= rVD
;
1413 rVD
= DomainDifference(rVD
, F
, MaxRays
);
1420 rVD
= DomainConstraintSimplify(rVD
, MaxRays
);
1422 Domain_Free(fVD
[nd
]);
1429 barvinok_count(rVD
, &c
, MaxRays
);
1430 if (value_zero_p(c
)) {
1439 evalue
* barvinok_enumerate_ev(Polyhedron
*P
, Polyhedron
* C
, unsigned MaxRays
)
1441 //P = unfringe(P, MaxRays);
1442 Polyhedron
*CEq
= NULL
, *rVD
, *pVD
, *CA
;
1444 Param_Polyhedron
*PP
= NULL
;
1445 Param_Domain
*D
, *next
;
1448 unsigned nparam
= C
->Dimension
;
1451 value_init(eres
->d
);
1452 value_set_si(eres
->d
, 0);
1455 value_init(factor
.d
);
1456 evalue_set_si(&factor
, 1, 1);
1458 CA
= align_context(C
, P
->Dimension
, MaxRays
);
1459 P
= DomainIntersection(P
, CA
, MaxRays
);
1460 Polyhedron_Free(CA
);
1462 if (C
->Dimension
== 0 || emptyQ(P
)) {
1464 eres
->x
.p
= new_enode(partition
, 2, C
->Dimension
);
1465 EVALUE_SET_DOMAIN(eres
->x
.p
->arr
[0],
1466 DomainConstraintSimplify(CEq
? CEq
: Polyhedron_Copy(C
), MaxRays
));
1467 value_set_si(eres
->x
.p
->arr
[1].d
, 1);
1468 value_init(eres
->x
.p
->arr
[1].x
.n
);
1470 value_set_si(eres
->x
.p
->arr
[1].x
.n
, 0);
1472 barvinok_count(P
, &eres
->x
.p
->arr
[1].x
.n
, MaxRays
);
1474 emul(&factor
, eres
);
1475 reduce_evalue(eres
);
1476 free_evalue_refs(&factor
);
1481 Param_Polyhedron_Free(PP
);
1485 for (r
= 0; r
< P
->NbRays
; ++r
)
1486 if (value_zero_p(P
->Ray
[r
][0]) ||
1487 value_zero_p(P
->Ray
[r
][P
->Dimension
+1])) {
1489 for (i
= P
->Dimension
- nparam
; i
< P
->Dimension
; ++i
)
1490 if (value_notzero_p(P
->Ray
[r
][i
+1]))
1492 if (i
>= P
->Dimension
)
1500 P
= remove_equalities_p(P
, P
->Dimension
-nparam
, &f
);
1504 if (P
->Dimension
== nparam
) {
1506 P
= Universe_Polyhedron(0);
1510 Polyhedron
*Q
= ParamPolyhedron_Reduce(P
, P
->Dimension
-nparam
, &factor
);
1513 if (Q
->Dimension
== nparam
) {
1515 P
= Universe_Polyhedron(0);
1520 Polyhedron
*oldP
= P
;
1521 PP
= Polyhedron2Param_SimplifiedDomain(&P
,C
,MaxRays
,&CEq
,&CT
);
1523 Polyhedron_Free(oldP
);
1525 if (isIdentity(CT
)) {
1529 assert(CT
->NbRows
!= CT
->NbColumns
);
1530 if (CT
->NbRows
== 1) // no more parameters
1532 nparam
= CT
->NbRows
- 1;
1535 unsigned dim
= P
->Dimension
- nparam
;
1536 Polyhedron
** vcone
= new (Polyhedron
*)[PP
->nbV
];
1537 int * npos
= new int[PP
->nbV
];
1538 int * nneg
= new int[PP
->nbV
];
1542 for (i
= 0, V
= PP
->V
; V
; ++i
, V
= V
->next
) {
1543 Polyhedron
*C
= supporting_cone_p(P
, V
);
1544 decompose(C
, &vcone
[i
], &npos
[i
], &nneg
[i
], MaxRays
);
1547 Vector
*c
= Vector_Alloc(dim
+2);
1550 for (nd
= 0, D
=PP
->D
; D
; ++nd
, D
=D
->next
);
1551 struct section
{ Polyhedron
*D
; evalue E
; };
1552 section
*s
= new section
[nd
];
1553 Polyhedron
**fVD
= new (Polyhedron
*)[nd
];
1555 for(nd
= 0, D
=PP
->D
; D
; D
=next
) {
1558 Polyhedron
*rVD
= reduce_domain(D
->Domain
, CT
, CEq
,
1563 pVD
= CT
? DomainImage(rVD
,CT
,MaxRays
) : rVD
;
1566 sign
.SetLength(ncone
);
1567 FORALL_PVertex_in_ParamPolyhedron(V
,D
,PP
) // _i is internal counter
1568 ncone
+= npos
[_i
] + nneg
[_i
];
1569 sign
.SetLength(ncone
);
1570 for (int k
= 0; k
< npos
[_i
]; ++k
)
1571 sign
[ncone
-nneg
[_i
]-k
-1] = 1;
1572 for (int k
= 0; k
< nneg
[_i
]; ++k
)
1573 sign
[ncone
-k
-1] = -1;
1574 END_FORALL_PVertex_in_ParamPolyhedron
;
1577 rays
.SetDims(ncone
* dim
, dim
);
1579 FORALL_PVertex_in_ParamPolyhedron(V
,D
,PP
) // _i is internal counter
1580 for (Polyhedron
*i
= vcone
[_i
]; i
; i
= i
->next
) {
1581 assert(i
->NbRays
-1 == dim
);
1582 add_rays(rays
, i
, &r
);
1584 END_FORALL_PVertex_in_ParamPolyhedron
;
1586 nonorthog(rays
, lambda
);
1589 den
.SetDims(ncone
,dim
);
1590 term_info
*num
= new term_info
[ncone
];
1593 FORALL_PVertex_in_ParamPolyhedron(V
,D
,PP
)
1594 for (Polyhedron
*i
= vcone
[_i
]; i
; i
= i
->next
) {
1595 lattice_point(V
, i
, lambda
, &num
[f
], pVD
);
1596 normalize(i
, lambda
, sign
[f
], num
[f
].constant
, den
[f
]);
1599 END_FORALL_PVertex_in_ParamPolyhedron
;
1600 ZZ min
= num
[0].constant
;
1601 for (int j
= 1; j
< ncone
; ++j
)
1602 if (num
[j
].constant
< min
)
1603 min
= num
[j
].constant
;
1604 for (int j
= 0; j
< ncone
; ++j
)
1605 num
[j
].constant
-= min
;
1607 value_init(s
[nd
].E
.d
);
1608 evalue_set_si(&s
[nd
].E
, 0, 1);
1611 FORALL_PVertex_in_ParamPolyhedron(V
,D
,PP
)
1612 for (Polyhedron
*i
= vcone
[_i
]; i
; i
= i
->next
) {
1613 dpoly
n(dim
, den
[f
][0], 1);
1614 for (int k
= 1; k
< dim
; ++k
) {
1615 dpoly
fact(dim
, den
[f
][k
], 1);
1618 if (num
[f
].E
!= NULL
) {
1619 ZZ
one(INIT_VAL
, 1);
1620 dpoly_n
d(dim
, num
[f
].constant
, one
);
1621 d
.div(n
, c
, sign
[f
]);
1623 multi_polynom(c
, num
[f
].E
, &EV
);
1624 eadd(&EV
, &s
[nd
].E
);
1625 free_evalue_refs(&EV
);
1626 free_evalue_refs(num
[f
].E
);
1628 } else if (num
[f
].pos
!= -1) {
1629 dpoly_n
d(dim
, num
[f
].constant
, num
[f
].coeff
);
1630 d
.div(n
, c
, sign
[f
]);
1632 uni_polynom(num
[f
].pos
, c
, &EV
);
1633 eadd(&EV
, &s
[nd
].E
);
1634 free_evalue_refs(&EV
);
1636 mpq_set_si(count
, 0, 1);
1637 dpoly
d(dim
, num
[f
].constant
);
1638 d
.div(n
, count
, sign
[f
]);
1641 evalue_set(&EV
, &count
[0]._mp_num
, &count
[0]._mp_den
);
1642 eadd(&EV
, &s
[nd
].E
);
1643 free_evalue_refs(&EV
);
1647 END_FORALL_PVertex_in_ParamPolyhedron
;
1653 addeliminatedparams_evalue(&s
[nd
].E
, CT
);
1661 evalue_set_si(eres
, 0, 1);
1663 eres
->x
.p
= new_enode(partition
, 2*nd
, C
->Dimension
);
1664 for (int j
= 0; j
< nd
; ++j
) {
1665 EVALUE_SET_DOMAIN(eres
->x
.p
->arr
[2*j
], s
[j
].D
);
1666 value_clear(eres
->x
.p
->arr
[2*j
+1].d
);
1667 eres
->x
.p
->arr
[2*j
+1] = s
[j
].E
;
1668 Domain_Free(fVD
[j
]);
1676 for (int j
= 0; j
< PP
->nbV
; ++j
)
1677 Domain_Free(vcone
[j
]);
1683 Polyhedron_Free(CEq
);
1688 Enumeration
* barvinok_enumerate(Polyhedron
*P
, Polyhedron
* C
, unsigned MaxRays
)
1690 evalue
*EP
= barvinok_enumerate_ev(P
, C
, MaxRays
);
1692 return partition2enumeration(EP
);
1695 static void SwapColumns(Value
**V
, int n
, int i
, int j
)
1697 for (int r
= 0; r
< n
; ++r
)
1698 value_swap(V
[r
][i
], V
[r
][j
]);
1701 static void SwapColumns(Polyhedron
*P
, int i
, int j
)
1703 SwapColumns(P
->Constraint
, P
->NbConstraints
, i
, j
);
1704 SwapColumns(P
->Ray
, P
->NbRays
, i
, j
);
1707 static bool SplitOnConstraint(Polyhedron
*P
, int i
, int l
, int u
,
1708 int nvar
, int len
, int exist
, int MaxRays
,
1709 Vector
*row
, Value
& f
, bool independent
,
1710 Polyhedron
**pos
, Polyhedron
**neg
)
1712 value_oppose(f
, P
->Constraint
[u
][nvar
+i
+1]);
1713 Vector_Combine(P
->Constraint
[l
]+1, P
->Constraint
[u
]+1,
1715 f
, P
->Constraint
[l
][nvar
+i
+1], len
-1);
1717 //printf("l: %d, u: %d\n", l, u);
1718 value_multiply(f
, f
, P
->Constraint
[l
][nvar
+i
+1]);
1719 value_substract(row
->p
[len
-1], row
->p
[len
-1], f
);
1720 value_set_si(f
, -1);
1721 Vector_Scale(row
->p
+1, row
->p
+1, f
, len
-1);
1722 value_decrement(row
->p
[len
-1], row
->p
[len
-1]);
1723 Vector_Gcd(row
->p
+1, len
- 2, &f
);
1724 if (value_notone_p(f
)) {
1725 Vector_AntiScale(row
->p
+1, row
->p
+1, f
, len
-2);
1726 mpz_fdiv_q(row
->p
[len
-1], row
->p
[len
-1], f
);
1728 *neg
= AddConstraints(row
->p
, 1, P
, MaxRays
);
1730 /* We found an independent, but useless constraint
1731 * Maybe we should detect this earlier and not
1732 * mark the variable as INDEPENDENT
1734 if (emptyQ((*neg
))) {
1735 Polyhedron_Free(*neg
);
1739 value_set_si(f
, -1);
1740 Vector_Scale(row
->p
+1, row
->p
+1, f
, len
-1);
1741 value_decrement(row
->p
[len
-1], row
->p
[len
-1]);
1742 *pos
= AddConstraints(row
->p
, 1, P
, MaxRays
);
1744 if (emptyQ((*pos
))) {
1745 Polyhedron_Free(*neg
);
1746 Polyhedron_Free(*pos
);
1754 * unimodularly transform P such that constraint r is transformed
1755 * into a constraint that involves only a single (the first)
1756 * existential variable
1759 static Polyhedron
*rotate_along(Polyhedron
*P
, int r
, int nvar
, int exist
,
1765 Vector
*row
= Vector_Alloc(exist
);
1766 Vector_Copy(P
->Constraint
[r
]+1+nvar
, row
->p
, exist
);
1767 Vector_Gcd(row
->p
, exist
, &g
);
1768 if (value_notone_p(g
))
1769 Vector_AntiScale(row
->p
, row
->p
, g
, exist
);
1772 Matrix
*M
= unimodular_complete(row
);
1773 Matrix
*M2
= Matrix_Alloc(P
->Dimension
+1, P
->Dimension
+1);
1774 for (r
= 0; r
< nvar
; ++r
)
1775 value_set_si(M2
->p
[r
][r
], 1);
1776 for ( ; r
< nvar
+exist
; ++r
)
1777 Vector_Copy(M
->p
[r
-nvar
], M2
->p
[r
]+nvar
, exist
);
1778 for ( ; r
< P
->Dimension
+1; ++r
)
1779 value_set_si(M2
->p
[r
][r
], 1);
1780 Polyhedron
*T
= Polyhedron_Image(P
, M2
, MaxRays
);
1789 static bool SplitOnVar(Polyhedron
*P
, int i
,
1790 int nvar
, int len
, int exist
, int MaxRays
,
1791 Vector
*row
, Value
& f
, bool independent
,
1792 Polyhedron
**pos
, Polyhedron
**neg
)
1796 for (int l
= P
->NbEq
; l
< P
->NbConstraints
; ++l
) {
1797 if (value_negz_p(P
->Constraint
[l
][nvar
+i
+1]))
1801 for (j
= 0; j
< exist
; ++j
)
1802 if (j
!= i
&& value_notzero_p(P
->Constraint
[l
][nvar
+j
+1]))
1808 for (int u
= P
->NbEq
; u
< P
->NbConstraints
; ++u
) {
1809 if (value_posz_p(P
->Constraint
[u
][nvar
+i
+1]))
1813 for (j
= 0; j
< exist
; ++j
)
1814 if (j
!= i
&& value_notzero_p(P
->Constraint
[u
][nvar
+j
+1]))
1820 if (SplitOnConstraint(P
, i
, l
, u
,
1821 nvar
, len
, exist
, MaxRays
,
1822 row
, f
, independent
,
1826 SwapColumns(*neg
, nvar
+1, nvar
+1+i
);
1836 static bool double_bound_pair(Polyhedron
*P
, int nvar
, int exist
,
1837 int i
, int l1
, int l2
,
1838 Polyhedron
**pos
, Polyhedron
**neg
)
1842 Vector
*row
= Vector_Alloc(P
->Dimension
+2);
1843 value_set_si(row
->p
[0], 1);
1844 value_oppose(f
, P
->Constraint
[l1
][nvar
+i
+1]);
1845 Vector_Combine(P
->Constraint
[l1
]+1, P
->Constraint
[l2
]+1,
1847 P
->Constraint
[l2
][nvar
+i
+1], f
,
1849 ConstraintSimplify(row
->p
, row
->p
, P
->Dimension
+2, &f
);
1850 *pos
= AddConstraints(row
->p
, 1, P
, 0);
1851 value_set_si(f
, -1);
1852 Vector_Scale(row
->p
+1, row
->p
+1, f
, P
->Dimension
+1);
1853 value_decrement(row
->p
[P
->Dimension
+1], row
->p
[P
->Dimension
+1]);
1854 *neg
= AddConstraints(row
->p
, 1, P
, 0);
1858 return !emptyQ((*pos
)) && !emptyQ((*neg
));
1861 static bool double_bound(Polyhedron
*P
, int nvar
, int exist
,
1862 Polyhedron
**pos
, Polyhedron
**neg
)
1864 for (int i
= 0; i
< exist
; ++i
) {
1866 for (l1
= P
->NbEq
; l1
< P
->NbConstraints
; ++l1
) {
1867 if (value_negz_p(P
->Constraint
[l1
][nvar
+i
+1]))
1869 for (l2
= l1
+ 1; l2
< P
->NbConstraints
; ++l2
) {
1870 if (value_negz_p(P
->Constraint
[l2
][nvar
+i
+1]))
1872 if (double_bound_pair(P
, nvar
, exist
, i
, l1
, l2
, pos
, neg
))
1876 for (l1
= P
->NbEq
; l1
< P
->NbConstraints
; ++l1
) {
1877 if (value_posz_p(P
->Constraint
[l1
][nvar
+i
+1]))
1879 if (l1
< P
->NbConstraints
)
1880 for (l2
= l1
+ 1; l2
< P
->NbConstraints
; ++l2
) {
1881 if (value_posz_p(P
->Constraint
[l2
][nvar
+i
+1]))
1883 if (double_bound_pair(P
, nvar
, exist
, i
, l1
, l2
, pos
, neg
))
1895 INDEPENDENT
= 1 << 2,
1898 static evalue
* enumerate_or(Polyhedron
*pos
, Polyhedron
*neg
,
1899 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
1902 fprintf(stderr
, "\nER: Or\n");
1903 #endif /* DEBUG_ER */
1906 barvinok_enumerate_e(neg
, exist
, nparam
, MaxRays
);
1908 barvinok_enumerate_e(pos
, exist
, nparam
, MaxRays
);
1911 evalue_copy(&E
, EP
);
1914 free_evalue_refs(EN
);
1916 evalue_set_si(EN
, -1, 1);
1920 free_evalue_refs(EN
);
1922 free_evalue_refs(&E
);
1923 Polyhedron_Free(neg
);
1924 Polyhedron_Free(pos
);
1931 static evalue
* enumerate_sum(Polyhedron
*P
,
1932 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
1934 int nvar
= P
->Dimension
- exist
- nparam
;
1935 int toswap
= nvar
< exist
? nvar
: exist
;
1936 for (int i
= 0; i
< toswap
; ++i
)
1937 SwapColumns(P
, 1 + i
, nvar
+exist
- i
);
1941 fprintf(stderr
, "\nER: Sum\n");
1942 #endif /* DEBUG_ER */
1944 evalue
*EP
= barvinok_enumerate_e(P
, exist
, nparam
, MaxRays
);
1946 for (int i
= 0; i
< /* nvar */ nparam
; ++i
) {
1947 Matrix
*C
= Matrix_Alloc(1, 1 + nparam
+ 1);
1948 value_set_si(C
->p
[0][0], 1);
1950 value_init(split
.d
);
1951 value_set_si(split
.d
, 0);
1952 split
.x
.p
= new_enode(partition
, 4, nparam
);
1953 value_set_si(C
->p
[0][1+i
], 1);
1954 Matrix
*C2
= Matrix_Copy(C
);
1955 EVALUE_SET_DOMAIN(split
.x
.p
->arr
[0],
1956 Constraints2Polyhedron(C2
, MaxRays
));
1958 evalue_set_si(&split
.x
.p
->arr
[1], 1, 1);
1959 value_set_si(C
->p
[0][1+i
], -1);
1960 value_set_si(C
->p
[0][1+nparam
], -1);
1961 EVALUE_SET_DOMAIN(split
.x
.p
->arr
[2],
1962 Constraints2Polyhedron(C
, MaxRays
));
1963 evalue_set_si(&split
.x
.p
->arr
[3], 1, 1);
1965 free_evalue_refs(&split
);
1969 evalue_range_reduction(EP
);
1971 evalue_frac2floor(EP
);
1973 evalue
*sum
= esum(EP
, nvar
);
1975 free_evalue_refs(EP
);
1979 evalue_range_reduction(EP
);
1984 static evalue
* split_sure(Polyhedron
*P
, Polyhedron
*S
,
1985 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
1987 int nvar
= P
->Dimension
- exist
- nparam
;
1989 Matrix
*M
= Matrix_Alloc(exist
, S
->Dimension
+2);
1990 for (int i
= 0; i
< exist
; ++i
)
1991 value_set_si(M
->p
[i
][nvar
+i
+1], 1);
1993 S
= DomainAddRays(S
, M
, MaxRays
);
1995 Polyhedron
*F
= DomainAddRays(P
, M
, MaxRays
);
1996 Polyhedron
*D
= DomainDifference(F
, S
, MaxRays
);
1998 D
= Disjoint_Domain(D
, 0, MaxRays
);
2003 M
= Matrix_Alloc(P
->Dimension
+1-exist
, P
->Dimension
+1);
2004 for (int j
= 0; j
< nvar
; ++j
)
2005 value_set_si(M
->p
[j
][j
], 1);
2006 for (int j
= 0; j
< nparam
+1; ++j
)
2007 value_set_si(M
->p
[nvar
+j
][nvar
+exist
+j
], 1);
2008 Polyhedron
*T
= Polyhedron_Image(S
, M
, MaxRays
);
2009 evalue
*EP
= barvinok_enumerate_e(T
, 0, nparam
, MaxRays
);
2014 for (Polyhedron
*Q
= D
; Q
; Q
= Q
->next
) {
2015 Polyhedron
*N
= Q
->next
;
2017 T
= DomainIntersection(P
, Q
, MaxRays
);
2018 evalue
*E
= barvinok_enumerate_e(T
, exist
, nparam
, MaxRays
);
2020 free_evalue_refs(E
);
2029 static evalue
* enumerate_sure(Polyhedron
*P
,
2030 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
2034 int nvar
= P
->Dimension
- exist
- nparam
;
2040 for (i
= 0; i
< exist
; ++i
) {
2041 Matrix
*M
= Matrix_Alloc(S
->NbConstraints
, S
->Dimension
+2);
2043 value_set_si(lcm
, 1);
2044 for (int j
= 0; j
< S
->NbConstraints
; ++j
) {
2045 if (value_negz_p(S
->Constraint
[j
][1+nvar
+i
]))
2047 if (value_one_p(S
->Constraint
[j
][1+nvar
+i
]))
2049 value_lcm(lcm
, S
->Constraint
[j
][1+nvar
+i
], &lcm
);
2052 for (int j
= 0; j
< S
->NbConstraints
; ++j
) {
2053 if (value_negz_p(S
->Constraint
[j
][1+nvar
+i
]))
2055 if (value_one_p(S
->Constraint
[j
][1+nvar
+i
]))
2057 value_division(f
, lcm
, S
->Constraint
[j
][1+nvar
+i
]);
2058 Vector_Scale(S
->Constraint
[j
], M
->p
[c
], f
, S
->Dimension
+2);
2059 value_substract(M
->p
[c
][S
->Dimension
+1],
2060 M
->p
[c
][S
->Dimension
+1],
2062 value_increment(M
->p
[c
][S
->Dimension
+1],
2063 M
->p
[c
][S
->Dimension
+1]);
2067 S
= AddConstraints(M
->p
[0], c
, S
, MaxRays
);
2082 fprintf(stderr
, "\nER: Sure\n");
2083 #endif /* DEBUG_ER */
2085 return split_sure(P
, S
, exist
, nparam
, MaxRays
);
2088 static evalue
* enumerate_sure2(Polyhedron
*P
,
2089 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
2091 int nvar
= P
->Dimension
- exist
- nparam
;
2093 for (r
= 0; r
< P
->NbRays
; ++r
)
2094 if (value_one_p(P
->Ray
[r
][0]) &&
2095 value_one_p(P
->Ray
[r
][P
->Dimension
+1]))
2101 Matrix
*M
= Matrix_Alloc(nvar
+ 1 + nparam
, P
->Dimension
+2);
2102 for (int i
= 0; i
< nvar
; ++i
)
2103 value_set_si(M
->p
[i
][1+i
], 1);
2104 for (int i
= 0; i
< nparam
; ++i
)
2105 value_set_si(M
->p
[i
+nvar
][1+nvar
+exist
+i
], 1);
2106 Vector_Copy(P
->Ray
[r
]+1+nvar
, M
->p
[nvar
+nparam
]+1+nvar
, exist
);
2107 value_set_si(M
->p
[nvar
+nparam
][0], 1);
2108 value_set_si(M
->p
[nvar
+nparam
][P
->Dimension
+1], 1);
2109 Polyhedron
* F
= Rays2Polyhedron(M
, MaxRays
);
2112 Polyhedron
*I
= DomainIntersection(F
, P
, MaxRays
);
2116 fprintf(stderr
, "\nER: Sure2\n");
2117 #endif /* DEBUG_ER */
2119 return split_sure(P
, I
, exist
, nparam
, MaxRays
);
2122 static evalue
* new_zero_ep()
2127 evalue_set_si(EP
, 0, 1);
2131 static evalue
* enumerate_vd(Polyhedron
**PA
,
2132 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
2134 Polyhedron
*P
= *PA
;
2135 int nvar
= P
->Dimension
- exist
- nparam
;
2136 Param_Polyhedron
*PP
= NULL
;
2137 Polyhedron
*C
= Universe_Polyhedron(nparam
);
2141 PP
= Polyhedron2Param_SimplifiedDomain(&PR
,C
,MaxRays
,&CEq
,&CT
);
2145 Param_Domain
*D
, *last
;
2148 for (nd
= 0, D
=PP
->D
; D
; D
=D
->next
, ++nd
)
2151 Polyhedron
**VD
= new (Polyhedron
*)[nd
];
2152 Polyhedron
**fVD
= new (Polyhedron
*)[nd
];
2153 for(nd
= 0, D
=PP
->D
; D
; D
=D
->next
) {
2154 Polyhedron
*rVD
= reduce_domain(D
->Domain
, CT
, CEq
,
2168 /* This doesn't seem to have any effect */
2170 Polyhedron
*CA
= align_context(VD
[0], P
->Dimension
, MaxRays
);
2172 P
= DomainIntersection(P
, CA
, MaxRays
);
2175 Polyhedron_Free(CA
);
2180 if (!EP
&& CT
->NbColumns
!= CT
->NbRows
) {
2181 Polyhedron
*CEqr
= DomainImage(CEq
, CT
, MaxRays
);
2182 Polyhedron
*CA
= align_context(CEqr
, PR
->Dimension
, MaxRays
);
2183 Polyhedron
*I
= DomainIntersection(PR
, CA
, MaxRays
);
2184 Polyhedron_Free(CEqr
);
2185 Polyhedron_Free(CA
);
2187 fprintf(stderr
, "\nER: Eliminate\n");
2188 #endif /* DEBUG_ER */
2189 nparam
-= CT
->NbColumns
- CT
->NbRows
;
2190 EP
= barvinok_enumerate_e(I
, exist
, nparam
, MaxRays
);
2191 addeliminatedparams_enum(EP
, CT
, CEq
, MaxRays
, nparam
);
2195 Polyhedron_Free(PR
);
2198 if (!EP
&& nd
> 1) {
2200 fprintf(stderr
, "\nER: VD\n");
2201 #endif /* DEBUG_ER */
2202 for (int i
= 0; i
< nd
; ++i
) {
2203 Polyhedron
*CA
= align_context(VD
[i
], P
->Dimension
, MaxRays
);
2204 Polyhedron
*I
= DomainIntersection(P
, CA
, MaxRays
);
2207 EP
= barvinok_enumerate_e(I
, exist
, nparam
, MaxRays
);
2209 evalue
*E
= barvinok_enumerate_e(I
, exist
, nparam
, MaxRays
);
2211 free_evalue_refs(E
);
2215 Polyhedron_Free(CA
);
2219 for (int i
= 0; i
< nd
; ++i
) {
2220 Polyhedron_Free(VD
[i
]);
2221 Polyhedron_Free(fVD
[i
]);
2227 if (!EP
&& nvar
== 0) {
2230 Param_Vertices
*V
, *V2
;
2231 Matrix
* M
= Matrix_Alloc(1, P
->Dimension
+2);
2233 FORALL_PVertex_in_ParamPolyhedron(V
, last
, PP
) {
2235 FORALL_PVertex_in_ParamPolyhedron(V2
, last
, PP
) {
2242 for (int i
= 0; i
< exist
; ++i
) {
2243 value_oppose(f
, V
->Vertex
->p
[i
][nparam
+1]);
2244 Vector_Combine(V
->Vertex
->p
[i
],
2246 M
->p
[0] + 1 + nvar
+ exist
,
2247 V2
->Vertex
->p
[i
][nparam
+1],
2251 for (j
= 0; j
< nparam
; ++j
)
2252 if (value_notzero_p(M
->p
[0][1+nvar
+exist
+j
]))
2256 ConstraintSimplify(M
->p
[0], M
->p
[0],
2257 P
->Dimension
+2, &f
);
2258 value_set_si(M
->p
[0][0], 0);
2259 Polyhedron
*para
= AddConstraints(M
->p
[0], 1, P
,
2262 Polyhedron_Free(para
);
2265 Polyhedron
*pos
, *neg
;
2266 value_set_si(M
->p
[0][0], 1);
2267 value_decrement(M
->p
[0][P
->Dimension
+1],
2268 M
->p
[0][P
->Dimension
+1]);
2269 neg
= AddConstraints(M
->p
[0], 1, P
, MaxRays
);
2270 value_set_si(f
, -1);
2271 Vector_Scale(M
->p
[0]+1, M
->p
[0]+1, f
,
2273 value_decrement(M
->p
[0][P
->Dimension
+1],
2274 M
->p
[0][P
->Dimension
+1]);
2275 value_decrement(M
->p
[0][P
->Dimension
+1],
2276 M
->p
[0][P
->Dimension
+1]);
2277 pos
= AddConstraints(M
->p
[0], 1, P
, MaxRays
);
2278 if (emptyQ(neg
) && emptyQ(pos
)) {
2279 Polyhedron_Free(para
);
2280 Polyhedron_Free(pos
);
2281 Polyhedron_Free(neg
);
2285 fprintf(stderr
, "\nER: Order\n");
2286 #endif /* DEBUG_ER */
2287 EP
= barvinok_enumerate_e(para
, exist
, nparam
, MaxRays
);
2290 E
= barvinok_enumerate_e(pos
, exist
, nparam
, MaxRays
);
2292 free_evalue_refs(E
);
2296 E
= barvinok_enumerate_e(neg
, exist
, nparam
, MaxRays
);
2298 free_evalue_refs(E
);
2301 Polyhedron_Free(para
);
2302 Polyhedron_Free(pos
);
2303 Polyhedron_Free(neg
);
2308 } END_FORALL_PVertex_in_ParamPolyhedron
;
2311 } END_FORALL_PVertex_in_ParamPolyhedron
;
2314 /* Search for vertex coordinate to split on */
2315 /* First look for one independent of the parameters */
2316 FORALL_PVertex_in_ParamPolyhedron(V
, last
, PP
) {
2317 for (int i
= 0; i
< exist
; ++i
) {
2319 for (j
= 0; j
< nparam
; ++j
)
2320 if (value_notzero_p(V
->Vertex
->p
[i
][j
]))
2324 value_set_si(M
->p
[0][0], 1);
2325 Vector_Set(M
->p
[0]+1, 0, nvar
+exist
);
2326 Vector_Copy(V
->Vertex
->p
[i
],
2327 M
->p
[0] + 1 + nvar
+ exist
, nparam
+1);
2328 value_oppose(M
->p
[0][1+nvar
+i
],
2329 V
->Vertex
->p
[i
][nparam
+1]);
2331 Polyhedron
*pos
, *neg
;
2332 value_set_si(M
->p
[0][0], 1);
2333 value_decrement(M
->p
[0][P
->Dimension
+1],
2334 M
->p
[0][P
->Dimension
+1]);
2335 neg
= AddConstraints(M
->p
[0], 1, P
, MaxRays
);
2336 value_set_si(f
, -1);
2337 Vector_Scale(M
->p
[0]+1, M
->p
[0]+1, f
,
2339 value_decrement(M
->p
[0][P
->Dimension
+1],
2340 M
->p
[0][P
->Dimension
+1]);
2341 value_decrement(M
->p
[0][P
->Dimension
+1],
2342 M
->p
[0][P
->Dimension
+1]);
2343 pos
= AddConstraints(M
->p
[0], 1, P
, MaxRays
);
2344 if (emptyQ(neg
) || emptyQ(pos
)) {
2345 Polyhedron_Free(pos
);
2346 Polyhedron_Free(neg
);
2349 Polyhedron_Free(pos
);
2350 value_increment(M
->p
[0][P
->Dimension
+1],
2351 M
->p
[0][P
->Dimension
+1]);
2352 pos
= AddConstraints(M
->p
[0], 1, P
, MaxRays
);
2354 fprintf(stderr
, "\nER: Vertex\n");
2355 #endif /* DEBUG_ER */
2356 EP
= enumerate_or(pos
, neg
, exist
, nparam
, MaxRays
);
2361 } END_FORALL_PVertex_in_ParamPolyhedron
;
2365 /* Search for vertex coordinate to split on */
2366 /* Now look for one that depends on the parameters */
2367 FORALL_PVertex_in_ParamPolyhedron(V
, last
, PP
) {
2368 for (int i
= 0; i
< exist
; ++i
) {
2369 value_set_si(M
->p
[0][0], 1);
2370 Vector_Set(M
->p
[0]+1, 0, nvar
+exist
);
2371 Vector_Copy(V
->Vertex
->p
[i
],
2372 M
->p
[0] + 1 + nvar
+ exist
, nparam
+1);
2373 value_oppose(M
->p
[0][1+nvar
+i
],
2374 V
->Vertex
->p
[i
][nparam
+1]);
2376 Polyhedron
*pos
, *neg
;
2377 value_set_si(M
->p
[0][0], 1);
2378 value_decrement(M
->p
[0][P
->Dimension
+1],
2379 M
->p
[0][P
->Dimension
+1]);
2380 neg
= AddConstraints(M
->p
[0], 1, P
, MaxRays
);
2381 value_set_si(f
, -1);
2382 Vector_Scale(M
->p
[0]+1, M
->p
[0]+1, f
,
2384 value_decrement(M
->p
[0][P
->Dimension
+1],
2385 M
->p
[0][P
->Dimension
+1]);
2386 value_decrement(M
->p
[0][P
->Dimension
+1],
2387 M
->p
[0][P
->Dimension
+1]);
2388 pos
= AddConstraints(M
->p
[0], 1, P
, MaxRays
);
2389 if (emptyQ(neg
) || emptyQ(pos
)) {
2390 Polyhedron_Free(pos
);
2391 Polyhedron_Free(neg
);
2394 Polyhedron_Free(pos
);
2395 value_increment(M
->p
[0][P
->Dimension
+1],
2396 M
->p
[0][P
->Dimension
+1]);
2397 pos
= AddConstraints(M
->p
[0], 1, P
, MaxRays
);
2399 fprintf(stderr
, "\nER: ParamVertex\n");
2400 #endif /* DEBUG_ER */
2401 EP
= enumerate_or(pos
, neg
, exist
, nparam
, MaxRays
);
2406 } END_FORALL_PVertex_in_ParamPolyhedron
;
2414 Polyhedron_Free(CEq
);
2418 Param_Polyhedron_Free(PP
);
2424 static evalue
* barvinok_enumerate_e_r(Polyhedron
*P
,
2425 unsigned exist
, unsigned nparam
, unsigned MaxRays
);
2428 static int er_level
= 0;
2430 evalue
* barvinok_enumerate_e(Polyhedron
*P
,
2431 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
2433 fprintf(stderr
, "\nER: level %i\n", er_level
);
2434 int nvar
= P
->Dimension
- exist
- nparam
;
2435 fprintf(stderr
, "%d %d %d\n", nvar
, exist
, nparam
);
2437 Polyhedron_Print(stderr
, P_VALUE_FMT
, P
);
2439 P
= DomainConstraintSimplify(Polyhedron_Copy(P
), MaxRays
);
2440 evalue
*EP
= barvinok_enumerate_e_r(P
, exist
, nparam
, MaxRays
);
2446 evalue
* barvinok_enumerate_e(Polyhedron
*P
,
2447 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
2449 P
= DomainConstraintSimplify(Polyhedron_Copy(P
), MaxRays
);
2450 evalue
*EP
= barvinok_enumerate_e_r(P
, exist
, nparam
, MaxRays
);
2456 static evalue
* barvinok_enumerate_e_r(Polyhedron
*P
,
2457 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
2460 Polyhedron
*U
= Universe_Polyhedron(nparam
);
2461 evalue
*EP
= barvinok_enumerate_ev(P
, U
, MaxRays
);
2462 //char *param_name[] = {"P", "Q", "R", "S", "T" };
2463 //print_evalue(stdout, EP, param_name);
2468 int nvar
= P
->Dimension
- exist
- nparam
;
2469 int len
= P
->Dimension
+ 2;
2472 return new_zero_ep();
2474 if (nvar
== 0 && nparam
== 0) {
2475 evalue
*EP
= new_zero_ep();
2476 barvinok_count(P
, &EP
->x
.n
, MaxRays
);
2477 if (value_pos_p(EP
->x
.n
))
2478 value_set_si(EP
->x
.n
, 1);
2483 for (r
= 0; r
< P
->NbRays
; ++r
)
2484 if (value_zero_p(P
->Ray
[r
][0]) ||
2485 value_zero_p(P
->Ray
[r
][P
->Dimension
+1])) {
2487 for (i
= 0; i
< nvar
; ++i
)
2488 if (value_notzero_p(P
->Ray
[r
][i
+1]))
2492 for (i
= nvar
+ exist
; i
< nvar
+ exist
+ nparam
; ++i
)
2493 if (value_notzero_p(P
->Ray
[r
][i
+1]))
2495 if (i
>= nvar
+ exist
+ nparam
)
2498 if (r
< P
->NbRays
) {
2499 evalue
*EP
= new_zero_ep();
2500 value_set_si(EP
->x
.n
, -1);
2505 for (r
= 0; r
< P
->NbEq
; ++r
)
2506 if ((first
= First_Non_Zero(P
->Constraint
[r
]+1+nvar
, exist
)) != -1)
2509 if (First_Non_Zero(P
->Constraint
[r
]+1+nvar
+first
+1,
2510 exist
-first
-1) != -1) {
2511 Polyhedron
*T
= rotate_along(P
, r
, nvar
, exist
, MaxRays
);
2513 fprintf(stderr
, "\nER: Equality\n");
2514 #endif /* DEBUG_ER */
2515 evalue
*EP
= barvinok_enumerate_e(T
, exist
-1, nparam
, MaxRays
);
2520 fprintf(stderr
, "\nER: Fixed\n");
2521 #endif /* DEBUG_ER */
2523 return barvinok_enumerate_e(P
, exist
-1, nparam
, MaxRays
);
2525 Polyhedron
*T
= Polyhedron_Copy(P
);
2526 SwapColumns(T
, nvar
+1, nvar
+1+first
);
2527 evalue
*EP
= barvinok_enumerate_e(T
, exist
-1, nparam
, MaxRays
);
2534 Vector
*row
= Vector_Alloc(len
);
2535 value_set_si(row
->p
[0], 1);
2540 enum constraint info
[exist
];
2541 for (int i
= 0; i
< exist
; ++i
) {
2543 for (int l
= P
->NbEq
; l
< P
->NbConstraints
; ++l
) {
2544 if (value_negz_p(P
->Constraint
[l
][nvar
+i
+1]))
2546 for (int u
= P
->NbEq
; u
< P
->NbConstraints
; ++u
) {
2547 if (value_posz_p(P
->Constraint
[u
][nvar
+i
+1]))
2549 value_oppose(f
, P
->Constraint
[u
][nvar
+i
+1]);
2550 Vector_Combine(P
->Constraint
[l
]+1, P
->Constraint
[u
]+1, row
->p
+1,
2551 f
, P
->Constraint
[l
][nvar
+i
+1], len
-1);
2552 if (!(info
[i
] & INDEPENDENT
)) {
2554 for (j
= 0; j
< exist
; ++j
)
2555 if (j
!= i
&& value_notzero_p(row
->p
[nvar
+j
+1]))
2558 //printf("independent: i: %d, l: %d, u: %d\n", i, l, u);
2559 info
[i
] = (constraint
)(info
[i
] | INDEPENDENT
);
2562 if (info
[i
] & ALL_POS
) {
2563 value_addto(row
->p
[len
-1], row
->p
[len
-1],
2564 P
->Constraint
[l
][nvar
+i
+1]);
2565 value_addto(row
->p
[len
-1], row
->p
[len
-1], f
);
2566 value_multiply(f
, f
, P
->Constraint
[l
][nvar
+i
+1]);
2567 value_substract(row
->p
[len
-1], row
->p
[len
-1], f
);
2568 value_decrement(row
->p
[len
-1], row
->p
[len
-1]);
2569 Vector_Gcd(row
->p
+1, len
- 2, &f
);
2570 if (value_notone_p(f
)) {
2571 Vector_AntiScale(row
->p
+1, row
->p
+1, f
, len
-2);
2572 mpz_fdiv_q(row
->p
[len
-1], row
->p
[len
-1], f
);
2574 value_set_si(f
, -1);
2575 Vector_Scale(row
->p
+1, row
->p
+1, f
, len
-1);
2576 value_decrement(row
->p
[len
-1], row
->p
[len
-1]);
2577 Polyhedron
*T
= AddConstraints(row
->p
, 1, P
, MaxRays
);
2579 //printf("not all_pos: i: %d, l: %d, u: %d\n", i, l, u);
2580 info
[i
] = (constraint
)(info
[i
] ^ ALL_POS
);
2582 //puts("pos remainder");
2583 //Polyhedron_Print(stdout, P_VALUE_FMT, T);
2586 if (!(info
[i
] & ONE_NEG
)) {
2588 for (j
= 0; j
< exist
; ++j
)
2590 value_notzero_p(P
->Constraint
[l
][nvar
+j
+1]))
2593 for (j
= 0; j
< exist
; ++j
)
2595 value_notzero_p(P
->Constraint
[u
][nvar
+j
+1]))
2598 /* recalculate constant */
2599 /* We actually recalculate the whole row for
2600 * now, because it may have already been scaled
2602 value_oppose(f
, P
->Constraint
[u
][nvar
+i
+1]);
2603 Vector_Combine(P
->Constraint
[l
]+1, P
->Constraint
[u
]+1,
2605 f
, P
->Constraint
[l
][nvar
+i
+1], len
-1);
2607 Vector_Combine(P->Constraint[l]+len-1,
2608 P->Constraint[u]+len-1, row->p+len-1,
2609 f, P->Constraint[l][nvar+i+1], 1);
2611 value_multiply(f
, f
, P
->Constraint
[l
][nvar
+i
+1]);
2612 value_substract(row
->p
[len
-1], row
->p
[len
-1], f
);
2613 value_set_si(f
, -1);
2614 Vector_Scale(row
->p
+1, row
->p
+1, f
, len
-1);
2615 value_decrement(row
->p
[len
-1], row
->p
[len
-1]);
2616 Vector_Gcd(row
->p
+1, len
- 2, &f
);
2617 if (value_notone_p(f
)) {
2618 Vector_AntiScale(row
->p
+1, row
->p
+1, f
, len
-2);
2619 mpz_fdiv_q(row
->p
[len
-1], row
->p
[len
-1], f
);
2621 value_set_si(f
, -1);
2622 Vector_Scale(row
->p
+1, row
->p
+1, f
, len
-1);
2623 value_decrement(row
->p
[len
-1], row
->p
[len
-1]);
2625 //Vector_Print(stdout, P_VALUE_FMT, row);
2626 Polyhedron
*T
= AddConstraints(row
->p
, 1, P
, MaxRays
);
2628 //printf("one_neg i: %d, l: %d, u: %d\n", i, l, u);
2629 info
[i
] = (constraint
)(info
[i
] | ONE_NEG
);
2631 //puts("neg remainder");
2632 //Polyhedron_Print(stdout, P_VALUE_FMT, T);
2636 if (!(info
[i
] & ALL_POS
) && (info
[i
] & ONE_NEG
))
2640 if (info
[i
] & ALL_POS
)
2647 for (int i = 0; i < exist; ++i)
2648 printf("%i: %i\n", i, info[i]);
2650 for (int i
= 0; i
< exist
; ++i
)
2651 if (info
[i
] & ALL_POS
) {
2653 fprintf(stderr
, "\nER: Positive\n");
2654 #endif /* DEBUG_ER */
2656 // Maybe we should chew off some of the fat here
2657 Matrix
*M
= Matrix_Alloc(P
->Dimension
, P
->Dimension
+1);
2658 for (int j
= 0; j
< P
->Dimension
; ++j
)
2659 value_set_si(M
->p
[j
][j
+ (j
>= i
+nvar
)], 1);
2660 Polyhedron
*T
= Polyhedron_Image(P
, M
, MaxRays
);
2662 evalue
*EP
= barvinok_enumerate_e(T
, exist
-1, nparam
, MaxRays
);
2668 for (int i
= 0; i
< exist
; ++i
)
2669 if (info
[i
] & ONE_NEG
) {
2671 fprintf(stderr
, "\nER: Negative\n");
2672 #endif /* DEBUG_ER */
2676 return barvinok_enumerate_e(P
, exist
-1, nparam
, MaxRays
);
2678 Polyhedron
*T
= Polyhedron_Copy(P
);
2679 SwapColumns(T
, nvar
+1, nvar
+1+i
);
2680 evalue
*EP
= barvinok_enumerate_e(T
, exist
-1, nparam
, MaxRays
);
2685 for (int i
= 0; i
< exist
; ++i
)
2686 if (info
[i
] & INDEPENDENT
) {
2687 Polyhedron
*pos
, *neg
;
2689 /* Find constraint again and split off negative part */
2691 if (SplitOnVar(P
, i
, nvar
, len
, exist
, MaxRays
,
2692 row
, f
, true, &pos
, &neg
)) {
2694 fprintf(stderr
, "\nER: Split\n");
2695 #endif /* DEBUG_ER */
2698 barvinok_enumerate_e(neg
, exist
-1, nparam
, MaxRays
);
2700 barvinok_enumerate_e(pos
, exist
, nparam
, MaxRays
);
2702 free_evalue_refs(E
);
2704 Polyhedron_Free(neg
);
2705 Polyhedron_Free(pos
);
2716 EP
= enumerate_sure(P
, exist
, nparam
, MaxRays
);
2720 EP
= enumerate_sure2(P
, exist
, nparam
, MaxRays
);
2724 F
= unfringe(P
, MaxRays
);
2725 if (!PolyhedronIncludes(F
, P
)) {
2727 fprintf(stderr
, "\nER: Fringed\n");
2728 #endif /* DEBUG_ER */
2729 EP
= barvinok_enumerate_e(F
, exist
, nparam
, MaxRays
);
2736 EP
= enumerate_vd(&P
, exist
, nparam
, MaxRays
);
2741 EP
= enumerate_sum(P
, exist
, nparam
, MaxRays
);
2748 Polyhedron
*pos
, *neg
;
2749 for (i
= 0; i
< exist
; ++i
)
2750 if (SplitOnVar(P
, i
, nvar
, len
, exist
, MaxRays
,
2751 row
, f
, false, &pos
, &neg
))
2756 EP
= enumerate_or(pos
, neg
, exist
, nparam
, MaxRays
);