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])) {
707 value_oppose(I
->Constraint
[0][2], I
->Constraint
[0][2]);
708 /* There should never be a remainder here */
709 if (value_pos_p(I
->Constraint
[0][1]))
710 mpz_fdiv_q(min
, I
->Constraint
[0][2], I
->Constraint
[0][1]);
712 mpz_fdiv_q(min
, I
->Constraint
[0][2], I
->Constraint
[0][1]);
713 value_assign(max
, min
);
714 } else for (i
= 0; i
< I
->NbConstraints
; ++i
) {
715 value_oppose(I
->Constraint
[i
][2], I
->Constraint
[i
][2]);
716 if (value_pos_p(I
->Constraint
[i
][1]))
717 mpz_cdiv_q(min
, I
->Constraint
[i
][2], I
->Constraint
[i
][1]);
719 mpz_fdiv_q(max
, I
->Constraint
[i
][2], I
->Constraint
[i
][1]);
724 mpz_cdiv_q(min
, min
, d
);
726 mpz_fdiv_q(min
, min
, d
);
727 mpz_fdiv_q(max
, max
, d
);
728 if (value_eq(min
, max
)) {
731 value_assign(*ret
, min
);
739 * Normalize linear expression coef modulo m
740 * Removes common factor and reduces coefficients
741 * Returns index of first non-zero coefficient or len
743 static int normal_mod(Value
*coef
, int len
, Value
*m
)
748 Vector_Gcd(coef
, len
, &gcd
);
750 Vector_AntiScale(coef
, coef
, gcd
, len
);
752 value_division(*m
, *m
, gcd
);
759 for (j
= 0; j
< len
; ++j
)
760 mpz_fdiv_r(coef
[j
], coef
[j
], *m
);
761 for (j
= 0; j
< len
; ++j
)
762 if (value_notzero_p(coef
[j
]))
769 static void mask(Matrix
*f
, evalue
*factor
)
771 int nr
= f
->NbRows
, nc
= f
->NbColumns
;
774 for (n
= 0; n
< nr
&& value_notzero_p(f
->p
[n
][nc
-1]); ++n
)
775 if (value_notone_p(f
->p
[n
][nc
-1]) &&
776 value_notmone_p(f
->p
[n
][nc
-1]))
787 for (n
= 0; n
< nr
; ++n
) {
788 value_assign(m
, f
->p
[n
][nc
-1]);
789 if (value_one_p(m
) || value_mone_p(m
))
792 int j
= normal_mod(f
->p
[n
], nc
-1, &m
);
794 free_evalue_refs(factor
);
795 value_init(factor
->d
);
796 evalue_set_si(factor
, 0, 1);
800 values2zz(f
->p
[n
], row
, nc
-1);
803 if (j
< (nc
-1)-1 && row
[j
] > g
/2) {
804 for (int k
= j
; k
< (nc
-1); ++k
)
810 value_set_si(EP
.d
, 0);
811 EP
.x
.p
= new_enode(relation
, 2, 0);
812 value_clear(EP
.x
.p
->arr
[1].d
);
813 EP
.x
.p
->arr
[1] = *factor
;
814 evalue
*ev
= &EP
.x
.p
->arr
[0];
815 value_set_si(ev
->d
, 0);
816 ev
->x
.p
= new_enode(modulo
, 3, VALUE_TO_INT(m
));
817 evalue_set_si(&ev
->x
.p
->arr
[1], 0, 1);
818 evalue_set_si(&ev
->x
.p
->arr
[2], 1, 1);
819 evalue
*E
= multi_monom(row
);
820 value_clear(ev
->x
.p
->arr
[0].d
);
821 ev
->x
.p
->arr
[0] = *E
;
832 static void mask(Matrix
*f
, evalue
*factor
)
834 int nr
= f
->NbRows
, nc
= f
->NbColumns
;
837 for (n
= 0; n
< nr
&& value_notzero_p(f
->p
[n
][nc
-1]); ++n
)
838 if (value_notone_p(f
->p
[n
][nc
-1]) &&
839 value_notmone_p(f
->p
[n
][nc
-1]))
847 unsigned np
= nc
- 2;
848 Vector
*lcm
= Vector_Alloc(np
);
849 Vector
*val
= Vector_Alloc(nc
);
850 Vector_Set(val
->p
, 0, nc
);
851 value_set_si(val
->p
[np
], 1);
852 Vector_Set(lcm
->p
, 1, np
);
853 for (n
= 0; n
< nr
; ++n
) {
854 if (value_one_p(f
->p
[n
][nc
-1]) ||
855 value_mone_p(f
->p
[n
][nc
-1]))
857 for (int j
= 0; j
< np
; ++j
)
858 if (value_notzero_p(f
->p
[n
][j
])) {
859 Gcd(f
->p
[n
][j
], f
->p
[n
][nc
-1], &tmp
);
860 value_division(tmp
, f
->p
[n
][nc
-1], tmp
);
861 value_lcm(tmp
, lcm
->p
[j
], &lcm
->p
[j
]);
866 mask_r(f
, nr
, lcm
, 0, val
, &EP
);
871 free_evalue_refs(&EP
);
882 static bool mod_needed(Polyhedron
*PD
, vec_ZZ
& num
, Value d
, evalue
*E
)
884 Value
*q
= fixed_quotient(PD
, num
, d
, false);
889 value_oppose(*q
, *q
);
892 value_set_si(EV
.d
, 1);
894 value_multiply(EV
.x
.n
, *q
, d
);
896 free_evalue_refs(&EV
);
902 static void ceil_mod(Value
*coef
, int len
, Value d
, ZZ
& f
, evalue
*EP
, Polyhedron
*PD
)
908 Vector_Scale(coef
, coef
, m
, len
);
911 int j
= normal_mod(coef
, len
, &m
);
919 values2zz(coef
, num
, len
);
926 evalue_set_si(&tmp
, 0, 1);
928 if (j
< len
-1 && num
[j
] > g
/2) {
929 for (int k
= j
; k
< len
-1; ++k
)
932 num
[len
-1] = g
- 1 - num
[len
-1];
933 value_assign(tmp
.d
, m
);
935 zz2value(t
, tmp
.x
.n
);
941 ZZ t
= num
[len
-1] * f
;
942 zz2value(t
, tmp
.x
.n
);
943 value_assign(tmp
.d
, m
);
946 evalue
*E
= multi_monom(num
);
950 if (PD
&& !mod_needed(PD
, num
, m
, E
)) {
953 value_assign(EV
.d
, m
);
957 value_set_si(EV
.d
, 0);
958 EV
.x
.p
= new_enode(modulo
, 3, VALUE_TO_INT(m
));
959 evalue_copy(&EV
.x
.p
->arr
[0], E
);
960 evalue_set_si(&EV
.x
.p
->arr
[1], 0, 1);
961 value_init(EV
.x
.p
->arr
[2].x
.n
);
962 zz2value(f
, EV
.x
.p
->arr
[2].x
.n
);
963 value_assign(EV
.x
.p
->arr
[2].d
, m
);
968 free_evalue_refs(&EV
);
973 free_evalue_refs(&tmp
);
979 evalue
* bv_ceil3(Value
*coef
, int len
, Value d
, Polyhedron
*P
)
981 Vector
*val
= Vector_Alloc(len
);
986 Vector_Scale(coef
, val
->p
, t
, len
);
987 value_absolute(t
, d
);
990 values2zz(val
->p
, num
, len
);
991 evalue
*EP
= multi_monom(num
);
996 value_set_si(tmp
.x
.n
, 1);
997 value_assign(tmp
.d
, t
);
1003 ceil_mod(val
->p
, len
, t
, one
, EP
, P
);
1006 /* copy EP to malloc'ed evalue */
1012 free_evalue_refs(&tmp
);
1019 evalue
* lattice_point(
1020 Polyhedron
*i
, vec_ZZ
& lambda
, Matrix
*W
, Value lcm
, Polyhedron
*PD
)
1022 unsigned nparam
= W
->NbColumns
- 1;
1024 Matrix
* Rays
= rays2(i
);
1025 Matrix
*T
= Transpose(Rays
);
1026 Matrix
*T2
= Matrix_Copy(T
);
1027 Matrix
*inv
= Matrix_Alloc(T2
->NbRows
, T2
->NbColumns
);
1028 int ok
= Matrix_Inverse(T2
, inv
);
1033 matrix2zz(W
, vertex
, W
->NbRows
, W
->NbColumns
);
1036 num
= lambda
* vertex
;
1038 evalue
*EP
= multi_monom(num
);
1042 value_init(tmp
.x
.n
);
1043 value_set_si(tmp
.x
.n
, 1);
1044 value_assign(tmp
.d
, lcm
);
1048 Matrix
*L
= Matrix_Alloc(inv
->NbRows
, W
->NbColumns
);
1049 Matrix_Product(inv
, W
, L
);
1052 matrix2zz(T
, RT
, T
->NbRows
, T
->NbColumns
);
1055 vec_ZZ p
= lambda
* RT
;
1057 for (int i
= 0; i
< L
->NbRows
; ++i
) {
1058 ceil_mod(L
->p
[i
], nparam
+1, lcm
, p
[i
], EP
, PD
);
1064 free_evalue_refs(&tmp
);
1068 evalue
* lattice_point(
1069 Polyhedron
*i
, vec_ZZ
& lambda
, Matrix
*W
, Value lcm
, Polyhedron
*PD
)
1071 Matrix
*T
= Transpose(W
);
1072 unsigned nparam
= T
->NbRows
- 1;
1074 evalue
*EP
= new evalue();
1076 evalue_set_si(EP
, 0, 1);
1079 Vector
*val
= Vector_Alloc(nparam
+1);
1080 value_set_si(val
->p
[nparam
], 1);
1081 ZZ
offset(INIT_VAL
, 0);
1083 vertex_period(i
, lambda
, T
, lcm
, 0, val
, EP
, &ev
, offset
);
1086 free_evalue_refs(&ev
);
1097 Param_Vertices
* V
, Polyhedron
*i
, vec_ZZ
& lambda
, term_info
* term
,
1100 unsigned nparam
= V
->Vertex
->NbColumns
- 2;
1101 unsigned dim
= i
->Dimension
;
1103 vertex
.SetDims(V
->Vertex
->NbRows
, nparam
+1);
1107 value_set_si(lcm
, 1);
1108 for (int j
= 0; j
< V
->Vertex
->NbRows
; ++j
) {
1109 value_lcm(lcm
, V
->Vertex
->p
[j
][nparam
+1], &lcm
);
1111 if (value_notone_p(lcm
)) {
1112 Matrix
* mv
= Matrix_Alloc(dim
, nparam
+1);
1113 for (int j
= 0 ; j
< dim
; ++j
) {
1114 value_division(tmp
, lcm
, V
->Vertex
->p
[j
][nparam
+1]);
1115 Vector_Scale(V
->Vertex
->p
[j
], mv
->p
[j
], tmp
, nparam
+1);
1118 term
->E
= lattice_point(i
, lambda
, mv
, lcm
, PD
);
1126 for (int i
= 0; i
< V
->Vertex
->NbRows
; ++i
) {
1127 assert(value_one_p(V
->Vertex
->p
[i
][nparam
+1])); // for now
1128 values2zz(V
->Vertex
->p
[i
], vertex
[i
], nparam
+1);
1132 num
= lambda
* vertex
;
1136 for (int j
= 0; j
< nparam
; ++j
)
1142 term
->E
= multi_monom(num
);
1146 term
->constant
= num
[nparam
];
1149 term
->coeff
= num
[p
];
1156 void normalize(Polyhedron
*i
, vec_ZZ
& lambda
, ZZ
& sign
, ZZ
& num
, vec_ZZ
& den
)
1158 unsigned dim
= i
->Dimension
;
1162 rays
.SetDims(dim
, dim
);
1163 add_rays(rays
, i
, &r
);
1164 den
= rays
* lambda
;
1167 for (int j
= 0; j
< den
.length(); ++j
) {
1171 den
[j
] = abs(den
[j
]);
1179 void barvinok_count(Polyhedron
*P
, Value
* result
, unsigned NbMaxCons
)
1181 Polyhedron
** vcone
;
1184 sign
.SetLength(ncone
);
1192 value_set_si(*result
, 0);
1196 for (; r
< P
->NbRays
; ++r
)
1197 if (value_zero_p(P
->Ray
[r
][P
->Dimension
+1]))
1199 if (P
->NbBid
!=0 || r
< P
->NbRays
) {
1200 value_set_si(*result
, -1);
1204 P
= remove_equalities(P
);
1207 value_set_si(*result
, 0);
1213 value_set_si(factor
, 1);
1214 Q
= Polyhedron_Reduce(P
, &factor
);
1221 if (P
->Dimension
== 0) {
1222 value_assign(*result
, factor
);
1225 value_clear(factor
);
1230 vcone
= new (Polyhedron
*)[P
->NbRays
];
1232 for (int j
= 0; j
< P
->NbRays
; ++j
) {
1234 Polyhedron
*C
= supporting_cone(P
, j
);
1235 decompose(C
, &vcone
[j
], &npos
, &nneg
, NbMaxCons
);
1236 ncone
+= npos
+ nneg
;
1237 sign
.SetLength(ncone
);
1238 for (int k
= 0; k
< npos
; ++k
)
1239 sign
[ncone
-nneg
-k
-1] = 1;
1240 for (int k
= 0; k
< nneg
; ++k
)
1241 sign
[ncone
-k
-1] = -1;
1245 rays
.SetDims(ncone
* dim
, dim
);
1247 for (int j
= 0; j
< P
->NbRays
; ++j
) {
1248 for (Polyhedron
*i
= vcone
[j
]; i
; i
= i
->next
) {
1249 assert(i
->NbRays
-1 == dim
);
1250 add_rays(rays
, i
, &r
);
1254 nonorthog(rays
, lambda
);
1258 num
.SetLength(ncone
);
1259 den
.SetDims(ncone
,dim
);
1262 for (int j
= 0; j
< P
->NbRays
; ++j
) {
1263 for (Polyhedron
*i
= vcone
[j
]; i
; i
= i
->next
) {
1264 lattice_point(P
->Ray
[j
]+1, i
, lambda
, num
[f
]);
1265 normalize(i
, lambda
, sign
[f
], num
[f
], den
[f
]);
1270 for (int j
= 1; j
< num
.length(); ++j
)
1273 for (int j
= 0; j
< num
.length(); ++j
)
1279 for (int j
= 0; j
< P
->NbRays
; ++j
) {
1280 for (Polyhedron
*i
= vcone
[j
]; i
; i
= i
->next
) {
1281 dpoly
d(dim
, num
[f
]);
1282 dpoly
n(dim
, den
[f
][0], 1);
1283 for (int k
= 1; k
< dim
; ++k
) {
1284 dpoly
fact(dim
, den
[f
][k
], 1);
1287 d
.div(n
, count
, sign
[f
]);
1291 assert(value_one_p(&count
[0]._mp_den
));
1292 value_multiply(*result
, &count
[0]._mp_num
, factor
);
1295 for (int j
= 0; j
< P
->NbRays
; ++j
)
1296 Domain_Free(vcone
[j
]);
1302 value_clear(factor
);
1305 static void uni_polynom(int param
, Vector
*c
, evalue
*EP
)
1307 unsigned dim
= c
->Size
-2;
1309 value_set_si(EP
->d
,0);
1310 EP
->x
.p
= new_enode(polynomial
, dim
+1, param
+1);
1311 for (int j
= 0; j
<= dim
; ++j
)
1312 evalue_set(&EP
->x
.p
->arr
[j
], c
->p
[j
], c
->p
[dim
+1]);
1315 static void multi_polynom(Vector
*c
, evalue
* X
, evalue
*EP
)
1317 unsigned dim
= c
->Size
-2;
1321 evalue_set(&EC
, c
->p
[dim
], c
->p
[dim
+1]);
1324 evalue_set(EP
, c
->p
[dim
], c
->p
[dim
+1]);
1326 for (int i
= dim
-1; i
>= 0; --i
) {
1328 value_assign(EC
.x
.n
, c
->p
[i
]);
1331 free_evalue_refs(&EC
);
1334 Polyhedron
*unfringe (Polyhedron
*P
, unsigned MaxRays
)
1336 int len
= P
->Dimension
+2;
1337 Polyhedron
*T
, *R
= P
;
1341 Polyhedron_Print(stdout
, P_VALUE_FMT
, P
);
1342 Vector
*row
= Vector_Alloc(len
);
1343 value_set_si(row
->p
[0], 1);
1345 R
= DomainConstraintSimplify(Polyhedron_Copy(P
), MaxRays
);
1347 Matrix
*M
= Matrix_Alloc(2, len
-1);
1348 value_set_si(M
->p
[1][len
-2], 1);
1349 for (int v
= 0; v
< P
->Dimension
; ++v
) {
1350 value_set_si(M
->p
[0][v
], 1);
1351 Polyhedron
*I
= Polyhedron_Image(P
, M
, 2+1);
1352 Polyhedron_Print(stdout
, P_VALUE_FMT
, I
);
1353 value_set_si(M
->p
[0][v
], 0);
1354 for (int r
= 0; r
< I
->NbConstraints
; ++r
) {
1355 if (value_zero_p(I
->Constraint
[r
][0]))
1357 if (value_zero_p(I
->Constraint
[r
][1]))
1359 if (value_one_p(I
->Constraint
[r
][1]))
1361 if (value_mone_p(I
->Constraint
[r
][1]))
1363 value_absolute(g
, I
->Constraint
[r
][1]);
1364 Vector_Set(row
->p
+1, 0, len
-2);
1365 value_division(row
->p
[1+v
], I
->Constraint
[r
][1], g
);
1366 mpz_fdiv_q(row
->p
[len
-1], I
->Constraint
[r
][2], g
);
1368 Vector_Print(stdout
, P_VALUE_FMT
, row
);
1370 R
= AddConstraints(row
->p
, 1, R
, MaxRays
);
1376 Polyhedron_Print(stdout
, P_VALUE_FMT
, R
);
1381 evalue
* barvinok_enumerate_ev(Polyhedron
*P
, Polyhedron
* C
, unsigned MaxRays
)
1383 //P = unfringe(P, MaxRays);
1384 Polyhedron
*CEq
= NULL
, *rVD
, *pVD
, *CA
;
1386 Param_Polyhedron
*PP
= NULL
;
1387 Param_Domain
*D
, *next
;
1390 unsigned nparam
= C
->Dimension
;
1391 evalue
*eres
= new evalue
;
1392 value_init(eres
->d
);
1393 value_set_si(eres
->d
, 0);
1396 value_init(factor
.d
);
1397 evalue_set_si(&factor
, 1, 1);
1399 CA
= align_context(C
, P
->Dimension
, MaxRays
);
1400 P
= DomainIntersection(P
, CA
, MaxRays
);
1401 Polyhedron_Free(CA
);
1403 if (C
->Dimension
== 0 || emptyQ(P
)) {
1405 eres
->x
.p
= new_enode(partition
, 2, -1);
1406 EVALUE_SET_DOMAIN(eres
->x
.p
->arr
[0],
1407 DomainConstraintSimplify(CEq
? CEq
: Polyhedron_Copy(C
), MaxRays
));
1408 value_set_si(eres
->x
.p
->arr
[1].d
, 1);
1409 value_init(eres
->x
.p
->arr
[1].x
.n
);
1411 value_set_si(eres
->x
.p
->arr
[1].x
.n
, 0);
1413 barvinok_count(P
, &eres
->x
.p
->arr
[1].x
.n
, MaxRays
);
1415 emul(&factor
, eres
);
1416 reduce_evalue(eres
);
1417 free_evalue_refs(&factor
);
1422 Param_Polyhedron_Free(PP
);
1426 for (r
= 0; r
< P
->NbRays
; ++r
)
1427 if (value_zero_p(P
->Ray
[r
][0]) ||
1428 value_zero_p(P
->Ray
[r
][P
->Dimension
+1])) {
1430 for (i
= P
->Dimension
- nparam
; i
< P
->Dimension
; ++i
)
1431 if (value_notzero_p(P
->Ray
[r
][i
+1]))
1433 if (i
>= P
->Dimension
)
1441 P
= remove_equalities_p(P
, P
->Dimension
-nparam
, &f
);
1444 if (P
->Dimension
== nparam
) {
1446 P
= Universe_Polyhedron(0);
1451 Polyhedron
*Q
= ParamPolyhedron_Reduce(P
, P
->Dimension
-nparam
, &factor
);
1454 if (Q
->Dimension
== nparam
) {
1456 P
= Universe_Polyhedron(0);
1461 Polyhedron
*oldP
= P
;
1462 PP
= Polyhedron2Param_SimplifiedDomain(&P
,C
,MaxRays
,&CEq
,&CT
);
1464 Polyhedron_Free(oldP
);
1466 if (isIdentity(CT
)) {
1470 assert(CT
->NbRows
!= CT
->NbColumns
);
1471 if (CT
->NbRows
== 1) // no more parameters
1473 nparam
= CT
->NbRows
- 1;
1476 unsigned dim
= P
->Dimension
- nparam
;
1477 Polyhedron
** vcone
= new (Polyhedron
*)[PP
->nbV
];
1478 int * npos
= new int[PP
->nbV
];
1479 int * nneg
= new int[PP
->nbV
];
1483 for (i
= 0, V
= PP
->V
; V
; ++i
, V
= V
->next
) {
1484 Polyhedron
*C
= supporting_cone_p(P
, V
);
1485 decompose(C
, &vcone
[i
], &npos
[i
], &nneg
[i
], MaxRays
);
1488 Vector
*c
= Vector_Alloc(dim
+2);
1491 for (nd
= 0, D
=PP
->D
; D
; ++nd
, D
=D
->next
);
1492 struct section
{ Polyhedron
* D
; evalue E
; };
1493 section
*s
= new section
[nd
];
1495 for(nd
= 0, D
=PP
->D
; D
; D
=next
) {
1497 D
->Domain
= DomainConstraintSimplify(D
->Domain
, MaxRays
);
1499 pVD
= rVD
= D
->Domain
;
1503 Dt
= CT
? DomainPreimage(D
->Domain
,CT
,MaxRays
) : D
->Domain
;
1504 rVD
= DomainIntersection(Dt
,CEq
,MaxRays
);
1505 rVD
= DomainConstraintSimplify(rVD
, MaxRays
);
1507 /* if rVD is empty or too small in geometric dimension */
1508 if(!rVD
|| emptyQ(rVD
) ||
1509 (rVD
->Dimension
-rVD
->NbEq
< Dt
->Dimension
-Dt
->NbEq
-CEq
->NbEq
)) {
1514 continue; /* empty validity domain */
1518 pVD
= CT
? DomainImage(rVD
,CT
,MaxRays
) : rVD
;
1519 for (Param_Domain
*D2
= D
->next
; D2
; D2
=D2
->next
) {
1520 Polyhedron
*T
= D2
->Domain
;
1521 D2
->Domain
= DomainDifference(D2
->Domain
, D
->Domain
, MaxRays
);
1522 //D2->Domain = DomainConstraintSimplify(D2->Domain, MaxRays);
1527 sign
.SetLength(ncone
);
1528 FORALL_PVertex_in_ParamPolyhedron(V
,D
,PP
) // _i is internal counter
1529 ncone
+= npos
[_i
] + nneg
[_i
];
1530 sign
.SetLength(ncone
);
1531 for (int k
= 0; k
< npos
[_i
]; ++k
)
1532 sign
[ncone
-nneg
[_i
]-k
-1] = 1;
1533 for (int k
= 0; k
< nneg
[_i
]; ++k
)
1534 sign
[ncone
-k
-1] = -1;
1535 END_FORALL_PVertex_in_ParamPolyhedron
;
1538 rays
.SetDims(ncone
* dim
, dim
);
1540 FORALL_PVertex_in_ParamPolyhedron(V
,D
,PP
) // _i is internal counter
1541 for (Polyhedron
*i
= vcone
[_i
]; i
; i
= i
->next
) {
1542 assert(i
->NbRays
-1 == dim
);
1543 add_rays(rays
, i
, &r
);
1545 END_FORALL_PVertex_in_ParamPolyhedron
;
1547 nonorthog(rays
, lambda
);
1550 den
.SetDims(ncone
,dim
);
1551 term_info
*num
= new term_info
[ncone
];
1554 FORALL_PVertex_in_ParamPolyhedron(V
,D
,PP
)
1555 for (Polyhedron
*i
= vcone
[_i
]; i
; i
= i
->next
) {
1556 lattice_point(V
, i
, lambda
, &num
[f
], pVD
);
1557 normalize(i
, lambda
, sign
[f
], num
[f
].constant
, den
[f
]);
1560 END_FORALL_PVertex_in_ParamPolyhedron
;
1561 ZZ min
= num
[0].constant
;
1562 for (int j
= 1; j
< ncone
; ++j
)
1563 if (num
[j
].constant
< min
)
1564 min
= num
[j
].constant
;
1565 for (int j
= 0; j
< ncone
; ++j
)
1566 num
[j
].constant
-= min
;
1568 value_init(s
[nd
].E
.d
);
1569 evalue_set_si(&s
[nd
].E
, 0, 1);
1572 FORALL_PVertex_in_ParamPolyhedron(V
,D
,PP
)
1573 for (Polyhedron
*i
= vcone
[_i
]; i
; i
= i
->next
) {
1574 dpoly
n(dim
, den
[f
][0], 1);
1575 for (int k
= 1; k
< dim
; ++k
) {
1576 dpoly
fact(dim
, den
[f
][k
], 1);
1579 if (num
[f
].E
!= NULL
) {
1580 ZZ
one(INIT_VAL
, 1);
1581 dpoly_n
d(dim
, num
[f
].constant
, one
);
1582 d
.div(n
, c
, sign
[f
]);
1584 multi_polynom(c
, num
[f
].E
, &EV
);
1585 eadd(&EV
, &s
[nd
].E
);
1586 free_evalue_refs(&EV
);
1587 free_evalue_refs(num
[f
].E
);
1589 } else if (num
[f
].pos
!= -1) {
1590 dpoly_n
d(dim
, num
[f
].constant
, num
[f
].coeff
);
1591 d
.div(n
, c
, sign
[f
]);
1593 uni_polynom(num
[f
].pos
, c
, &EV
);
1594 eadd(&EV
, &s
[nd
].E
);
1595 free_evalue_refs(&EV
);
1597 mpq_set_si(count
, 0, 1);
1598 dpoly
d(dim
, num
[f
].constant
);
1599 d
.div(n
, count
, sign
[f
]);
1602 evalue_set(&EV
, &count
[0]._mp_num
, &count
[0]._mp_den
);
1603 eadd(&EV
, &s
[nd
].E
);
1604 free_evalue_refs(&EV
);
1608 END_FORALL_PVertex_in_ParamPolyhedron
;
1614 addeliminatedparams_evalue(&s
[nd
].E
, CT
);
1621 eres
->x
.p
= new_enode(partition
, 2*nd
, -1);
1622 for (int j
= 0; j
< nd
; ++j
) {
1623 EVALUE_SET_DOMAIN(eres
->x
.p
->arr
[2*j
], s
[j
].D
);
1624 value_clear(eres
->x
.p
->arr
[2*j
+1].d
);
1625 eres
->x
.p
->arr
[2*j
+1] = s
[j
].E
;
1631 for (int j
= 0; j
< PP
->nbV
; ++j
)
1632 Domain_Free(vcone
[j
]);
1638 Polyhedron_Free(CEq
);
1643 Enumeration
* barvinok_enumerate(Polyhedron
*P
, Polyhedron
* C
, unsigned MaxRays
)
1645 Enumeration
*en
, *res
= NULL
;
1646 evalue
*EP
= barvinok_enumerate_ev(P
, C
, MaxRays
);
1647 for (int i
= 0; i
< EP
->x
.p
->size
/2; ++i
) {
1648 en
= (Enumeration
*)malloc(sizeof(Enumeration
));
1651 res
->ValidityDomain
= EVALUE_DOMAIN(EP
->x
.p
->arr
[2*i
]);
1652 value_clear(EP
->x
.p
->arr
[2*i
].d
);
1653 res
->EP
= EP
->x
.p
->arr
[2*i
+1];
1661 static void SwapColumns(Value
**V
, int n
, int i
, int j
)
1663 for (int r
= 0; r
< n
; ++r
)
1664 value_swap(V
[r
][i
], V
[r
][j
]);
1667 static void SwapColumns(Polyhedron
*P
, int i
, int j
)
1669 SwapColumns(P
->Constraint
, P
->NbConstraints
, i
, j
);
1670 SwapColumns(P
->Ray
, P
->NbRays
, i
, j
);
1676 INDEPENDENT
= 1 << 2,
1679 evalue
* barvinok_enumerate_e(Polyhedron
*P
,
1680 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
1682 //Polyhedron_Print(stderr, P_VALUE_FMT, P);
1684 Polyhedron
*U
= Universe_Polyhedron(nparam
);
1685 evalue
*EP
= barvinok_enumerate_ev(P
, U
, MaxRays
);
1686 //char *param_name[] = {"P", "Q", "R", "S", "T" };
1687 //print_evalue(stdout, EP, param_name);
1692 int nvar
= P
->Dimension
- exist
- nparam
;
1693 int len
= P
->Dimension
+ 2;
1695 //printf("%d %d %d\n", nvar, exist, nparam);
1699 for (r
= 0; r
< P
->NbEq
; ++r
)
1700 if ((first
= First_Non_Zero(P
->Constraint
[r
]+1+nvar
, exist
)) != -1)
1703 if (First_Non_Zero(P
->Constraint
[r
]+1+nvar
+first
+1,
1704 exist
-first
-1) != -1) {
1708 Vector
*row
= Vector_Alloc(exist
);
1709 Vector_Copy(P
->Constraint
[r
]+1+nvar
, row
->p
, exist
);
1710 Vector_Gcd(row
->p
, exist
, &g
);
1711 if (value_notone_p(g
))
1712 Vector_AntiScale(row
->p
, row
->p
, g
, exist
);
1715 Matrix
*M
= unimodular_complete(row
);
1716 Matrix
*M2
= Matrix_Alloc(P
->Dimension
+1, P
->Dimension
+1);
1717 for (r
= 0; r
< nvar
; ++r
)
1718 value_set_si(M2
->p
[r
][r
], 1);
1719 for ( ; r
< nvar
+exist
; ++r
)
1720 Vector_Copy(M
->p
[r
-nvar
], M2
->p
[r
]+nvar
, exist
);
1721 for ( ; r
< P
->Dimension
+1; ++r
)
1722 value_set_si(M2
->p
[r
][r
], 1);
1723 Polyhedron
*T
= Polyhedron_Image(P
, M2
, MaxRays
);
1724 evalue
*EP
= barvinok_enumerate_e(T
, exist
-1, nparam
, MaxRays
);
1732 return barvinok_enumerate_e(P
, exist
-1, nparam
, MaxRays
);
1734 Polyhedron
*T
= Polyhedron_Copy(P
);
1735 SwapColumns(T
, nvar
+1, nvar
+1+first
);
1736 evalue
*EP
= barvinok_enumerate_e(T
, exist
-1, nparam
, MaxRays
);
1743 Vector
*row
= Vector_Alloc(len
);
1744 value_set_si(row
->p
[0], 1);
1749 enum constraint info
[exist
];
1750 for (int i
= 0; i
< exist
; ++i
) {
1752 for (int l
= P
->NbEq
; l
< P
->NbConstraints
; ++l
) {
1753 if (value_negz_p(P
->Constraint
[l
][nvar
+i
+1]))
1755 for (int u
= P
->NbEq
; u
< P
->NbConstraints
; ++u
) {
1756 if (value_posz_p(P
->Constraint
[u
][nvar
+i
+1]))
1758 value_oppose(f
, P
->Constraint
[u
][nvar
+i
+1]);
1759 Vector_Combine(P
->Constraint
[l
]+1, P
->Constraint
[u
]+1, row
->p
+1,
1760 f
, P
->Constraint
[l
][nvar
+i
+1], len
-1);
1761 if (!(info
[i
] & INDEPENDENT
)) {
1763 for (j
= 0; j
< exist
; ++j
)
1764 if (j
!= i
&& value_notzero_p(row
->p
[nvar
+j
+1]))
1767 //printf("independent: i: %d, l: %d, u: %d\n", i, l, u);
1768 info
[i
] = (constraint
)(info
[i
] | INDEPENDENT
);
1771 if (info
[i
] & ALL_POS
) {
1772 value_addto(row
->p
[len
-1], row
->p
[len
-1],
1773 P
->Constraint
[l
][nvar
+i
+1]);
1774 value_addto(row
->p
[len
-1], row
->p
[len
-1], f
);
1775 value_multiply(f
, f
, P
->Constraint
[l
][nvar
+i
+1]);
1776 value_substract(row
->p
[len
-1], row
->p
[len
-1], f
);
1777 value_decrement(row
->p
[len
-1], row
->p
[len
-1]);
1778 Vector_Gcd(row
->p
+1, len
- 2, &f
);
1779 if (value_notone_p(f
)) {
1780 Vector_AntiScale(row
->p
+1, row
->p
+1, f
, len
-2);
1781 mpz_fdiv_q(row
->p
[len
-1], row
->p
[len
-1], f
);
1783 value_set_si(f
, -1);
1784 Vector_Scale(row
->p
+1, row
->p
+1, f
, len
-1);
1785 value_decrement(row
->p
[len
-1], row
->p
[len
-1]);
1786 Polyhedron
*T
= AddConstraints(row
->p
, 1, P
, MaxRays
);
1788 //printf("not all_pos: i: %d, l: %d, u: %d\n", i, l, u);
1789 info
[i
] = (constraint
)(info
[i
] ^ ALL_POS
);
1791 //puts("pos remainder");
1792 //Polyhedron_Print(stdout, P_VALUE_FMT, T);
1795 if (!(info
[i
] & ONE_NEG
)) {
1797 for (j
= 0; j
< exist
; ++j
)
1799 (value_notzero_p(P
->Constraint
[l
][nvar
+j
+1]) ||
1800 value_notzero_p(P
->Constraint
[u
][nvar
+j
+1])))
1803 /* recalculate constant */
1804 /* We actually recalculate the whole row for
1805 * now, because it may have already been scaled
1807 value_oppose(f
, P
->Constraint
[u
][nvar
+i
+1]);
1808 Vector_Combine(P
->Constraint
[l
]+1, P
->Constraint
[u
]+1,
1810 f
, P
->Constraint
[l
][nvar
+i
+1], len
-1);
1812 Vector_Combine(P->Constraint[l]+len-1,
1813 P->Constraint[u]+len-1, row->p+len-1,
1814 f, P->Constraint[l][nvar+i+1], 1);
1816 value_multiply(f
, f
, P
->Constraint
[l
][nvar
+i
+1]);
1817 value_substract(row
->p
[len
-1], row
->p
[len
-1], f
);
1818 value_set_si(f
, -1);
1819 Vector_Scale(row
->p
+1, row
->p
+1, f
, len
-1);
1820 value_decrement(row
->p
[len
-1], row
->p
[len
-1]);
1821 Vector_Gcd(row
->p
+1, len
- 2, &f
);
1822 if (value_notone_p(f
)) {
1823 Vector_AntiScale(row
->p
+1, row
->p
+1, f
, len
-2);
1824 mpz_fdiv_q(row
->p
[len
-1], row
->p
[len
-1], f
);
1826 value_set_si(f
, -1);
1827 Vector_Scale(row
->p
+1, row
->p
+1, f
, len
-1);
1828 value_decrement(row
->p
[len
-1], row
->p
[len
-1]);
1830 //Vector_Print(stdout, P_VALUE_FMT, row);
1831 Polyhedron
*T
= AddConstraints(row
->p
, 1, P
, MaxRays
);
1833 //printf("one_neg i: %d, l: %d, u: %d\n", i, l, u);
1834 info
[i
] = (constraint
)(info
[i
] | ONE_NEG
);
1836 //puts("neg remainder");
1837 //Polyhedron_Print(stdout, P_VALUE_FMT, T);
1841 if (!(info
[i
] & ALL_POS
) && (info
[i
] & ONE_NEG
))
1845 if (info
[i
] & ALL_POS
)
1852 for (int i = 0; i < exist; ++i)
1853 printf("%i: %i\n", i, info[i]);
1855 for (int i
= 0; i
< exist
; ++i
)
1856 if (info
[i
] & ALL_POS
) {
1858 // Maybe we should chew off some of the fat here
1859 Matrix
*M
= Matrix_Alloc(P
->Dimension
, P
->Dimension
+1);
1860 for (int j
= 0; j
< P
->Dimension
; ++j
)
1861 value_set_si(M
->p
[j
][j
+ (j
>= i
+nvar
)], 1);
1862 Polyhedron
*T
= Polyhedron_Image(P
, M
, MaxRays
);
1864 evalue
*EP
= barvinok_enumerate_e(T
, exist
-1, nparam
, MaxRays
);
1868 for (int i
= 0; i
< exist
; ++i
)
1869 if (info
[i
] & ONE_NEG
) {
1872 return barvinok_enumerate_e(P
, exist
-1, nparam
, MaxRays
);
1874 Polyhedron
*T
= Polyhedron_Copy(P
);
1875 SwapColumns(T
, nvar
+1, nvar
+1+i
);
1876 evalue
*EP
= barvinok_enumerate_e(T
, exist
-1, nparam
, MaxRays
);
1881 for (int i
= 0; i
< exist
; ++i
)
1882 if (info
[i
] & INDEPENDENT
) {
1883 /* Find constraint again and split off negative part */
1885 for (int l
= P
->NbEq
; l
< P
->NbConstraints
; ++l
) {
1886 if (value_negz_p(P
->Constraint
[l
][nvar
+i
+1]))
1888 for (int u
= P
->NbEq
; u
< P
->NbConstraints
; ++u
) {
1889 if (value_posz_p(P
->Constraint
[u
][nvar
+i
+1]))
1891 value_oppose(f
, P
->Constraint
[u
][nvar
+i
+1]);
1892 Vector_Combine(P
->Constraint
[l
]+1, P
->Constraint
[u
]+1,
1894 f
, P
->Constraint
[l
][nvar
+i
+1], len
-1);
1897 for (j
= 0; j
< exist
; ++j
)
1898 if (j
!= i
&& value_notzero_p(row
->p
[nvar
+j
+1]))
1903 //printf("l: %d, u: %d\n", l, u);
1904 value_multiply(f
, f
, P
->Constraint
[l
][nvar
+i
+1]);
1905 value_substract(row
->p
[len
-1], row
->p
[len
-1], f
);
1906 value_set_si(f
, -1);
1907 Vector_Scale(row
->p
+1, row
->p
+1, f
, len
-1);
1908 value_decrement(row
->p
[len
-1], row
->p
[len
-1]);
1909 Vector_Gcd(row
->p
+1, len
- 2, &f
);
1910 if (value_notone_p(f
)) {
1911 Vector_AntiScale(row
->p
+1, row
->p
+1, f
, len
-2);
1912 mpz_fdiv_q(row
->p
[len
-1], row
->p
[len
-1], f
);
1914 Polyhedron
*neg
= AddConstraints(row
->p
, 1, P
, MaxRays
);
1915 value_set_si(f
, -1);
1916 Vector_Scale(row
->p
+1, row
->p
+1, f
, len
-1);
1917 value_decrement(row
->p
[len
-1], row
->p
[len
-1]);
1918 Polyhedron
*pos
= AddConstraints(row
->p
, 1, P
, MaxRays
);
1920 assert(i
== 0); // for now
1922 barvinok_enumerate_e(neg
, exist
-1, nparam
, MaxRays
);
1924 barvinok_enumerate_e(pos
, exist
, nparam
, MaxRays
);
1926 free_evalue_refs(E
);
1928 Polyhedron_Free(neg
);
1929 Polyhedron_Free(pos
);
1934 assert(0); // can't happen