8 #include <NTL/mat_ZZ.h>
12 #include <polylib/polylibgmp.h>
13 #include "ev_operations.h"
27 using std::ostringstream
;
29 #define ALLOC(p) (((long *) (p))[0])
30 #define SIZE(p) (((long *) (p))[1])
31 #define DATA(p) ((mp_limb_t *) (((long *) (p)) + 2))
33 static void value2zz(Value v
, ZZ
& z
)
35 int sa
= v
[0]._mp_size
;
36 int abs_sa
= sa
< 0 ? -sa
: sa
;
38 _ntl_gsetlength(&z
.rep
, abs_sa
);
39 mp_limb_t
* adata
= DATA(z
.rep
);
40 for (int i
= 0; i
< abs_sa
; ++i
)
41 adata
[i
] = v
[0]._mp_d
[i
];
45 static void zz2value(ZZ
& z
, Value
& v
)
53 int abs_sa
= sa
< 0 ? -sa
: sa
;
55 mp_limb_t
* adata
= DATA(z
.rep
);
56 _mpz_realloc(v
, abs_sa
);
57 for (int i
= 0; i
< abs_sa
; ++i
)
58 v
[0]._mp_d
[i
] = adata
[i
];
63 #define ALLOC(t,p) p = (t*)malloc(sizeof(*p))
66 * We just ignore the last column and row
67 * If the final element is not equal to one
68 * then the result will actually be a multiple of the input
70 static void matrix2zz(Matrix
*M
, mat_ZZ
& m
, unsigned nr
, unsigned nc
)
74 for (int i
= 0; i
< nr
; ++i
) {
75 // assert(value_one_p(M->p[i][M->NbColumns - 1]));
76 for (int j
= 0; j
< nc
; ++j
) {
77 value2zz(M
->p
[i
][j
], m
[i
][j
]);
82 static void values2zz(Value
*p
, vec_ZZ
& v
, int len
)
86 for (int i
= 0; i
< len
; ++i
) {
93 static void zz2values(vec_ZZ
& v
, Value
*p
)
95 for (int i
= 0; i
< v
.length(); ++i
)
99 static void rays(mat_ZZ
& r
, Polyhedron
*C
)
101 unsigned dim
= C
->NbRays
- 1; /* don't count zero vertex */
102 assert(C
->NbRays
- 1 == C
->Dimension
);
107 for (i
= 0, c
= 0; i
< dim
; ++i
)
108 if (value_zero_p(C
->Ray
[i
][dim
+1])) {
109 for (int j
= 0; j
< dim
; ++j
) {
110 value2zz(C
->Ray
[i
][j
+1], tmp
);
117 static Matrix
* rays(Polyhedron
*C
)
119 unsigned dim
= C
->NbRays
- 1; /* don't count zero vertex */
120 assert(C
->NbRays
- 1 == C
->Dimension
);
122 Matrix
*M
= Matrix_Alloc(dim
+1, dim
+1);
126 for (i
= 0, c
= 0; i
<= dim
&& c
< dim
; ++i
)
127 if (value_zero_p(C
->Ray
[i
][dim
+1])) {
128 Vector_Copy(C
->Ray
[i
] + 1, M
->p
[c
], dim
);
129 value_set_si(M
->p
[c
++][dim
], 0);
132 value_set_si(M
->p
[dim
][dim
], 1);
137 static Matrix
* rays2(Polyhedron
*C
)
139 unsigned dim
= C
->NbRays
- 1; /* don't count zero vertex */
140 assert(C
->NbRays
- 1 == C
->Dimension
);
142 Matrix
*M
= Matrix_Alloc(dim
, dim
);
146 for (i
= 0, c
= 0; i
<= dim
&& c
< dim
; ++i
)
147 if (value_zero_p(C
->Ray
[i
][dim
+1]))
148 Vector_Copy(C
->Ray
[i
] + 1, M
->p
[c
++], dim
);
155 * Returns the largest absolute value in the vector
157 static ZZ
max(vec_ZZ
& v
)
160 for (int i
= 1; i
< v
.length(); ++i
)
170 Rays
= Matrix_Copy(M
);
173 cone(Polyhedron
*C
) {
174 Cone
= Polyhedron_Copy(C
);
180 matrix2zz(Rays
, A
, Rays
->NbRows
- 1, Rays
->NbColumns
- 1);
181 det
= determinant(A
);
188 Vector
* short_vector(vec_ZZ
& lambda
) {
189 Matrix
*M
= Matrix_Copy(Rays
);
190 Matrix
*inv
= Matrix_Alloc(M
->NbRows
, M
->NbColumns
);
191 int ok
= Matrix_Inverse(M
, inv
);
198 matrix2zz(inv
, B
, inv
->NbRows
- 1, inv
->NbColumns
- 1);
199 long r
= LLL(det2
, B
, U
);
203 for (int i
= 1; i
< B
.NumRows(); ++i
) {
215 Vector
*z
= Vector_Alloc(U
[index
].length()+1);
217 zz2values(U
[index
], z
->p
);
218 value_set_si(z
->p
[U
[index
].length()], 0);
222 Polyhedron
*C
= poly();
224 for (i
= 0; i
< C
->NbConstraints
; ++i
) {
225 Inner_Product(z
->p
, C
->Constraint
[i
]+1, z
->Size
-1, &tmp
);
226 if (value_pos_p(tmp
))
229 if (i
== C
->NbConstraints
) {
230 value_set_si(tmp
, -1);
231 Vector_Scale(z
->p
, z
->p
, tmp
, z
->Size
-1);
238 Polyhedron_Free(Cone
);
244 Matrix
*M
= Matrix_Alloc(Rays
->NbRows
+1, Rays
->NbColumns
+1);
245 for (int i
= 0; i
< Rays
->NbRows
; ++i
) {
246 Vector_Copy(Rays
->p
[i
], M
->p
[i
]+1, Rays
->NbColumns
);
247 value_set_si(M
->p
[i
][0], 1);
249 Vector_Set(M
->p
[Rays
->NbRows
]+1, 0, Rays
->NbColumns
-1);
250 value_set_si(M
->p
[Rays
->NbRows
][0], 1);
251 value_set_si(M
->p
[Rays
->NbRows
][Rays
->NbColumns
], 1);
252 Cone
= Rays2Polyhedron(M
, M
->NbRows
+1);
253 assert(Cone
->NbConstraints
== Cone
->NbRays
);
267 dpoly(int d
, ZZ
& degree
, int offset
= 0) {
268 coeff
.SetLength(d
+1);
270 int min
= d
+ offset
;
271 if (degree
>= 0 && degree
< ZZ(INIT_VAL
, min
))
272 min
= to_int(degree
);
274 ZZ c
= ZZ(INIT_VAL
, 1);
277 for (int i
= 1; i
<= min
; ++i
) {
278 c
*= (degree
-i
+ 1);
283 void operator *= (dpoly
& f
) {
284 assert(coeff
.length() == f
.coeff
.length());
286 coeff
= f
.coeff
[0] * coeff
;
287 for (int i
= 1; i
< coeff
.length(); ++i
)
288 for (int j
= 0; i
+j
< coeff
.length(); ++j
)
289 coeff
[i
+j
] += f
.coeff
[i
] * old
[j
];
291 void div(dpoly
& d
, mpq_t count
, ZZ
& sign
) {
292 int len
= coeff
.length();
295 mpq_t
* c
= new mpq_t
[coeff
.length()];
298 for (int i
= 0; i
< len
; ++i
) {
300 zz2value(coeff
[i
], tmp
);
301 mpq_set_z(c
[i
], tmp
);
303 for (int j
= 1; j
<= i
; ++j
) {
304 zz2value(d
.coeff
[j
], tmp
);
305 mpq_set_z(qtmp
, tmp
);
306 mpq_mul(qtmp
, qtmp
, c
[i
-j
]);
307 mpq_sub(c
[i
], c
[i
], qtmp
);
310 zz2value(d
.coeff
[0], tmp
);
311 mpq_set_z(qtmp
, tmp
);
312 mpq_div(c
[i
], c
[i
], qtmp
);
315 mpq_sub(count
, count
, c
[len
-1]);
317 mpq_add(count
, count
, c
[len
-1]);
321 for (int i
= 0; i
< len
; ++i
)
333 dpoly_n(int d
, ZZ
& degree_0
, ZZ
& degree_1
, int offset
= 0) {
337 zz2value(degree_0
, d0
);
338 zz2value(degree_1
, d1
);
339 coeff
= Matrix_Alloc(d
+1, d
+1+1);
340 value_set_si(coeff
->p
[0][0], 1);
341 value_set_si(coeff
->p
[0][d
+1], 1);
342 for (int i
= 1; i
<= d
; ++i
) {
343 value_multiply(coeff
->p
[i
][0], coeff
->p
[i
-1][0], d0
);
344 Vector_Combine(coeff
->p
[i
-1], coeff
->p
[i
-1]+1, coeff
->p
[i
]+1,
346 value_set_si(coeff
->p
[i
][d
+1], i
);
347 value_multiply(coeff
->p
[i
][d
+1], coeff
->p
[i
][d
+1], coeff
->p
[i
-1][d
+1]);
348 value_decrement(d0
, d0
);
353 void div(dpoly
& d
, Vector
*count
, ZZ
& sign
) {
354 int len
= coeff
->NbRows
;
355 Matrix
* c
= Matrix_Alloc(coeff
->NbRows
, coeff
->NbColumns
);
358 for (int i
= 0; i
< len
; ++i
) {
359 Vector_Copy(coeff
->p
[i
], c
->p
[i
], len
+1);
360 for (int j
= 1; j
<= i
; ++j
) {
361 zz2value(d
.coeff
[j
], tmp
);
362 value_multiply(tmp
, tmp
, c
->p
[i
][len
]);
363 value_oppose(tmp
, tmp
);
364 Vector_Combine(c
->p
[i
], c
->p
[i
-j
], c
->p
[i
],
365 c
->p
[i
-j
][len
], tmp
, len
);
366 value_multiply(c
->p
[i
][len
], c
->p
[i
][len
], c
->p
[i
-j
][len
]);
368 zz2value(d
.coeff
[0], tmp
);
369 value_multiply(c
->p
[i
][len
], c
->p
[i
][len
], tmp
);
372 value_set_si(tmp
, -1);
373 Vector_Scale(c
->p
[len
-1], count
->p
, tmp
, len
);
374 value_assign(count
->p
[len
], c
->p
[len
-1][len
]);
376 Vector_Copy(c
->p
[len
-1], count
->p
, len
+1);
377 Vector_Normalize(count
->p
, len
+1);
384 * Barvinok's Decomposition of a simplicial cone
386 * Returns two lists of polyhedra
388 void barvinok_decompose(Polyhedron
*C
, Polyhedron
**ppos
, Polyhedron
**pneg
)
390 Polyhedron
*pos
= *ppos
, *neg
= *pneg
;
391 vector
<cone
*> nonuni
;
392 cone
* c
= new cone(C
);
399 Polyhedron
*p
= Polyhedron_Copy(c
->Cone
);
405 while (!nonuni
.empty()) {
408 Vector
* v
= c
->short_vector(lambda
);
409 for (int i
= 0; i
< c
->Rays
->NbRows
- 1; ++i
) {
412 Matrix
* M
= Matrix_Copy(c
->Rays
);
413 Vector_Copy(v
->p
, M
->p
[i
], v
->Size
);
414 cone
* pc
= new cone(M
);
415 assert (pc
->det
!= 0);
416 if (abs(pc
->det
) > 1) {
417 assert(abs(pc
->det
) < abs(c
->det
));
418 nonuni
.push_back(pc
);
420 Polyhedron
*p
= pc
->poly();
422 if (sign(pc
->det
) == s
) {
441 * Returns a single list of npos "positive" cones followed by nneg
443 * The input cone is freed
445 void decompose(Polyhedron
*cone
, Polyhedron
**parts
, int *npos
, int *nneg
, unsigned MaxRays
)
447 Polyhedron_Polarize(cone
);
448 if (cone
->NbRays
- 1 != cone
->Dimension
) {
449 Polyhedron
*tmp
= cone
;
450 cone
= triangularize_cone(cone
, MaxRays
);
451 Polyhedron_Free(tmp
);
453 Polyhedron
*polpos
= NULL
, *polneg
= NULL
;
454 *npos
= 0; *nneg
= 0;
455 for (Polyhedron
*Polar
= cone
; Polar
; Polar
= Polar
->next
)
456 barvinok_decompose(Polar
, &polpos
, &polneg
);
459 for (Polyhedron
*i
= polpos
; i
; i
= i
->next
) {
460 Polyhedron_Polarize(i
);
464 for (Polyhedron
*i
= polneg
; i
; i
= i
->next
) {
465 Polyhedron_Polarize(i
);
476 const int MAX_TRY
=10;
478 * Searches for a vector that is not othogonal to any
479 * of the rays in rays.
481 static void nonorthog(mat_ZZ
& rays
, vec_ZZ
& lambda
)
483 int dim
= rays
.NumCols();
485 lambda
.SetLength(dim
);
486 for (int i
= 2; !found
&& i
<= 50*dim
; i
+=4) {
487 for (int j
= 0; j
< MAX_TRY
; ++j
) {
488 for (int k
= 0; k
< dim
; ++k
) {
489 int r
= random_int(i
)+2;
490 int v
= (2*(r
%2)-1) * (r
>> 1);
494 for (; k
< rays
.NumRows(); ++k
)
495 if (lambda
* rays
[k
] == 0)
497 if (k
== rays
.NumRows()) {
506 static void add_rays(mat_ZZ
& rays
, Polyhedron
*i
, int *r
)
508 unsigned dim
= i
->Dimension
;
509 for (int k
= 0; k
< i
->NbRays
; ++k
) {
510 if (!value_zero_p(i
->Ray
[k
][dim
+1]))
512 values2zz(i
->Ray
[k
]+1, rays
[(*r
)++], dim
);
516 void lattice_point(Value
* values
, Polyhedron
*i
, vec_ZZ
& lambda
, ZZ
& num
)
519 unsigned dim
= i
->Dimension
;
520 if(!value_one_p(values
[dim
])) {
521 Matrix
* Rays
= rays(i
);
522 Matrix
*inv
= Matrix_Alloc(Rays
->NbRows
, Rays
->NbColumns
);
523 int ok
= Matrix_Inverse(Rays
, inv
);
527 Vector
*lambda
= Vector_Alloc(dim
+1);
528 Vector_Matrix_Product(values
, inv
, lambda
->p
);
530 for (int j
= 0; j
< dim
; ++j
)
531 mpz_cdiv_q(lambda
->p
[j
], lambda
->p
[j
], lambda
->p
[dim
]);
532 value_set_si(lambda
->p
[dim
], 1);
533 Vector
*A
= Vector_Alloc(dim
+1);
534 Vector_Matrix_Product(lambda
->p
, Rays
, A
->p
);
537 values2zz(A
->p
, vertex
, dim
);
540 values2zz(values
, vertex
, dim
);
542 num
= vertex
* lambda
;
545 static evalue
*term(int param
, ZZ
& c
, Value
*den
= NULL
)
547 evalue
*EP
= new evalue();
549 value_set_si(EP
->d
,0);
550 EP
->x
.p
= new_enode(polynomial
, 2, param
+ 1);
551 evalue_set_si(&EP
->x
.p
->arr
[0], 0, 1);
552 value_init(EP
->x
.p
->arr
[1].x
.n
);
554 value_set_si(EP
->x
.p
->arr
[1].d
, 1);
556 value_assign(EP
->x
.p
->arr
[1].d
, *den
);
557 zz2value(c
, EP
->x
.p
->arr
[1].x
.n
);
561 static void vertex_period(
562 Polyhedron
*i
, vec_ZZ
& lambda
, Matrix
*T
,
563 Value lcm
, int p
, Vector
*val
,
564 evalue
*E
, evalue
* ev
,
567 unsigned nparam
= T
->NbRows
- 1;
568 unsigned dim
= i
->Dimension
;
574 Vector
* values
= Vector_Alloc(dim
+ 1);
575 Vector_Matrix_Product(val
->p
, T
, values
->p
);
576 value_assign(values
->p
[dim
], lcm
);
577 lattice_point(values
->p
, i
, lambda
, num
);
582 zz2value(num
, ev
->x
.n
);
583 value_assign(ev
->d
, lcm
);
590 values2zz(T
->p
[p
], vertex
, dim
);
591 nump
= vertex
* lambda
;
592 if (First_Non_Zero(val
->p
, p
) == -1) {
593 value_assign(tmp
, lcm
);
594 evalue
*ET
= term(p
, nump
, &tmp
);
596 free_evalue_refs(ET
);
600 value_assign(tmp
, lcm
);
601 if (First_Non_Zero(T
->p
[p
], dim
) != -1)
602 Vector_Gcd(T
->p
[p
], dim
, &tmp
);
604 if (value_lt(tmp
, lcm
)) {
607 value_division(tmp
, lcm
, tmp
);
608 value_set_si(ev
->d
, 0);
609 ev
->x
.p
= new_enode(periodic
, VALUE_TO_INT(tmp
), p
+1);
610 value2zz(tmp
, count
);
612 value_decrement(tmp
, tmp
);
614 ZZ new_offset
= offset
- count
* nump
;
615 value_assign(val
->p
[p
], tmp
);
616 vertex_period(i
, lambda
, T
, lcm
, p
+1, val
, E
,
617 &ev
->x
.p
->arr
[VALUE_TO_INT(tmp
)], new_offset
);
618 } while (value_pos_p(tmp
));
620 vertex_period(i
, lambda
, T
, lcm
, p
+1, val
, E
, ev
, offset
);
624 static void mask_r(Matrix
*f
, int nr
, Vector
*lcm
, int p
, Vector
*val
, evalue
*ev
)
626 unsigned nparam
= lcm
->Size
;
629 Vector
* prod
= Vector_Alloc(f
->NbRows
);
630 Matrix_Vector_Product(f
, val
->p
, prod
->p
);
632 for (int i
= 0; i
< nr
; ++i
) {
633 value_modulus(prod
->p
[i
], prod
->p
[i
], f
->p
[i
][nparam
+1]);
634 isint
&= value_zero_p(prod
->p
[i
]);
636 value_set_si(ev
->d
, 1);
638 value_set_si(ev
->x
.n
, isint
);
645 if (value_one_p(lcm
->p
[p
]))
646 mask_r(f
, nr
, lcm
, p
+1, val
, ev
);
648 value_assign(tmp
, lcm
->p
[p
]);
649 value_set_si(ev
->d
, 0);
650 ev
->x
.p
= new_enode(periodic
, VALUE_TO_INT(tmp
), p
+1);
652 value_decrement(tmp
, tmp
);
653 value_assign(val
->p
[p
], tmp
);
654 mask_r(f
, nr
, lcm
, p
+1, val
, &ev
->x
.p
->arr
[VALUE_TO_INT(tmp
)]);
655 } while (value_pos_p(tmp
));
660 static evalue
*multi_monom(vec_ZZ
& p
)
662 evalue
*X
= new evalue();
665 unsigned nparam
= p
.length()-1;
666 zz2value(p
[nparam
], X
->x
.n
);
667 value_set_si(X
->d
, 1);
668 for (int i
= 0; i
< nparam
; ++i
) {
671 evalue
*T
= term(i
, p
[i
]);
680 * Check whether mapping polyhedron P on the affine combination
681 * num yields a range that has a fixed quotient on integer
683 * If zero is true, then we are only interested in the quotient
684 * for the cases where the remainder is zero.
685 * Returns NULL if false and a newly allocated value if true.
687 static Value
*fixed_quotient(Polyhedron
*P
, vec_ZZ
& num
, Value d
, bool zero
)
690 int len
= num
.length();
691 Matrix
*T
= Matrix_Alloc(2, len
);
692 zz2values(num
, T
->p
[0]);
693 value_set_si(T
->p
[1][len
-1], 1);
694 Polyhedron
*I
= Polyhedron_Image(P
, T
, P
->NbConstraints
);
698 for (i
= 0; i
< I
->NbRays
; ++i
)
699 if (value_zero_p(I
->Ray
[i
][2])) {
707 int bounded
= line_minmax(I
, &min
, &max
);
711 mpz_cdiv_q(min
, min
, d
);
713 mpz_fdiv_q(min
, min
, d
);
714 mpz_fdiv_q(max
, max
, d
);
716 if (value_eq(min
, max
)) {
719 value_assign(*ret
, min
);
727 * Normalize linear expression coef modulo m
728 * Removes common factor and reduces coefficients
729 * Returns index of first non-zero coefficient or len
731 static int normal_mod(Value
*coef
, int len
, Value
*m
)
736 Vector_Gcd(coef
, len
, &gcd
);
738 Vector_AntiScale(coef
, coef
, gcd
, len
);
740 value_division(*m
, *m
, gcd
);
747 for (j
= 0; j
< len
; ++j
)
748 mpz_fdiv_r(coef
[j
], coef
[j
], *m
);
749 for (j
= 0; j
< len
; ++j
)
750 if (value_notzero_p(coef
[j
]))
757 static void mask(Matrix
*f
, evalue
*factor
)
759 int nr
= f
->NbRows
, nc
= f
->NbColumns
;
762 for (n
= 0; n
< nr
&& value_notzero_p(f
->p
[n
][nc
-1]); ++n
)
763 if (value_notone_p(f
->p
[n
][nc
-1]) &&
764 value_notmone_p(f
->p
[n
][nc
-1]))
778 value_set_si(EV
.x
.n
, 1);
780 for (n
= 0; n
< nr
; ++n
) {
781 value_assign(m
, f
->p
[n
][nc
-1]);
782 if (value_one_p(m
) || value_mone_p(m
))
785 int j
= normal_mod(f
->p
[n
], nc
-1, &m
);
787 free_evalue_refs(factor
);
788 value_init(factor
->d
);
789 evalue_set_si(factor
, 0, 1);
793 values2zz(f
->p
[n
], row
, nc
-1);
796 if (j
< (nc
-1)-1 && row
[j
] > g
/2) {
797 for (int k
= j
; k
< (nc
-1); ++k
)
803 value_set_si(EP
.d
, 0);
804 EP
.x
.p
= new_enode(relation
, 2, 0);
805 value_clear(EP
.x
.p
->arr
[1].d
);
806 EP
.x
.p
->arr
[1] = *factor
;
807 evalue
*ev
= &EP
.x
.p
->arr
[0];
808 value_set_si(ev
->d
, 0);
809 ev
->x
.p
= new_enode(fractional
, 3, -1);
810 evalue_set_si(&ev
->x
.p
->arr
[1], 0, 1);
811 evalue_set_si(&ev
->x
.p
->arr
[2], 1, 1);
812 evalue
*E
= multi_monom(row
);
813 value_assign(EV
.d
, m
);
815 value_clear(ev
->x
.p
->arr
[0].d
);
816 ev
->x
.p
->arr
[0] = *E
;
822 free_evalue_refs(&EV
);
828 static void mask(Matrix
*f
, evalue
*factor
)
830 int nr
= f
->NbRows
, nc
= f
->NbColumns
;
833 for (n
= 0; n
< nr
&& value_notzero_p(f
->p
[n
][nc
-1]); ++n
)
834 if (value_notone_p(f
->p
[n
][nc
-1]) &&
835 value_notmone_p(f
->p
[n
][nc
-1]))
843 unsigned np
= nc
- 2;
844 Vector
*lcm
= Vector_Alloc(np
);
845 Vector
*val
= Vector_Alloc(nc
);
846 Vector_Set(val
->p
, 0, nc
);
847 value_set_si(val
->p
[np
], 1);
848 Vector_Set(lcm
->p
, 1, np
);
849 for (n
= 0; n
< nr
; ++n
) {
850 if (value_one_p(f
->p
[n
][nc
-1]) ||
851 value_mone_p(f
->p
[n
][nc
-1]))
853 for (int j
= 0; j
< np
; ++j
)
854 if (value_notzero_p(f
->p
[n
][j
])) {
855 Gcd(f
->p
[n
][j
], f
->p
[n
][nc
-1], &tmp
);
856 value_division(tmp
, f
->p
[n
][nc
-1], tmp
);
857 value_lcm(tmp
, lcm
->p
[j
], &lcm
->p
[j
]);
862 mask_r(f
, nr
, lcm
, 0, val
, &EP
);
867 free_evalue_refs(&EP
);
878 static bool mod_needed(Polyhedron
*PD
, vec_ZZ
& num
, Value d
, evalue
*E
)
880 Value
*q
= fixed_quotient(PD
, num
, d
, false);
885 value_oppose(*q
, *q
);
888 value_set_si(EV
.d
, 1);
890 value_multiply(EV
.x
.n
, *q
, d
);
892 free_evalue_refs(&EV
);
898 static void ceil_mod(Value
*coef
, int len
, Value d
, ZZ
& f
, evalue
*EP
, Polyhedron
*PD
)
904 Vector_Scale(coef
, coef
, m
, len
);
907 int j
= normal_mod(coef
, len
, &m
);
915 values2zz(coef
, num
, len
);
922 evalue_set_si(&tmp
, 0, 1);
926 while (j
< len
-1 && (num
[j
] == g
/2 || num
[j
] == 0))
928 if ((j
< len
-1 && num
[j
] > g
/2) || (j
== len
-1 && num
[j
] >= (g
+1)/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
);
958 value_set_si(EV
.x
.n
, 1);
959 value_assign(EV
.d
, m
);
962 value_set_si(EV
.d
, 0);
963 EV
.x
.p
= new_enode(fractional
, 3, -1);
964 evalue_copy(&EV
.x
.p
->arr
[0], E
);
965 evalue_set_si(&EV
.x
.p
->arr
[1], 0, 1);
966 value_init(EV
.x
.p
->arr
[2].x
.n
);
967 zz2value(f
, EV
.x
.p
->arr
[2].x
.n
);
968 value_set_si(EV
.x
.p
->arr
[2].d
, 1);
973 free_evalue_refs(&EV
);
978 free_evalue_refs(&tmp
);
984 evalue
* bv_ceil3(Value
*coef
, int len
, Value d
, Polyhedron
*P
)
986 Vector
*val
= Vector_Alloc(len
);
991 Vector_Scale(coef
, val
->p
, t
, len
);
992 value_absolute(t
, d
);
995 values2zz(val
->p
, num
, len
);
996 evalue
*EP
= multi_monom(num
);
1000 value_init(tmp
.x
.n
);
1001 value_set_si(tmp
.x
.n
, 1);
1002 value_assign(tmp
.d
, t
);
1008 ceil_mod(val
->p
, len
, t
, one
, EP
, P
);
1011 /* copy EP to malloc'ed evalue */
1017 free_evalue_refs(&tmp
);
1024 evalue
* lattice_point(
1025 Polyhedron
*i
, vec_ZZ
& lambda
, Matrix
*W
, Value lcm
, Polyhedron
*PD
)
1027 unsigned nparam
= W
->NbColumns
- 1;
1029 Matrix
* Rays
= rays2(i
);
1030 Matrix
*T
= Transpose(Rays
);
1031 Matrix
*T2
= Matrix_Copy(T
);
1032 Matrix
*inv
= Matrix_Alloc(T2
->NbRows
, T2
->NbColumns
);
1033 int ok
= Matrix_Inverse(T2
, inv
);
1038 matrix2zz(W
, vertex
, W
->NbRows
, W
->NbColumns
);
1041 num
= lambda
* vertex
;
1043 evalue
*EP
= multi_monom(num
);
1047 value_init(tmp
.x
.n
);
1048 value_set_si(tmp
.x
.n
, 1);
1049 value_assign(tmp
.d
, lcm
);
1053 Matrix
*L
= Matrix_Alloc(inv
->NbRows
, W
->NbColumns
);
1054 Matrix_Product(inv
, W
, L
);
1057 matrix2zz(T
, RT
, T
->NbRows
, T
->NbColumns
);
1060 vec_ZZ p
= lambda
* RT
;
1062 for (int i
= 0; i
< L
->NbRows
; ++i
) {
1063 ceil_mod(L
->p
[i
], nparam
+1, lcm
, p
[i
], EP
, PD
);
1069 free_evalue_refs(&tmp
);
1073 evalue
* lattice_point(
1074 Polyhedron
*i
, vec_ZZ
& lambda
, Matrix
*W
, Value lcm
, Polyhedron
*PD
)
1076 Matrix
*T
= Transpose(W
);
1077 unsigned nparam
= T
->NbRows
- 1;
1079 evalue
*EP
= new evalue();
1081 evalue_set_si(EP
, 0, 1);
1084 Vector
*val
= Vector_Alloc(nparam
+1);
1085 value_set_si(val
->p
[nparam
], 1);
1086 ZZ
offset(INIT_VAL
, 0);
1088 vertex_period(i
, lambda
, T
, lcm
, 0, val
, EP
, &ev
, offset
);
1091 free_evalue_refs(&ev
);
1102 Param_Vertices
* V
, Polyhedron
*i
, vec_ZZ
& lambda
, term_info
* term
,
1105 unsigned nparam
= V
->Vertex
->NbColumns
- 2;
1106 unsigned dim
= i
->Dimension
;
1108 vertex
.SetDims(V
->Vertex
->NbRows
, nparam
+1);
1112 value_set_si(lcm
, 1);
1113 for (int j
= 0; j
< V
->Vertex
->NbRows
; ++j
) {
1114 value_lcm(lcm
, V
->Vertex
->p
[j
][nparam
+1], &lcm
);
1116 if (value_notone_p(lcm
)) {
1117 Matrix
* mv
= Matrix_Alloc(dim
, nparam
+1);
1118 for (int j
= 0 ; j
< dim
; ++j
) {
1119 value_division(tmp
, lcm
, V
->Vertex
->p
[j
][nparam
+1]);
1120 Vector_Scale(V
->Vertex
->p
[j
], mv
->p
[j
], tmp
, nparam
+1);
1123 term
->E
= lattice_point(i
, lambda
, mv
, lcm
, PD
);
1131 for (int i
= 0; i
< V
->Vertex
->NbRows
; ++i
) {
1132 assert(value_one_p(V
->Vertex
->p
[i
][nparam
+1])); // for now
1133 values2zz(V
->Vertex
->p
[i
], vertex
[i
], nparam
+1);
1137 num
= lambda
* vertex
;
1141 for (int j
= 0; j
< nparam
; ++j
)
1147 term
->E
= multi_monom(num
);
1151 term
->constant
= num
[nparam
];
1154 term
->coeff
= num
[p
];
1161 void normalize(Polyhedron
*i
, vec_ZZ
& lambda
, ZZ
& sign
, ZZ
& num
, vec_ZZ
& den
)
1163 unsigned dim
= i
->Dimension
;
1167 rays
.SetDims(dim
, dim
);
1168 add_rays(rays
, i
, &r
);
1169 den
= rays
* lambda
;
1172 for (int j
= 0; j
< den
.length(); ++j
) {
1176 den
[j
] = abs(den
[j
]);
1184 typedef Polyhedron
* Polyhedron_p
;
1186 void barvinok_count(Polyhedron
*P
, Value
* result
, unsigned NbMaxCons
)
1188 Polyhedron
** vcone
;
1191 sign
.SetLength(ncone
);
1199 value_set_si(*result
, 0);
1203 for (; r
< P
->NbRays
; ++r
)
1204 if (value_zero_p(P
->Ray
[r
][P
->Dimension
+1]))
1206 if (P
->NbBid
!=0 || r
< P
->NbRays
) {
1207 value_set_si(*result
, -1);
1211 P
= remove_equalities(P
);
1214 value_set_si(*result
, 0);
1220 value_set_si(factor
, 1);
1221 Q
= Polyhedron_Reduce(P
, &factor
);
1228 if (P
->Dimension
== 0) {
1229 value_assign(*result
, factor
);
1232 value_clear(factor
);
1237 vcone
= new Polyhedron_p
[P
->NbRays
];
1239 for (int j
= 0; j
< P
->NbRays
; ++j
) {
1241 Polyhedron
*C
= supporting_cone(P
, j
);
1242 decompose(C
, &vcone
[j
], &npos
, &nneg
, NbMaxCons
);
1243 ncone
+= npos
+ nneg
;
1244 sign
.SetLength(ncone
);
1245 for (int k
= 0; k
< npos
; ++k
)
1246 sign
[ncone
-nneg
-k
-1] = 1;
1247 for (int k
= 0; k
< nneg
; ++k
)
1248 sign
[ncone
-k
-1] = -1;
1252 rays
.SetDims(ncone
* dim
, dim
);
1254 for (int j
= 0; j
< P
->NbRays
; ++j
) {
1255 for (Polyhedron
*i
= vcone
[j
]; i
; i
= i
->next
) {
1256 assert(i
->NbRays
-1 == dim
);
1257 add_rays(rays
, i
, &r
);
1261 nonorthog(rays
, lambda
);
1271 for (int j
= 0; j
< P
->NbRays
; ++j
) {
1272 for (Polyhedron
*i
= vcone
[j
]; i
; i
= i
->next
) {
1273 lattice_point(P
->Ray
[j
]+1, i
, vertex
);
1274 num
= vertex
* lambda
;
1275 normalize(i
, lambda
, sign
[f
], num
, den
);
1278 dpoly
n(dim
, den
[0], 1);
1279 for (int k
= 1; k
< dim
; ++k
) {
1280 dpoly
fact(dim
, den
[k
], 1);
1283 d
.div(n
, count
, sign
[f
]);
1287 Domain_Free(vcone
[j
]);
1290 assert(value_one_p(&count
[0]._mp_den
));
1291 value_multiply(*result
, &count
[0]._mp_num
, factor
);
1298 value_clear(factor
);
1301 static void uni_polynom(int param
, Vector
*c
, evalue
*EP
)
1303 unsigned dim
= c
->Size
-2;
1305 value_set_si(EP
->d
,0);
1306 EP
->x
.p
= new_enode(polynomial
, dim
+1, param
+1);
1307 for (int j
= 0; j
<= dim
; ++j
)
1308 evalue_set(&EP
->x
.p
->arr
[j
], c
->p
[j
], c
->p
[dim
+1]);
1311 static void multi_polynom(Vector
*c
, evalue
* X
, evalue
*EP
)
1313 unsigned dim
= c
->Size
-2;
1317 evalue_set(&EC
, c
->p
[dim
], c
->p
[dim
+1]);
1320 evalue_set(EP
, c
->p
[dim
], c
->p
[dim
+1]);
1322 for (int i
= dim
-1; i
>= 0; --i
) {
1324 value_assign(EC
.x
.n
, c
->p
[i
]);
1327 free_evalue_refs(&EC
);
1330 Polyhedron
*unfringe (Polyhedron
*P
, unsigned MaxRays
)
1332 int len
= P
->Dimension
+2;
1333 Polyhedron
*T
, *R
= P
;
1336 Vector
*row
= Vector_Alloc(len
);
1337 value_set_si(row
->p
[0], 1);
1339 R
= DomainConstraintSimplify(Polyhedron_Copy(P
), MaxRays
);
1341 Matrix
*M
= Matrix_Alloc(2, len
-1);
1342 value_set_si(M
->p
[1][len
-2], 1);
1343 for (int v
= 0; v
< P
->Dimension
; ++v
) {
1344 value_set_si(M
->p
[0][v
], 1);
1345 Polyhedron
*I
= Polyhedron_Image(P
, M
, 2+1);
1346 value_set_si(M
->p
[0][v
], 0);
1347 for (int r
= 0; r
< I
->NbConstraints
; ++r
) {
1348 if (value_zero_p(I
->Constraint
[r
][0]))
1350 if (value_zero_p(I
->Constraint
[r
][1]))
1352 if (value_one_p(I
->Constraint
[r
][1]))
1354 if (value_mone_p(I
->Constraint
[r
][1]))
1356 value_absolute(g
, I
->Constraint
[r
][1]);
1357 Vector_Set(row
->p
+1, 0, len
-2);
1358 value_division(row
->p
[1+v
], I
->Constraint
[r
][1], g
);
1359 mpz_fdiv_q(row
->p
[len
-1], I
->Constraint
[r
][2], g
);
1361 R
= AddConstraints(row
->p
, 1, R
, MaxRays
);
1373 static Polyhedron
*reduce_domain(Polyhedron
*D
, Matrix
*CT
, Polyhedron
*CEq
,
1374 Polyhedron
**fVD
, int nd
, unsigned MaxRays
)
1379 Dt
= CT
? DomainPreimage(D
, CT
, MaxRays
) : D
;
1380 Polyhedron
*rVD
= DomainIntersection(Dt
, CEq
, MaxRays
);
1382 /* if rVD is empty or too small in geometric dimension */
1383 if(!rVD
|| emptyQ(rVD
) ||
1384 (rVD
->Dimension
-rVD
->NbEq
< Dt
->Dimension
-Dt
->NbEq
-CEq
->NbEq
)) {
1389 return 0; /* empty validity domain */
1395 fVD
[nd
] = Domain_Copy(rVD
);
1396 for (int i
= 0 ; i
< nd
; ++i
) {
1397 Polyhedron
*I
= DomainIntersection(fVD
[nd
], fVD
[i
], MaxRays
);
1402 Polyhedron
*F
= DomainSimplify(I
, fVD
[nd
], MaxRays
);
1404 Polyhedron
*T
= rVD
;
1405 rVD
= DomainDifference(rVD
, F
, MaxRays
);
1412 rVD
= DomainConstraintSimplify(rVD
, MaxRays
);
1414 Domain_Free(fVD
[nd
]);
1421 barvinok_count(rVD
, &c
, MaxRays
);
1422 if (value_zero_p(c
)) {
1431 evalue
* barvinok_enumerate_ev(Polyhedron
*P
, Polyhedron
* C
, unsigned MaxRays
)
1433 //P = unfringe(P, MaxRays);
1434 Polyhedron
*CEq
= NULL
, *rVD
, *pVD
, *CA
;
1436 Param_Polyhedron
*PP
= NULL
;
1437 Param_Domain
*D
, *next
;
1440 unsigned nparam
= C
->Dimension
;
1442 ALLOC(evalue
, eres
);
1443 value_init(eres
->d
);
1444 value_set_si(eres
->d
, 0);
1447 value_init(factor
.d
);
1448 evalue_set_si(&factor
, 1, 1);
1450 CA
= align_context(C
, P
->Dimension
, MaxRays
);
1451 P
= DomainIntersection(P
, CA
, MaxRays
);
1452 Polyhedron_Free(CA
);
1454 if (C
->Dimension
== 0 || emptyQ(P
)) {
1456 eres
->x
.p
= new_enode(partition
, 2, C
->Dimension
);
1457 EVALUE_SET_DOMAIN(eres
->x
.p
->arr
[0],
1458 DomainConstraintSimplify(CEq
? CEq
: Polyhedron_Copy(C
), MaxRays
));
1459 value_set_si(eres
->x
.p
->arr
[1].d
, 1);
1460 value_init(eres
->x
.p
->arr
[1].x
.n
);
1462 value_set_si(eres
->x
.p
->arr
[1].x
.n
, 0);
1464 barvinok_count(P
, &eres
->x
.p
->arr
[1].x
.n
, MaxRays
);
1466 emul(&factor
, eres
);
1467 reduce_evalue(eres
);
1468 free_evalue_refs(&factor
);
1473 Param_Polyhedron_Free(PP
);
1477 for (r
= 0; r
< P
->NbRays
; ++r
)
1478 if (value_zero_p(P
->Ray
[r
][0]) ||
1479 value_zero_p(P
->Ray
[r
][P
->Dimension
+1])) {
1481 for (i
= P
->Dimension
- nparam
; i
< P
->Dimension
; ++i
)
1482 if (value_notzero_p(P
->Ray
[r
][i
+1]))
1484 if (i
>= P
->Dimension
)
1492 P
= remove_equalities_p(P
, P
->Dimension
-nparam
, &f
);
1496 if (P
->Dimension
== nparam
) {
1498 P
= Universe_Polyhedron(0);
1502 Polyhedron
*Q
= ParamPolyhedron_Reduce(P
, P
->Dimension
-nparam
, &factor
);
1505 if (Q
->Dimension
== nparam
) {
1507 P
= Universe_Polyhedron(0);
1512 Polyhedron
*oldP
= P
;
1513 PP
= Polyhedron2Param_SimplifiedDomain(&P
,C
,MaxRays
,&CEq
,&CT
);
1515 Polyhedron_Free(oldP
);
1517 if (isIdentity(CT
)) {
1521 assert(CT
->NbRows
!= CT
->NbColumns
);
1522 if (CT
->NbRows
== 1) // no more parameters
1524 nparam
= CT
->NbRows
- 1;
1527 unsigned dim
= P
->Dimension
- nparam
;
1528 Polyhedron
** vcone
= new Polyhedron_p
[PP
->nbV
];
1529 int * npos
= new int[PP
->nbV
];
1530 int * nneg
= new int[PP
->nbV
];
1534 for (i
= 0, V
= PP
->V
; V
; ++i
, V
= V
->next
) {
1535 Polyhedron
*C
= supporting_cone_p(P
, V
);
1536 decompose(C
, &vcone
[i
], &npos
[i
], &nneg
[i
], MaxRays
);
1539 Vector
*c
= Vector_Alloc(dim
+2);
1542 for (nd
= 0, D
=PP
->D
; D
; ++nd
, D
=D
->next
);
1543 struct section
{ Polyhedron
*D
; evalue E
; };
1544 section
*s
= new section
[nd
];
1545 Polyhedron
**fVD
= new Polyhedron_p
[nd
];
1547 for(nd
= 0, D
=PP
->D
; D
; D
=next
) {
1550 Polyhedron
*rVD
= reduce_domain(D
->Domain
, CT
, CEq
,
1555 pVD
= CT
? DomainImage(rVD
,CT
,MaxRays
) : rVD
;
1558 sign
.SetLength(ncone
);
1559 FORALL_PVertex_in_ParamPolyhedron(V
,D
,PP
) // _i is internal counter
1560 ncone
+= npos
[_i
] + nneg
[_i
];
1561 sign
.SetLength(ncone
);
1562 for (int k
= 0; k
< npos
[_i
]; ++k
)
1563 sign
[ncone
-nneg
[_i
]-k
-1] = 1;
1564 for (int k
= 0; k
< nneg
[_i
]; ++k
)
1565 sign
[ncone
-k
-1] = -1;
1566 END_FORALL_PVertex_in_ParamPolyhedron
;
1569 rays
.SetDims(ncone
* dim
, dim
);
1571 FORALL_PVertex_in_ParamPolyhedron(V
,D
,PP
) // _i is internal counter
1572 for (Polyhedron
*i
= vcone
[_i
]; i
; i
= i
->next
) {
1573 assert(i
->NbRays
-1 == dim
);
1574 add_rays(rays
, i
, &r
);
1576 END_FORALL_PVertex_in_ParamPolyhedron
;
1578 nonorthog(rays
, lambda
);
1581 den
.SetDims(ncone
,dim
);
1582 term_info
*num
= new term_info
[ncone
];
1585 FORALL_PVertex_in_ParamPolyhedron(V
,D
,PP
)
1586 for (Polyhedron
*i
= vcone
[_i
]; i
; i
= i
->next
) {
1587 lattice_point(V
, i
, lambda
, &num
[f
], pVD
);
1588 normalize(i
, lambda
, sign
[f
], num
[f
].constant
, den
[f
]);
1591 END_FORALL_PVertex_in_ParamPolyhedron
;
1592 ZZ min
= num
[0].constant
;
1593 for (int j
= 1; j
< ncone
; ++j
)
1594 if (num
[j
].constant
< min
)
1595 min
= num
[j
].constant
;
1596 for (int j
= 0; j
< ncone
; ++j
)
1597 num
[j
].constant
-= min
;
1599 value_init(s
[nd
].E
.d
);
1600 evalue_set_si(&s
[nd
].E
, 0, 1);
1603 FORALL_PVertex_in_ParamPolyhedron(V
,D
,PP
)
1604 for (Polyhedron
*i
= vcone
[_i
]; i
; i
= i
->next
) {
1605 dpoly
n(dim
, den
[f
][0], 1);
1606 for (int k
= 1; k
< dim
; ++k
) {
1607 dpoly
fact(dim
, den
[f
][k
], 1);
1610 if (num
[f
].E
!= NULL
) {
1611 ZZ
one(INIT_VAL
, 1);
1612 dpoly_n
d(dim
, num
[f
].constant
, one
);
1613 d
.div(n
, c
, sign
[f
]);
1615 multi_polynom(c
, num
[f
].E
, &EV
);
1616 eadd(&EV
, &s
[nd
].E
);
1617 free_evalue_refs(&EV
);
1618 free_evalue_refs(num
[f
].E
);
1620 } else if (num
[f
].pos
!= -1) {
1621 dpoly_n
d(dim
, num
[f
].constant
, num
[f
].coeff
);
1622 d
.div(n
, c
, sign
[f
]);
1624 uni_polynom(num
[f
].pos
, c
, &EV
);
1625 eadd(&EV
, &s
[nd
].E
);
1626 free_evalue_refs(&EV
);
1628 mpq_set_si(count
, 0, 1);
1629 dpoly
d(dim
, num
[f
].constant
);
1630 d
.div(n
, count
, sign
[f
]);
1633 evalue_set(&EV
, &count
[0]._mp_num
, &count
[0]._mp_den
);
1634 eadd(&EV
, &s
[nd
].E
);
1635 free_evalue_refs(&EV
);
1639 END_FORALL_PVertex_in_ParamPolyhedron
;
1645 addeliminatedparams_evalue(&s
[nd
].E
, CT
);
1653 evalue_set_si(eres
, 0, 1);
1655 eres
->x
.p
= new_enode(partition
, 2*nd
, C
->Dimension
);
1656 for (int j
= 0; j
< nd
; ++j
) {
1657 EVALUE_SET_DOMAIN(eres
->x
.p
->arr
[2*j
], s
[j
].D
);
1658 value_clear(eres
->x
.p
->arr
[2*j
+1].d
);
1659 eres
->x
.p
->arr
[2*j
+1] = s
[j
].E
;
1660 Domain_Free(fVD
[j
]);
1668 for (int j
= 0; j
< PP
->nbV
; ++j
)
1669 Domain_Free(vcone
[j
]);
1675 Polyhedron_Free(CEq
);
1680 Enumeration
* barvinok_enumerate(Polyhedron
*P
, Polyhedron
* C
, unsigned MaxRays
)
1682 evalue
*EP
= barvinok_enumerate_ev(P
, C
, MaxRays
);
1684 return partition2enumeration(EP
);
1687 static void SwapColumns(Value
**V
, int n
, int i
, int j
)
1689 for (int r
= 0; r
< n
; ++r
)
1690 value_swap(V
[r
][i
], V
[r
][j
]);
1693 static void SwapColumns(Polyhedron
*P
, int i
, int j
)
1695 SwapColumns(P
->Constraint
, P
->NbConstraints
, i
, j
);
1696 SwapColumns(P
->Ray
, P
->NbRays
, i
, j
);
1699 static void negative_test_constraint(Value
*l
, Value
*u
, Value
*c
, int pos
,
1702 value_oppose(*v
, u
[pos
+1]);
1703 Vector_Combine(l
+1, u
+1, c
+1, *v
, l
[pos
+1], len
-1);
1704 value_multiply(*v
, *v
, l
[pos
+1]);
1705 value_substract(c
[len
-1], c
[len
-1], *v
);
1706 value_set_si(*v
, -1);
1707 Vector_Scale(c
+1, c
+1, *v
, len
-1);
1708 value_decrement(c
[len
-1], c
[len
-1]);
1709 ConstraintSimplify(c
, c
, len
, v
);
1712 static void oppose_constraint(Value
*c
, int len
, Value
*v
)
1714 value_set_si(*v
, -1);
1715 Vector_Scale(c
+1, c
+1, *v
, len
-1);
1716 value_decrement(c
[len
-1], c
[len
-1]);
1719 static bool SplitOnConstraint(Polyhedron
*P
, int i
, int l
, int u
,
1720 int nvar
, int len
, int exist
, int MaxRays
,
1721 Vector
*row
, Value
& f
, bool independent
,
1722 Polyhedron
**pos
, Polyhedron
**neg
)
1724 negative_test_constraint(P
->Constraint
[l
], P
->Constraint
[u
],
1725 row
->p
, nvar
+i
, len
, &f
);
1726 *neg
= AddConstraints(row
->p
, 1, P
, MaxRays
);
1728 /* We found an independent, but useless constraint
1729 * Maybe we should detect this earlier and not
1730 * mark the variable as INDEPENDENT
1732 if (emptyQ((*neg
))) {
1733 Polyhedron_Free(*neg
);
1737 oppose_constraint(row
->p
, len
, &f
);
1738 *pos
= AddConstraints(row
->p
, 1, P
, MaxRays
);
1740 if (emptyQ((*pos
))) {
1741 Polyhedron_Free(*neg
);
1742 Polyhedron_Free(*pos
);
1750 * unimodularly transform P such that constraint r is transformed
1751 * into a constraint that involves only a single (the first)
1752 * existential variable
1755 static Polyhedron
*rotate_along(Polyhedron
*P
, int r
, int nvar
, int exist
,
1761 Vector
*row
= Vector_Alloc(exist
);
1762 Vector_Copy(P
->Constraint
[r
]+1+nvar
, row
->p
, exist
);
1763 Vector_Gcd(row
->p
, exist
, &g
);
1764 if (value_notone_p(g
))
1765 Vector_AntiScale(row
->p
, row
->p
, g
, exist
);
1768 Matrix
*M
= unimodular_complete(row
);
1769 Matrix
*M2
= Matrix_Alloc(P
->Dimension
+1, P
->Dimension
+1);
1770 for (r
= 0; r
< nvar
; ++r
)
1771 value_set_si(M2
->p
[r
][r
], 1);
1772 for ( ; r
< nvar
+exist
; ++r
)
1773 Vector_Copy(M
->p
[r
-nvar
], M2
->p
[r
]+nvar
, exist
);
1774 for ( ; r
< P
->Dimension
+1; ++r
)
1775 value_set_si(M2
->p
[r
][r
], 1);
1776 Polyhedron
*T
= Polyhedron_Image(P
, M2
, MaxRays
);
1785 static bool SplitOnVar(Polyhedron
*P
, int i
,
1786 int nvar
, int len
, int exist
, int MaxRays
,
1787 Vector
*row
, Value
& f
, bool independent
,
1788 Polyhedron
**pos
, Polyhedron
**neg
)
1792 for (int l
= P
->NbEq
; l
< P
->NbConstraints
; ++l
) {
1793 if (value_negz_p(P
->Constraint
[l
][nvar
+i
+1]))
1797 for (j
= 0; j
< exist
; ++j
)
1798 if (j
!= i
&& value_notzero_p(P
->Constraint
[l
][nvar
+j
+1]))
1804 for (int u
= P
->NbEq
; u
< P
->NbConstraints
; ++u
) {
1805 if (value_posz_p(P
->Constraint
[u
][nvar
+i
+1]))
1809 for (j
= 0; j
< exist
; ++j
)
1810 if (j
!= i
&& value_notzero_p(P
->Constraint
[u
][nvar
+j
+1]))
1816 if (SplitOnConstraint(P
, i
, l
, u
,
1817 nvar
, len
, exist
, MaxRays
,
1818 row
, f
, independent
,
1822 SwapColumns(*neg
, nvar
+1, nvar
+1+i
);
1832 static bool double_bound_pair(Polyhedron
*P
, int nvar
, int exist
,
1833 int i
, int l1
, int l2
,
1834 Polyhedron
**pos
, Polyhedron
**neg
)
1838 Vector
*row
= Vector_Alloc(P
->Dimension
+2);
1839 value_set_si(row
->p
[0], 1);
1840 value_oppose(f
, P
->Constraint
[l1
][nvar
+i
+1]);
1841 Vector_Combine(P
->Constraint
[l1
]+1, P
->Constraint
[l2
]+1,
1843 P
->Constraint
[l2
][nvar
+i
+1], f
,
1845 ConstraintSimplify(row
->p
, row
->p
, P
->Dimension
+2, &f
);
1846 *pos
= AddConstraints(row
->p
, 1, P
, 0);
1847 value_set_si(f
, -1);
1848 Vector_Scale(row
->p
+1, row
->p
+1, f
, P
->Dimension
+1);
1849 value_decrement(row
->p
[P
->Dimension
+1], row
->p
[P
->Dimension
+1]);
1850 *neg
= AddConstraints(row
->p
, 1, P
, 0);
1854 return !emptyQ((*pos
)) && !emptyQ((*neg
));
1857 static bool double_bound(Polyhedron
*P
, int nvar
, int exist
,
1858 Polyhedron
**pos
, Polyhedron
**neg
)
1860 for (int i
= 0; i
< exist
; ++i
) {
1862 for (l1
= P
->NbEq
; l1
< P
->NbConstraints
; ++l1
) {
1863 if (value_negz_p(P
->Constraint
[l1
][nvar
+i
+1]))
1865 for (l2
= l1
+ 1; l2
< P
->NbConstraints
; ++l2
) {
1866 if (value_negz_p(P
->Constraint
[l2
][nvar
+i
+1]))
1868 if (double_bound_pair(P
, nvar
, exist
, i
, l1
, l2
, pos
, neg
))
1872 for (l1
= P
->NbEq
; l1
< P
->NbConstraints
; ++l1
) {
1873 if (value_posz_p(P
->Constraint
[l1
][nvar
+i
+1]))
1875 if (l1
< P
->NbConstraints
)
1876 for (l2
= l1
+ 1; l2
< P
->NbConstraints
; ++l2
) {
1877 if (value_posz_p(P
->Constraint
[l2
][nvar
+i
+1]))
1879 if (double_bound_pair(P
, nvar
, exist
, i
, l1
, l2
, pos
, neg
))
1891 INDEPENDENT
= 1 << 2
1894 static evalue
* enumerate_or(Polyhedron
*D
,
1895 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
1898 fprintf(stderr
, "\nER: Or\n");
1899 #endif /* DEBUG_ER */
1901 Polyhedron
*N
= D
->next
;
1904 barvinok_enumerate_e(D
, exist
, nparam
, MaxRays
);
1907 for (D
= N
; D
; D
= N
) {
1912 barvinok_enumerate_e(D
, exist
, nparam
, MaxRays
);
1915 free_evalue_refs(EN
);
1925 static evalue
* enumerate_sum(Polyhedron
*P
,
1926 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
1928 int nvar
= P
->Dimension
- exist
- nparam
;
1929 int toswap
= nvar
< exist
? nvar
: exist
;
1930 for (int i
= 0; i
< toswap
; ++i
)
1931 SwapColumns(P
, 1 + i
, nvar
+exist
- i
);
1935 fprintf(stderr
, "\nER: Sum\n");
1936 #endif /* DEBUG_ER */
1938 evalue
*EP
= barvinok_enumerate_e(P
, exist
, nparam
, MaxRays
);
1940 for (int i
= 0; i
< /* nvar */ nparam
; ++i
) {
1941 Matrix
*C
= Matrix_Alloc(1, 1 + nparam
+ 1);
1942 value_set_si(C
->p
[0][0], 1);
1944 value_init(split
.d
);
1945 value_set_si(split
.d
, 0);
1946 split
.x
.p
= new_enode(partition
, 4, nparam
);
1947 value_set_si(C
->p
[0][1+i
], 1);
1948 Matrix
*C2
= Matrix_Copy(C
);
1949 EVALUE_SET_DOMAIN(split
.x
.p
->arr
[0],
1950 Constraints2Polyhedron(C2
, MaxRays
));
1952 evalue_set_si(&split
.x
.p
->arr
[1], 1, 1);
1953 value_set_si(C
->p
[0][1+i
], -1);
1954 value_set_si(C
->p
[0][1+nparam
], -1);
1955 EVALUE_SET_DOMAIN(split
.x
.p
->arr
[2],
1956 Constraints2Polyhedron(C
, MaxRays
));
1957 evalue_set_si(&split
.x
.p
->arr
[3], 1, 1);
1959 free_evalue_refs(&split
);
1963 evalue_range_reduction(EP
);
1965 evalue_frac2floor(EP
);
1967 evalue
*sum
= esum(EP
, nvar
);
1969 free_evalue_refs(EP
);
1973 evalue_range_reduction(EP
);
1978 static evalue
* split_sure(Polyhedron
*P
, Polyhedron
*S
,
1979 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
1981 int nvar
= P
->Dimension
- exist
- nparam
;
1983 Matrix
*M
= Matrix_Alloc(exist
, S
->Dimension
+2);
1984 for (int i
= 0; i
< exist
; ++i
)
1985 value_set_si(M
->p
[i
][nvar
+i
+1], 1);
1987 S
= DomainAddRays(S
, M
, MaxRays
);
1989 Polyhedron
*F
= DomainAddRays(P
, M
, MaxRays
);
1990 Polyhedron
*D
= DomainDifference(F
, S
, MaxRays
);
1992 D
= Disjoint_Domain(D
, 0, MaxRays
);
1997 M
= Matrix_Alloc(P
->Dimension
+1-exist
, P
->Dimension
+1);
1998 for (int j
= 0; j
< nvar
; ++j
)
1999 value_set_si(M
->p
[j
][j
], 1);
2000 for (int j
= 0; j
< nparam
+1; ++j
)
2001 value_set_si(M
->p
[nvar
+j
][nvar
+exist
+j
], 1);
2002 Polyhedron
*T
= Polyhedron_Image(S
, M
, MaxRays
);
2003 evalue
*EP
= barvinok_enumerate_e(T
, 0, nparam
, MaxRays
);
2008 for (Polyhedron
*Q
= D
; Q
; Q
= Q
->next
) {
2009 Polyhedron
*N
= Q
->next
;
2011 T
= DomainIntersection(P
, Q
, MaxRays
);
2012 evalue
*E
= barvinok_enumerate_e(T
, exist
, nparam
, MaxRays
);
2014 free_evalue_refs(E
);
2023 static evalue
* enumerate_sure(Polyhedron
*P
,
2024 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
2028 int nvar
= P
->Dimension
- exist
- nparam
;
2034 for (i
= 0; i
< exist
; ++i
) {
2035 Matrix
*M
= Matrix_Alloc(S
->NbConstraints
, S
->Dimension
+2);
2037 value_set_si(lcm
, 1);
2038 for (int j
= 0; j
< S
->NbConstraints
; ++j
) {
2039 if (value_negz_p(S
->Constraint
[j
][1+nvar
+i
]))
2041 if (value_one_p(S
->Constraint
[j
][1+nvar
+i
]))
2043 value_lcm(lcm
, S
->Constraint
[j
][1+nvar
+i
], &lcm
);
2046 for (int j
= 0; j
< S
->NbConstraints
; ++j
) {
2047 if (value_negz_p(S
->Constraint
[j
][1+nvar
+i
]))
2049 if (value_one_p(S
->Constraint
[j
][1+nvar
+i
]))
2051 value_division(f
, lcm
, S
->Constraint
[j
][1+nvar
+i
]);
2052 Vector_Scale(S
->Constraint
[j
], M
->p
[c
], f
, S
->Dimension
+2);
2053 value_substract(M
->p
[c
][S
->Dimension
+1],
2054 M
->p
[c
][S
->Dimension
+1],
2056 value_increment(M
->p
[c
][S
->Dimension
+1],
2057 M
->p
[c
][S
->Dimension
+1]);
2061 S
= AddConstraints(M
->p
[0], c
, S
, MaxRays
);
2076 fprintf(stderr
, "\nER: Sure\n");
2077 #endif /* DEBUG_ER */
2079 return split_sure(P
, S
, exist
, nparam
, MaxRays
);
2082 static evalue
* enumerate_sure2(Polyhedron
*P
,
2083 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
2085 int nvar
= P
->Dimension
- exist
- nparam
;
2087 for (r
= 0; r
< P
->NbRays
; ++r
)
2088 if (value_one_p(P
->Ray
[r
][0]) &&
2089 value_one_p(P
->Ray
[r
][P
->Dimension
+1]))
2095 Matrix
*M
= Matrix_Alloc(nvar
+ 1 + nparam
, P
->Dimension
+2);
2096 for (int i
= 0; i
< nvar
; ++i
)
2097 value_set_si(M
->p
[i
][1+i
], 1);
2098 for (int i
= 0; i
< nparam
; ++i
)
2099 value_set_si(M
->p
[i
+nvar
][1+nvar
+exist
+i
], 1);
2100 Vector_Copy(P
->Ray
[r
]+1+nvar
, M
->p
[nvar
+nparam
]+1+nvar
, exist
);
2101 value_set_si(M
->p
[nvar
+nparam
][0], 1);
2102 value_set_si(M
->p
[nvar
+nparam
][P
->Dimension
+1], 1);
2103 Polyhedron
* F
= Rays2Polyhedron(M
, MaxRays
);
2106 Polyhedron
*I
= DomainIntersection(F
, P
, MaxRays
);
2110 fprintf(stderr
, "\nER: Sure2\n");
2111 #endif /* DEBUG_ER */
2113 return split_sure(P
, I
, exist
, nparam
, MaxRays
);
2116 static evalue
* enumerate_cyclic(Polyhedron
*P
,
2117 unsigned exist
, unsigned nparam
,
2118 evalue
* EP
, int r
, int p
, unsigned MaxRays
)
2120 int nvar
= P
->Dimension
- exist
- nparam
;
2122 /* If EP in its fractional maps only contains references
2123 * to the remainder parameter with appropriate coefficients
2124 * then we could in principle avoid adding existentially
2125 * quantified variables to the validity domains.
2126 * We'd have to replace the remainder by m { p/m }
2127 * and multiply with an appropriate factor that is one
2128 * only in the appropriate range.
2129 * This last multiplication can be avoided if EP
2130 * has a single validity domain with no (further)
2131 * constraints on the remainder parameter
2134 Matrix
*CT
= Matrix_Alloc(nparam
+1, nparam
+3);
2135 Matrix
*M
= Matrix_Alloc(1, 1+nparam
+3);
2136 for (int j
= 0; j
< nparam
; ++j
)
2138 value_set_si(CT
->p
[j
][j
], 1);
2139 value_set_si(CT
->p
[p
][nparam
+1], 1);
2140 value_set_si(CT
->p
[nparam
][nparam
+2], 1);
2141 value_set_si(M
->p
[0][1+p
], -1);
2142 value_absolute(M
->p
[0][1+nparam
], P
->Ray
[0][1+nvar
+exist
+p
]);
2143 value_set_si(M
->p
[0][1+nparam
+1], 1);
2144 Polyhedron
*CEq
= Constraints2Polyhedron(M
, 1);
2146 addeliminatedparams_enum(EP
, CT
, CEq
, MaxRays
, nparam
);
2147 Polyhedron_Free(CEq
);
2153 static void enumerate_vd_add_ray(evalue
*EP
, Matrix
*Rays
, unsigned MaxRays
)
2155 if (value_notzero_p(EP
->d
))
2158 assert(EP
->x
.p
->type
== partition
);
2159 assert(EP
->x
.p
->pos
== EVALUE_DOMAIN(EP
->x
.p
->arr
[0])->Dimension
);
2160 for (int i
= 0; i
< EP
->x
.p
->size
/2; ++i
) {
2161 Polyhedron
*D
= EVALUE_DOMAIN(EP
->x
.p
->arr
[2*i
]);
2162 Polyhedron
*N
= DomainAddRays(D
, Rays
, MaxRays
);
2163 EVALUE_SET_DOMAIN(EP
->x
.p
->arr
[2*i
], N
);
2168 static evalue
* enumerate_line(Polyhedron
*P
,
2169 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
2175 fprintf(stderr
, "\nER: Line\n");
2176 #endif /* DEBUG_ER */
2178 int nvar
= P
->Dimension
- exist
- nparam
;
2180 for (i
= 0; i
< nparam
; ++i
)
2181 if (value_notzero_p(P
->Ray
[0][1+nvar
+exist
+i
]))
2184 for (j
= i
+1; j
< nparam
; ++j
)
2185 if (value_notzero_p(P
->Ray
[0][1+nvar
+exist
+i
]))
2187 assert(j
>= nparam
); // for now
2189 Matrix
*M
= Matrix_Alloc(2, P
->Dimension
+2);
2190 value_set_si(M
->p
[0][0], 1);
2191 value_set_si(M
->p
[0][1+nvar
+exist
+i
], 1);
2192 value_set_si(M
->p
[1][0], 1);
2193 value_set_si(M
->p
[1][1+nvar
+exist
+i
], -1);
2194 value_absolute(M
->p
[1][1+P
->Dimension
], P
->Ray
[0][1+nvar
+exist
+i
]);
2195 value_decrement(M
->p
[1][1+P
->Dimension
], M
->p
[1][1+P
->Dimension
]);
2196 Polyhedron
*S
= AddConstraints(M
->p
[0], 2, P
, MaxRays
);
2197 evalue
*EP
= barvinok_enumerate_e(S
, exist
, nparam
, MaxRays
);
2201 return enumerate_cyclic(P
, exist
, nparam
, EP
, 0, i
, MaxRays
);
2204 static int single_param_pos(Polyhedron
*P
, unsigned exist
, unsigned nparam
,
2207 int nvar
= P
->Dimension
- exist
- nparam
;
2208 if (First_Non_Zero(P
->Ray
[r
]+1, nvar
) != -1)
2210 int i
= First_Non_Zero(P
->Ray
[r
]+1+nvar
+exist
, nparam
);
2213 if (First_Non_Zero(P
->Ray
[r
]+1+nvar
+exist
+1, nparam
-i
-1) != -1)
2218 static evalue
* enumerate_remove_ray(Polyhedron
*P
, int r
,
2219 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
2222 fprintf(stderr
, "\nER: RedundantRay\n");
2223 #endif /* DEBUG_ER */
2227 value_set_si(one
, 1);
2228 int len
= P
->NbRays
-1;
2229 Matrix
*M
= Matrix_Alloc(2 * len
, P
->Dimension
+2);
2230 Vector_Copy(P
->Ray
[0], M
->p
[0], r
* (P
->Dimension
+2));
2231 Vector_Copy(P
->Ray
[r
+1], M
->p
[r
], (len
-r
) * (P
->Dimension
+2));
2232 for (int j
= 0; j
< P
->NbRays
; ++j
) {
2235 Vector_Combine(P
->Ray
[j
], P
->Ray
[r
], M
->p
[len
+j
-(j
>r
)],
2236 one
, P
->Ray
[j
][P
->Dimension
+1], P
->Dimension
+2);
2239 P
= Rays2Polyhedron(M
, MaxRays
);
2241 evalue
*EP
= barvinok_enumerate_e(P
, exist
, nparam
, MaxRays
);
2248 static evalue
* enumerate_redundant_ray(Polyhedron
*P
,
2249 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
2251 assert(P
->NbBid
== 0);
2252 int nvar
= P
->Dimension
- exist
- nparam
;
2256 for (int r
= 0; r
< P
->NbRays
; ++r
) {
2257 if (value_notzero_p(P
->Ray
[r
][P
->Dimension
+1]))
2259 int i1
= single_param_pos(P
, exist
, nparam
, r
);
2262 for (int r2
= r
+1; r2
< P
->NbRays
; ++r2
) {
2263 if (value_notzero_p(P
->Ray
[r2
][P
->Dimension
+1]))
2265 int i2
= single_param_pos(P
, exist
, nparam
, r2
);
2271 value_division(m
, P
->Ray
[r
][1+nvar
+exist
+i1
],
2272 P
->Ray
[r2
][1+nvar
+exist
+i1
]);
2273 value_multiply(m
, m
, P
->Ray
[r2
][1+nvar
+exist
+i1
]);
2274 /* r2 divides r => r redundant */
2275 if (value_eq(m
, P
->Ray
[r
][1+nvar
+exist
+i1
])) {
2277 return enumerate_remove_ray(P
, r
, exist
, nparam
, MaxRays
);
2280 value_division(m
, P
->Ray
[r2
][1+nvar
+exist
+i1
],
2281 P
->Ray
[r
][1+nvar
+exist
+i1
]);
2282 value_multiply(m
, m
, P
->Ray
[r
][1+nvar
+exist
+i1
]);
2283 /* r divides r2 => r2 redundant */
2284 if (value_eq(m
, P
->Ray
[r2
][1+nvar
+exist
+i1
])) {
2286 return enumerate_remove_ray(P
, r2
, exist
, nparam
, MaxRays
);
2294 static Polyhedron
*upper_bound(Polyhedron
*P
,
2295 int pos
, Value
*max
, Polyhedron
**R
)
2304 for (Polyhedron
*Q
= P
; Q
; Q
= N
) {
2306 for (r
= 0; r
< P
->NbRays
; ++r
) {
2307 if (value_zero_p(P
->Ray
[r
][P
->Dimension
+1]) &&
2308 value_pos_p(P
->Ray
[r
][1+pos
]))
2311 if (r
< P
->NbRays
) {
2319 for (r
= 0; r
< P
->NbRays
; ++r
) {
2320 if (value_zero_p(P
->Ray
[r
][P
->Dimension
+1]))
2322 mpz_fdiv_q(v
, P
->Ray
[r
][1+pos
], P
->Ray
[r
][1+P
->Dimension
]);
2323 if ((!Q
->next
&& r
== 0) || value_gt(v
, *max
))
2324 value_assign(*max
, v
);
2331 static evalue
* enumerate_ray(Polyhedron
*P
,
2332 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
2334 assert(P
->NbBid
== 0);
2335 int nvar
= P
->Dimension
- exist
- nparam
;
2338 for (r
= 0; r
< P
->NbRays
; ++r
)
2339 if (value_zero_p(P
->Ray
[r
][P
->Dimension
+1]))
2345 for (r2
= r
+1; r2
< P
->NbRays
; ++r2
)
2346 if (value_zero_p(P
->Ray
[r2
][P
->Dimension
+1]))
2348 if (r2
< P
->NbRays
) {
2350 return enumerate_sum(P
, exist
, nparam
, MaxRays
);
2354 fprintf(stderr
, "\nER: Ray\n");
2355 #endif /* DEBUG_ER */
2361 value_set_si(one
, 1);
2362 int i
= single_param_pos(P
, exist
, nparam
, r
);
2363 assert(i
!= -1); // for now;
2365 Matrix
*M
= Matrix_Alloc(P
->NbRays
, P
->Dimension
+2);
2366 for (int j
= 0; j
< P
->NbRays
; ++j
) {
2367 Vector_Combine(P
->Ray
[j
], P
->Ray
[r
], M
->p
[j
],
2368 one
, P
->Ray
[j
][P
->Dimension
+1], P
->Dimension
+2);
2370 Polyhedron
*S
= Rays2Polyhedron(M
, MaxRays
);
2372 Polyhedron
*D
= DomainDifference(P
, S
, MaxRays
);
2374 // Polyhedron_Print(stderr, P_VALUE_FMT, D);
2375 assert(value_pos_p(P
->Ray
[r
][1+nvar
+exist
+i
])); // for now
2377 D
= upper_bound(D
, nvar
+exist
+i
, &m
, &R
);
2381 M
= Matrix_Alloc(2, P
->Dimension
+2);
2382 value_set_si(M
->p
[0][0], 1);
2383 value_set_si(M
->p
[1][0], 1);
2384 value_set_si(M
->p
[0][1+nvar
+exist
+i
], -1);
2385 value_set_si(M
->p
[1][1+nvar
+exist
+i
], 1);
2386 value_assign(M
->p
[0][1+P
->Dimension
], m
);
2387 value_oppose(M
->p
[1][1+P
->Dimension
], m
);
2388 value_addto(M
->p
[1][1+P
->Dimension
], M
->p
[1][1+P
->Dimension
],
2389 P
->Ray
[r
][1+nvar
+exist
+i
]);
2390 value_decrement(M
->p
[1][1+P
->Dimension
], M
->p
[1][1+P
->Dimension
]);
2391 // Matrix_Print(stderr, P_VALUE_FMT, M);
2392 D
= AddConstraints(M
->p
[0], 2, P
, MaxRays
);
2393 // Polyhedron_Print(stderr, P_VALUE_FMT, D);
2394 value_substract(M
->p
[0][1+P
->Dimension
], M
->p
[0][1+P
->Dimension
],
2395 P
->Ray
[r
][1+nvar
+exist
+i
]);
2396 // Matrix_Print(stderr, P_VALUE_FMT, M);
2397 S
= AddConstraints(M
->p
[0], 1, P
, MaxRays
);
2398 // Polyhedron_Print(stderr, P_VALUE_FMT, S);
2401 evalue
*EP
= barvinok_enumerate_e(D
, exist
, nparam
, MaxRays
);
2406 if (value_notone_p(P
->Ray
[r
][1+nvar
+exist
+i
]))
2407 EP
= enumerate_cyclic(P
, exist
, nparam
, EP
, r
, i
, MaxRays
);
2409 M
= Matrix_Alloc(1, nparam
+2);
2410 value_set_si(M
->p
[0][0], 1);
2411 value_set_si(M
->p
[0][1+i
], 1);
2412 enumerate_vd_add_ray(EP
, M
, MaxRays
);
2417 evalue
*E
= barvinok_enumerate_e(S
, exist
, nparam
, MaxRays
);
2419 free_evalue_refs(E
);
2426 evalue
*ER
= enumerate_or(R
, exist
, nparam
, MaxRays
);
2428 free_evalue_refs(ER
);
2435 static evalue
* new_zero_ep()
2440 evalue_set_si(EP
, 0, 1);
2444 static evalue
* enumerate_vd(Polyhedron
**PA
,
2445 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
2447 Polyhedron
*P
= *PA
;
2448 int nvar
= P
->Dimension
- exist
- nparam
;
2449 Param_Polyhedron
*PP
= NULL
;
2450 Polyhedron
*C
= Universe_Polyhedron(nparam
);
2454 PP
= Polyhedron2Param_SimplifiedDomain(&PR
,C
,MaxRays
,&CEq
,&CT
);
2458 Param_Domain
*D
, *last
;
2461 for (nd
= 0, D
=PP
->D
; D
; D
=D
->next
, ++nd
)
2464 Polyhedron
**VD
= new Polyhedron_p
[nd
];
2465 Polyhedron
**fVD
= new Polyhedron_p
[nd
];
2466 for(nd
= 0, D
=PP
->D
; D
; D
=D
->next
) {
2467 Polyhedron
*rVD
= reduce_domain(D
->Domain
, CT
, CEq
,
2481 /* This doesn't seem to have any effect */
2483 Polyhedron
*CA
= align_context(VD
[0], P
->Dimension
, MaxRays
);
2485 P
= DomainIntersection(P
, CA
, MaxRays
);
2488 Polyhedron_Free(CA
);
2493 if (!EP
&& CT
->NbColumns
!= CT
->NbRows
) {
2494 Polyhedron
*CEqr
= DomainImage(CEq
, CT
, MaxRays
);
2495 Polyhedron
*CA
= align_context(CEqr
, PR
->Dimension
, MaxRays
);
2496 Polyhedron
*I
= DomainIntersection(PR
, CA
, MaxRays
);
2497 Polyhedron_Free(CEqr
);
2498 Polyhedron_Free(CA
);
2500 fprintf(stderr
, "\nER: Eliminate\n");
2501 #endif /* DEBUG_ER */
2502 nparam
-= CT
->NbColumns
- CT
->NbRows
;
2503 EP
= barvinok_enumerate_e(I
, exist
, nparam
, MaxRays
);
2504 nparam
+= CT
->NbColumns
- CT
->NbRows
;
2505 addeliminatedparams_enum(EP
, CT
, CEq
, MaxRays
, nparam
);
2509 Polyhedron_Free(PR
);
2512 if (!EP
&& nd
> 1) {
2514 fprintf(stderr
, "\nER: VD\n");
2515 #endif /* DEBUG_ER */
2516 for (int i
= 0; i
< nd
; ++i
) {
2517 Polyhedron
*CA
= align_context(VD
[i
], P
->Dimension
, MaxRays
);
2518 Polyhedron
*I
= DomainIntersection(P
, CA
, MaxRays
);
2521 EP
= barvinok_enumerate_e(I
, exist
, nparam
, MaxRays
);
2523 evalue
*E
= barvinok_enumerate_e(I
, exist
, nparam
, MaxRays
);
2525 free_evalue_refs(E
);
2529 Polyhedron_Free(CA
);
2533 for (int i
= 0; i
< nd
; ++i
) {
2534 Polyhedron_Free(VD
[i
]);
2535 Polyhedron_Free(fVD
[i
]);
2541 if (!EP
&& nvar
== 0) {
2544 Param_Vertices
*V
, *V2
;
2545 Matrix
* M
= Matrix_Alloc(1, P
->Dimension
+2);
2547 FORALL_PVertex_in_ParamPolyhedron(V
, last
, PP
) {
2549 FORALL_PVertex_in_ParamPolyhedron(V2
, last
, PP
) {
2556 for (int i
= 0; i
< exist
; ++i
) {
2557 value_oppose(f
, V
->Vertex
->p
[i
][nparam
+1]);
2558 Vector_Combine(V
->Vertex
->p
[i
],
2560 M
->p
[0] + 1 + nvar
+ exist
,
2561 V2
->Vertex
->p
[i
][nparam
+1],
2565 for (j
= 0; j
< nparam
; ++j
)
2566 if (value_notzero_p(M
->p
[0][1+nvar
+exist
+j
]))
2570 ConstraintSimplify(M
->p
[0], M
->p
[0],
2571 P
->Dimension
+2, &f
);
2572 value_set_si(M
->p
[0][0], 0);
2573 Polyhedron
*para
= AddConstraints(M
->p
[0], 1, P
,
2576 Polyhedron_Free(para
);
2579 Polyhedron
*pos
, *neg
;
2580 value_set_si(M
->p
[0][0], 1);
2581 value_decrement(M
->p
[0][P
->Dimension
+1],
2582 M
->p
[0][P
->Dimension
+1]);
2583 neg
= AddConstraints(M
->p
[0], 1, P
, MaxRays
);
2584 value_set_si(f
, -1);
2585 Vector_Scale(M
->p
[0]+1, M
->p
[0]+1, f
,
2587 value_decrement(M
->p
[0][P
->Dimension
+1],
2588 M
->p
[0][P
->Dimension
+1]);
2589 value_decrement(M
->p
[0][P
->Dimension
+1],
2590 M
->p
[0][P
->Dimension
+1]);
2591 pos
= AddConstraints(M
->p
[0], 1, P
, MaxRays
);
2592 if (emptyQ(neg
) && emptyQ(pos
)) {
2593 Polyhedron_Free(para
);
2594 Polyhedron_Free(pos
);
2595 Polyhedron_Free(neg
);
2599 fprintf(stderr
, "\nER: Order\n");
2600 #endif /* DEBUG_ER */
2601 EP
= barvinok_enumerate_e(para
, exist
, nparam
, MaxRays
);
2604 E
= barvinok_enumerate_e(pos
, exist
, nparam
, MaxRays
);
2606 free_evalue_refs(E
);
2610 E
= barvinok_enumerate_e(neg
, exist
, nparam
, MaxRays
);
2612 free_evalue_refs(E
);
2615 Polyhedron_Free(para
);
2616 Polyhedron_Free(pos
);
2617 Polyhedron_Free(neg
);
2622 } END_FORALL_PVertex_in_ParamPolyhedron
;
2625 } END_FORALL_PVertex_in_ParamPolyhedron
;
2628 /* Search for vertex coordinate to split on */
2629 /* First look for one independent of the parameters */
2630 FORALL_PVertex_in_ParamPolyhedron(V
, last
, PP
) {
2631 for (int i
= 0; i
< exist
; ++i
) {
2633 for (j
= 0; j
< nparam
; ++j
)
2634 if (value_notzero_p(V
->Vertex
->p
[i
][j
]))
2638 value_set_si(M
->p
[0][0], 1);
2639 Vector_Set(M
->p
[0]+1, 0, nvar
+exist
);
2640 Vector_Copy(V
->Vertex
->p
[i
],
2641 M
->p
[0] + 1 + nvar
+ exist
, nparam
+1);
2642 value_oppose(M
->p
[0][1+nvar
+i
],
2643 V
->Vertex
->p
[i
][nparam
+1]);
2645 Polyhedron
*pos
, *neg
;
2646 value_set_si(M
->p
[0][0], 1);
2647 value_decrement(M
->p
[0][P
->Dimension
+1],
2648 M
->p
[0][P
->Dimension
+1]);
2649 neg
= AddConstraints(M
->p
[0], 1, P
, MaxRays
);
2650 value_set_si(f
, -1);
2651 Vector_Scale(M
->p
[0]+1, M
->p
[0]+1, f
,
2653 value_decrement(M
->p
[0][P
->Dimension
+1],
2654 M
->p
[0][P
->Dimension
+1]);
2655 value_decrement(M
->p
[0][P
->Dimension
+1],
2656 M
->p
[0][P
->Dimension
+1]);
2657 pos
= AddConstraints(M
->p
[0], 1, P
, MaxRays
);
2658 if (emptyQ(neg
) || emptyQ(pos
)) {
2659 Polyhedron_Free(pos
);
2660 Polyhedron_Free(neg
);
2663 Polyhedron_Free(pos
);
2664 value_increment(M
->p
[0][P
->Dimension
+1],
2665 M
->p
[0][P
->Dimension
+1]);
2666 pos
= AddConstraints(M
->p
[0], 1, P
, MaxRays
);
2668 fprintf(stderr
, "\nER: Vertex\n");
2669 #endif /* DEBUG_ER */
2671 EP
= enumerate_or(pos
, exist
, nparam
, MaxRays
);
2676 } END_FORALL_PVertex_in_ParamPolyhedron
;
2680 /* Search for vertex coordinate to split on */
2681 /* Now look for one that depends on the parameters */
2682 FORALL_PVertex_in_ParamPolyhedron(V
, last
, PP
) {
2683 for (int i
= 0; i
< exist
; ++i
) {
2684 value_set_si(M
->p
[0][0], 1);
2685 Vector_Set(M
->p
[0]+1, 0, nvar
+exist
);
2686 Vector_Copy(V
->Vertex
->p
[i
],
2687 M
->p
[0] + 1 + nvar
+ exist
, nparam
+1);
2688 value_oppose(M
->p
[0][1+nvar
+i
],
2689 V
->Vertex
->p
[i
][nparam
+1]);
2691 Polyhedron
*pos
, *neg
;
2692 value_set_si(M
->p
[0][0], 1);
2693 value_decrement(M
->p
[0][P
->Dimension
+1],
2694 M
->p
[0][P
->Dimension
+1]);
2695 neg
= AddConstraints(M
->p
[0], 1, P
, MaxRays
);
2696 value_set_si(f
, -1);
2697 Vector_Scale(M
->p
[0]+1, M
->p
[0]+1, f
,
2699 value_decrement(M
->p
[0][P
->Dimension
+1],
2700 M
->p
[0][P
->Dimension
+1]);
2701 value_decrement(M
->p
[0][P
->Dimension
+1],
2702 M
->p
[0][P
->Dimension
+1]);
2703 pos
= AddConstraints(M
->p
[0], 1, P
, MaxRays
);
2704 if (emptyQ(neg
) || emptyQ(pos
)) {
2705 Polyhedron_Free(pos
);
2706 Polyhedron_Free(neg
);
2709 Polyhedron_Free(pos
);
2710 value_increment(M
->p
[0][P
->Dimension
+1],
2711 M
->p
[0][P
->Dimension
+1]);
2712 pos
= AddConstraints(M
->p
[0], 1, P
, MaxRays
);
2714 fprintf(stderr
, "\nER: ParamVertex\n");
2715 #endif /* DEBUG_ER */
2717 EP
= enumerate_or(pos
, exist
, nparam
, MaxRays
);
2722 } END_FORALL_PVertex_in_ParamPolyhedron
;
2730 Polyhedron_Free(CEq
);
2734 Param_Polyhedron_Free(PP
);
2741 evalue
*barvinok_enumerate_pip(Polyhedron
*P
,
2742 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
2747 evalue
*barvinok_enumerate_pip(Polyhedron
*P
,
2748 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
2750 int nvar
= P
->Dimension
- exist
- nparam
;
2751 evalue
*EP
= new_zero_ep();
2752 Polyhedron
*Q
, *N
, *T
= 0;
2758 fprintf(stderr
, "\nER: PIP\n");
2759 #endif /* DEBUG_ER */
2761 for (int i
= 0; i
< P
->Dimension
; ++i
) {
2764 bool posray
= false;
2765 bool negray
= false;
2766 value_set_si(min
, 0);
2767 for (int j
= 0; j
< P
->NbRays
; ++j
) {
2768 if (value_pos_p(P
->Ray
[j
][1+i
])) {
2770 if (value_zero_p(P
->Ray
[j
][1+P
->Dimension
]))
2772 } else if (value_neg_p(P
->Ray
[j
][1+i
])) {
2774 if (value_zero_p(P
->Ray
[j
][1+P
->Dimension
]))
2778 P
->Ray
[j
][1+i
], P
->Ray
[j
][1+P
->Dimension
]);
2779 if (value_lt(tmp
, min
))
2780 value_assign(min
, tmp
);
2785 assert(!(posray
&& negray
));
2786 assert(!negray
); // for now
2787 Polyhedron
*O
= T
? T
: P
;
2788 /* shift by a safe amount */
2789 Matrix
*M
= Matrix_Alloc(O
->NbRays
, O
->Dimension
+2);
2790 Vector_Copy(O
->Ray
[0], M
->p
[0], O
->NbRays
* (O
->Dimension
+2));
2791 for (int j
= 0; j
< P
->NbRays
; ++j
) {
2792 if (value_notzero_p(M
->p
[j
][1+P
->Dimension
])) {
2793 value_multiply(tmp
, min
, M
->p
[j
][1+P
->Dimension
]);
2794 value_substract(M
->p
[j
][1+i
], M
->p
[j
][1+i
], tmp
);
2799 T
= Rays2Polyhedron(M
, MaxRays
);
2802 /* negating a parameter requires that we substitute in the
2803 * sign again afterwards.
2806 assert(i
< nvar
+exist
);
2808 T
= Polyhedron_Copy(P
);
2809 for (int j
= 0; j
< T
->NbRays
; ++j
)
2810 value_oppose(T
->Ray
[j
][1+i
], T
->Ray
[j
][1+i
]);
2811 for (int j
= 0; j
< T
->NbConstraints
; ++j
)
2812 value_oppose(T
->Constraint
[j
][1+i
], T
->Constraint
[j
][1+i
]);
2818 Polyhedron
*D
= pip_lexmin(T
? T
: P
, exist
, nparam
);
2819 for (Q
= D
; Q
; Q
= N
) {
2823 exist
= Q
->Dimension
- nvar
- nparam
;
2824 E
= barvinok_enumerate_e(Q
, exist
, nparam
, MaxRays
);
2827 free_evalue_refs(E
);
2839 static bool is_single(Value
*row
, int pos
, int len
)
2841 return First_Non_Zero(row
, pos
) == -1 &&
2842 First_Non_Zero(row
+pos
+1, len
-pos
-1) == -1;
2845 static evalue
* barvinok_enumerate_e_r(Polyhedron
*P
,
2846 unsigned exist
, unsigned nparam
, unsigned MaxRays
);
2849 static int er_level
= 0;
2851 evalue
* barvinok_enumerate_e(Polyhedron
*P
,
2852 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
2854 fprintf(stderr
, "\nER: level %i\n", er_level
);
2855 int nvar
= P
->Dimension
- exist
- nparam
;
2856 fprintf(stderr
, "%d %d %d\n", nvar
, exist
, nparam
);
2858 Polyhedron_Print(stderr
, P_VALUE_FMT
, P
);
2860 P
= DomainConstraintSimplify(Polyhedron_Copy(P
), MaxRays
);
2861 evalue
*EP
= barvinok_enumerate_e_r(P
, exist
, nparam
, MaxRays
);
2867 evalue
* barvinok_enumerate_e(Polyhedron
*P
,
2868 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
2870 P
= DomainConstraintSimplify(Polyhedron_Copy(P
), MaxRays
);
2871 evalue
*EP
= barvinok_enumerate_e_r(P
, exist
, nparam
, MaxRays
);
2877 static evalue
* barvinok_enumerate_e_r(Polyhedron
*P
,
2878 unsigned exist
, unsigned nparam
, unsigned MaxRays
)
2881 Polyhedron
*U
= Universe_Polyhedron(nparam
);
2882 evalue
*EP
= barvinok_enumerate_ev(P
, U
, MaxRays
);
2883 //char *param_name[] = {"P", "Q", "R", "S", "T" };
2884 //print_evalue(stdout, EP, param_name);
2889 int nvar
= P
->Dimension
- exist
- nparam
;
2890 int len
= P
->Dimension
+ 2;
2893 return new_zero_ep();
2895 if (nvar
== 0 && nparam
== 0) {
2896 evalue
*EP
= new_zero_ep();
2897 barvinok_count(P
, &EP
->x
.n
, MaxRays
);
2898 if (value_pos_p(EP
->x
.n
))
2899 value_set_si(EP
->x
.n
, 1);
2904 for (r
= 0; r
< P
->NbRays
; ++r
)
2905 if (value_zero_p(P
->Ray
[r
][0]) ||
2906 value_zero_p(P
->Ray
[r
][P
->Dimension
+1])) {
2908 for (i
= 0; i
< nvar
; ++i
)
2909 if (value_notzero_p(P
->Ray
[r
][i
+1]))
2913 for (i
= nvar
+ exist
; i
< nvar
+ exist
+ nparam
; ++i
)
2914 if (value_notzero_p(P
->Ray
[r
][i
+1]))
2916 if (i
>= nvar
+ exist
+ nparam
)
2919 if (r
< P
->NbRays
) {
2920 evalue
*EP
= new_zero_ep();
2921 value_set_si(EP
->x
.n
, -1);
2926 for (r
= 0; r
< P
->NbEq
; ++r
)
2927 if ((first
= First_Non_Zero(P
->Constraint
[r
]+1+nvar
, exist
)) != -1)
2930 if (First_Non_Zero(P
->Constraint
[r
]+1+nvar
+first
+1,
2931 exist
-first
-1) != -1) {
2932 Polyhedron
*T
= rotate_along(P
, r
, nvar
, exist
, MaxRays
);
2934 fprintf(stderr
, "\nER: Equality\n");
2935 #endif /* DEBUG_ER */
2936 evalue
*EP
= barvinok_enumerate_e(T
, exist
-1, nparam
, MaxRays
);
2941 fprintf(stderr
, "\nER: Fixed\n");
2942 #endif /* DEBUG_ER */
2944 return barvinok_enumerate_e(P
, exist
-1, nparam
, MaxRays
);
2946 Polyhedron
*T
= Polyhedron_Copy(P
);
2947 SwapColumns(T
, nvar
+1, nvar
+1+first
);
2948 evalue
*EP
= barvinok_enumerate_e(T
, exist
-1, nparam
, MaxRays
);
2955 Vector
*row
= Vector_Alloc(len
);
2956 value_set_si(row
->p
[0], 1);
2961 enum constraint
* info
= new constraint
[exist
];
2962 for (int i
= 0; i
< exist
; ++i
) {
2964 for (int l
= P
->NbEq
; l
< P
->NbConstraints
; ++l
) {
2965 if (value_negz_p(P
->Constraint
[l
][nvar
+i
+1]))
2967 bool l_parallel
= is_single(P
->Constraint
[l
]+nvar
+1, i
, exist
);
2968 for (int u
= P
->NbEq
; u
< P
->NbConstraints
; ++u
) {
2969 if (value_posz_p(P
->Constraint
[u
][nvar
+i
+1]))
2971 bool lu_parallel
= l_parallel
||
2972 is_single(P
->Constraint
[u
]+nvar
+1, i
, exist
);
2973 value_oppose(f
, P
->Constraint
[u
][nvar
+i
+1]);
2974 Vector_Combine(P
->Constraint
[l
]+1, P
->Constraint
[u
]+1, row
->p
+1,
2975 f
, P
->Constraint
[l
][nvar
+i
+1], len
-1);
2976 if (!(info
[i
] & INDEPENDENT
)) {
2978 for (j
= 0; j
< exist
; ++j
)
2979 if (j
!= i
&& value_notzero_p(row
->p
[nvar
+j
+1]))
2982 //printf("independent: i: %d, l: %d, u: %d\n", i, l, u);
2983 info
[i
] = (constraint
)(info
[i
] | INDEPENDENT
);
2986 if (info
[i
] & ALL_POS
) {
2987 value_addto(row
->p
[len
-1], row
->p
[len
-1],
2988 P
->Constraint
[l
][nvar
+i
+1]);
2989 value_addto(row
->p
[len
-1], row
->p
[len
-1], f
);
2990 value_multiply(f
, f
, P
->Constraint
[l
][nvar
+i
+1]);
2991 value_substract(row
->p
[len
-1], row
->p
[len
-1], f
);
2992 value_decrement(row
->p
[len
-1], row
->p
[len
-1]);
2993 ConstraintSimplify(row
->p
, row
->p
, len
, &f
);
2994 value_set_si(f
, -1);
2995 Vector_Scale(row
->p
+1, row
->p
+1, f
, len
-1);
2996 value_decrement(row
->p
[len
-1], row
->p
[len
-1]);
2997 Polyhedron
*T
= AddConstraints(row
->p
, 1, P
, MaxRays
);
2999 //printf("not all_pos: i: %d, l: %d, u: %d\n", i, l, u);
3000 info
[i
] = (constraint
)(info
[i
] ^ ALL_POS
);
3002 //puts("pos remainder");
3003 //Polyhedron_Print(stdout, P_VALUE_FMT, T);
3006 if (!(info
[i
] & ONE_NEG
)) {
3008 negative_test_constraint(P
->Constraint
[l
],
3010 row
->p
, nvar
+i
, len
, &f
);
3011 oppose_constraint(row
->p
, len
, &f
);
3012 Polyhedron
*T
= AddConstraints(row
->p
, 1, P
, MaxRays
);
3014 //printf("one_neg i: %d, l: %d, u: %d\n", i, l, u);
3015 info
[i
] = (constraint
)(info
[i
] | ONE_NEG
);
3017 //puts("neg remainder");
3018 //Polyhedron_Print(stdout, P_VALUE_FMT, T);
3022 if (!(info
[i
] & ALL_POS
) && (info
[i
] & ONE_NEG
))
3026 if (info
[i
] & ALL_POS
)
3033 for (int i = 0; i < exist; ++i)
3034 printf("%i: %i\n", i, info[i]);
3036 for (int i
= 0; i
< exist
; ++i
)
3037 if (info
[i
] & ALL_POS
) {
3039 fprintf(stderr
, "\nER: Positive\n");
3040 #endif /* DEBUG_ER */
3042 // Maybe we should chew off some of the fat here
3043 Matrix
*M
= Matrix_Alloc(P
->Dimension
, P
->Dimension
+1);
3044 for (int j
= 0; j
< P
->Dimension
; ++j
)
3045 value_set_si(M
->p
[j
][j
+ (j
>= i
+nvar
)], 1);
3046 Polyhedron
*T
= Polyhedron_Image(P
, M
, MaxRays
);
3048 evalue
*EP
= barvinok_enumerate_e(T
, exist
-1, nparam
, MaxRays
);
3055 for (int i
= 0; i
< exist
; ++i
)
3056 if (info
[i
] & ONE_NEG
) {
3058 fprintf(stderr
, "\nER: Negative\n");
3059 #endif /* DEBUG_ER */
3064 return barvinok_enumerate_e(P
, exist
-1, nparam
, MaxRays
);
3066 Polyhedron
*T
= Polyhedron_Copy(P
);
3067 SwapColumns(T
, nvar
+1, nvar
+1+i
);
3068 evalue
*EP
= barvinok_enumerate_e(T
, exist
-1, nparam
, MaxRays
);
3073 for (int i
= 0; i
< exist
; ++i
)
3074 if (info
[i
] & INDEPENDENT
) {
3075 Polyhedron
*pos
, *neg
;
3077 /* Find constraint again and split off negative part */
3079 if (SplitOnVar(P
, i
, nvar
, len
, exist
, MaxRays
,
3080 row
, f
, true, &pos
, &neg
)) {
3082 fprintf(stderr
, "\nER: Split\n");
3083 #endif /* DEBUG_ER */
3086 barvinok_enumerate_e(neg
, exist
-1, nparam
, MaxRays
);
3088 barvinok_enumerate_e(pos
, exist
, nparam
, MaxRays
);
3090 free_evalue_refs(E
);
3092 Polyhedron_Free(neg
);
3093 Polyhedron_Free(pos
);
3107 EP
= enumerate_line(P
, exist
, nparam
, MaxRays
);
3111 EP
= barvinok_enumerate_pip(P
, exist
, nparam
, MaxRays
);
3115 EP
= enumerate_redundant_ray(P
, exist
, nparam
, MaxRays
);
3119 EP
= enumerate_sure(P
, exist
, nparam
, MaxRays
);
3123 EP
= enumerate_ray(P
, exist
, nparam
, MaxRays
);
3127 EP
= enumerate_sure2(P
, exist
, nparam
, MaxRays
);
3131 F
= unfringe(P
, MaxRays
);
3132 if (!PolyhedronIncludes(F
, P
)) {
3134 fprintf(stderr
, "\nER: Fringed\n");
3135 #endif /* DEBUG_ER */
3136 EP
= barvinok_enumerate_e(F
, exist
, nparam
, MaxRays
);
3143 EP
= enumerate_vd(&P
, exist
, nparam
, MaxRays
);
3148 EP
= enumerate_sum(P
, exist
, nparam
, MaxRays
);
3155 Polyhedron
*pos
, *neg
;
3156 for (i
= 0; i
< exist
; ++i
)
3157 if (SplitOnVar(P
, i
, nvar
, len
, exist
, MaxRays
,
3158 row
, f
, false, &pos
, &neg
))
3164 EP
= enumerate_or(pos
, exist
, nparam
, MaxRays
);