remove_parameter_equalities: stop "removing" equalities if polyhedron is empty
[barvinok.git] / euler.cc
blob9e7c2eb946ce59fddae1d395ceaf42e5f461ce33
1 /*
2 * Sum polynomial over integer points in polytope using local
3 * Euler-Maclaurin formula by Berline and Vergne.
4 */
6 #include <assert.h>
7 #include <barvinok/options.h>
8 #include <barvinok/util.h>
9 #include "bernoulli.h"
10 #include "conversion.h"
11 #include "decomposer.h"
12 #include "euler.h"
13 #include "lattice_point.h"
14 #include "param_util.h"
15 #include "reduce_domain.h"
17 /* Compute total degree in first nvar variables */
18 static unsigned total_degree(const evalue *e, unsigned nvar)
20 int i;
21 unsigned max_degree;
23 if (value_notzero_p(e->d))
24 return 0;
25 assert(e->x.p->type == polynomial);
26 if (e->x.p->pos-1 >= nvar)
27 return 0;
29 max_degree = total_degree(&e->x.p->arr[0], nvar);
30 for (i = 1; i < e->x.p->size; ++i) {
31 unsigned degree = i + total_degree(&e->x.p->arr[i], nvar);
32 if (degree > max_degree)
33 max_degree = degree;
35 return max_degree;
38 struct fact {
39 Value *fact;
40 unsigned size;
41 unsigned n;
43 struct fact fact;
45 Value *factorial(unsigned n)
47 if (n < fact.n)
48 return &fact.fact[n];
50 if (n >= fact.size) {
51 int size = 3*(n + 5)/2;
53 fact.fact = (Value *)realloc(fact.fact, size*sizeof(Value));
54 fact.size = size;
56 for (int i = fact.n; i <= n; ++i) {
57 value_init(fact.fact[i]);
58 if (!i)
59 value_set_si(fact.fact[0], 1);
60 else
61 mpz_mul_ui(fact.fact[i], fact.fact[i-1], i);
63 fact.n = n+1;
64 return &fact.fact[n];
67 struct binom {
68 Vector **binom;
69 unsigned size;
70 unsigned n;
72 struct binom binom;
74 Value *binomial(unsigned n, unsigned k)
76 if (n < binom.n)
77 return &binom.binom[n]->p[k];
79 if (n >= binom.size) {
80 int size = 3*(n + 5)/2;
82 binom.binom = (Vector **)realloc(binom.binom, size*sizeof(Vector *));
83 binom.size = size;
85 for (int i = binom.n; i <= n; ++i) {
86 binom.binom[i] = Vector_Alloc(i+1);
87 if (!i)
88 value_set_si(binom.binom[0]->p[0], 1);
89 else {
90 value_set_si(binom.binom[i]->p[0], 1);
91 value_set_si(binom.binom[i]->p[i], 1);
92 for (int j = 1; j < i; ++j)
93 value_addto(binom.binom[i]->p[j],
94 binom.binom[i-1]->p[j-1], binom.binom[i-1]->p[j]);
97 binom.n = n+1;
98 return &binom.binom[n]->p[k];
101 struct power {
102 int n;
103 Vector *powers;
105 power(Value v, int max) {
106 powers = Vector_Alloc(max+1);
107 assert(powers);
108 value_set_si(powers->p[0], 1);
109 if (max >= 1)
110 value_assign(powers->p[1], v);
111 n = 2;
113 ~power() {
114 Vector_Free(powers);
116 Value *operator[](int exp) {
117 assert(exp >= 0);
118 assert(exp < powers->Size);
119 for (; n <= exp; ++n)
120 value_multiply(powers->p[n], powers->p[n-1], powers->p[1]);
121 return &powers->p[exp];
126 * Computes the coefficients of
128 * \mu(-t + R_+)(\xi) = \sum_{n=0)^\infty -b(n+1, t)/(n+1)! \xi^n
130 * where b(n, t) is the Bernoulli polynomial of degree n in t
131 * and t(p) is an expression (a fractional part) of the parameters p
132 * such that 0 <= t(p) < 1 for all values of the parameters.
133 * The coefficients are computed on demand up to (and including)
134 * the maximal degree max_degree.
136 struct mu_1d {
137 unsigned max_degree;
138 evalue **coefficients;
139 const evalue *t;
141 mu_1d(unsigned max_degree, evalue *t) : max_degree(max_degree), t(t) {
142 coefficients = new evalue *[max_degree+1];
143 for (int i = 0; i < max_degree+1; ++i)
144 coefficients[i] = NULL;
146 ~mu_1d() {
147 for (int i = 0; i < max_degree+1; ++i)
148 if (coefficients[i])
149 evalue_free(coefficients[i]);
150 delete [] coefficients;
152 void compute_coefficient(unsigned n);
153 const evalue *coefficient(unsigned n) {
154 if (!coefficients[n])
155 compute_coefficient(n);
156 return coefficients[n];
160 void mu_1d::compute_coefficient(unsigned n)
162 struct poly_list *bernoulli = bernoulli_compute(n+1);
163 evalue *c = evalue_polynomial(bernoulli->poly[n+1], t);
164 evalue_negate(c);
165 evalue_div(c, *factorial(n+1));
167 coefficients[n] = c;
171 * Computes the coefficients of
173 * \mu(a)(y) = \sum_{n_1} \sum_{n_2} c_{n_1,n_2} y^{n_1} y^{n_2}
175 * with c_{n1,n2} given by
177 * b(n1+1,t1)/(n1+1)! b(n2+1,t2)/(n2+1)!
178 * - b(n1+n2+2,t2)/(n1+n2+2)! (-c1)^{n1+1}
179 * - b(n1+n2+2,t1)/(n1+n2+2)! (-c2)^{n2+1}
181 * where b(n, t) is the Bernoulli polynomial of degree n in t,
182 * t1(p) and t2(p) are expressions (fractional parts) of the parameters p
183 * such that 0 <= t1(p), t2(p) < 1 for all values of the parameters
184 * and c1 = cn/c1d and c2 = cn/c2d.
185 * The coefficients are computed on demand up to (and including)
186 * the maximal degree (n1,n2) = (max_degree,max_degree).
188 * bernoulli_t[i][j] contains b(j+1, t_i)/(j+1)!
190 struct mu_2d {
191 unsigned max_degree;
192 evalue ***coefficients;
193 /* bernoulli_t[i][n] stores b(n+1, t_i)/(n+1)! */
194 evalue **bernoulli_t[2];
195 /* stores the powers of -cn */
196 struct power *power_cn;
197 struct power *power_c1d;
198 struct power *power_c2d;
199 const evalue *t[2];
201 mu_2d(unsigned max_degree, evalue *t1, evalue *t2,
202 Value cn, Value c1d, Value c2d) : max_degree(max_degree) {
203 t[0] = t1;
204 t[1] = t2;
205 coefficients = new evalue **[max_degree+1];
206 for (int i = 0; i < max_degree+1; ++i) {
207 coefficients[i] = new evalue *[max_degree+1];
208 for (int j = 0; j < max_degree+1; ++j)
209 coefficients[i][j] = NULL;
211 for (int i = 0; i < 2; ++i) {
212 bernoulli_t[i] = new evalue *[max_degree+2];
213 for (int j = 0; j < max_degree+2; ++j)
214 bernoulli_t[i][j] = NULL;
216 value_oppose(cn, cn);
217 power_cn = new struct power(cn, max_degree+1);
218 value_oppose(cn, cn);
219 power_c1d = new struct power(c1d, max_degree+1);
220 power_c2d = new struct power(c2d, max_degree+1);
222 ~mu_2d() {
223 for (int i = 0; i < max_degree+1; ++i) {
224 for (int j = 0; j < max_degree+1; ++j)
225 if (coefficients[i][j])
226 evalue_free(coefficients[i][j]);
227 delete [] coefficients[i];
229 delete [] coefficients;
230 for (int i = 0; i < 2; ++i)
231 for (int j = 0; j < max_degree+2; ++j)
232 if (bernoulli_t[i][j])
233 evalue_free(bernoulli_t[i][j]);
234 for (int i = 0; i < 2; ++i)
235 delete [] bernoulli_t[i];
236 delete power_cn;
237 delete power_c1d;
238 delete power_c2d;
240 const evalue *get_bernoulli(int n, int i);
242 void compute_coefficient(unsigned n1, unsigned n2);
243 const evalue *coefficient(unsigned n1, unsigned n2) {
244 if (!coefficients[n1][n2])
245 compute_coefficient(n1, n2);
246 return coefficients[n1][n2];
251 * Returns b(n, t_i)/n!
253 const evalue *mu_2d::get_bernoulli(int n, int i)
255 if (!bernoulli_t[i][n-1]) {
256 struct poly_list *bernoulli = bernoulli_compute(n);
257 bernoulli_t[i][n-1] = evalue_polynomial(bernoulli->poly[n], t[i]);
258 evalue_div(bernoulli_t[i][n-1], *factorial(n));
260 return bernoulli_t[i][n-1];
263 void mu_2d::compute_coefficient(unsigned n1, unsigned n2)
265 evalue *c = evalue_dup(get_bernoulli(n1+1, 0));
266 emul(get_bernoulli(n2+1, 1), c);
268 if (value_notzero_p(*(*power_cn)[1])) {
269 evalue *t;
270 Value neg_power;
272 value_init(neg_power);
274 t = evalue_dup(get_bernoulli(n1+n2+2, 1));
275 value_multiply(neg_power,
276 *(*power_cn)[n1+1], *binomial(n1+n2+1, n1+1));
277 value_oppose(neg_power, neg_power);
278 evalue_mul_div(t, neg_power, *(*power_c1d)[n1+1]);
279 eadd(t, c);
280 evalue_free(t);
282 t = evalue_dup(get_bernoulli(n1+n2+2, 0));
283 value_multiply(neg_power,
284 *(*power_cn)[n2+1], *binomial(n1+n2+1, n2+1));
285 value_oppose(neg_power, neg_power);
286 evalue_mul_div(t, neg_power, *(*power_c2d)[n2+1]);
287 eadd(t, c);
288 evalue_free(t);
290 value_clear(neg_power);
293 coefficients[n1][n2] = c;
296 /* Later: avoid recomputation of mu of faces that appear in
297 * more than one domain.
299 struct summator_2d : public signed_cone_consumer, public vertex_decomposer {
300 const evalue *polynomial;
301 Value *inner;
302 unsigned nparam;
304 /* substitutions to use when result is 0-dimensional (only parameters) */
305 evalue **subs_0d;
306 /* substitutions to use when result is 1-dimensional */
307 evalue **subs_1d;
308 evalue *sum;
310 summator_2d(evalue *e, Param_Polyhedron *PP, Value *inner,
311 unsigned nparam) :
312 polynomial(e), vertex_decomposer(PP, *this),
313 inner(inner), nparam(nparam) {
314 sum = evalue_zero();
316 subs_0d = new evalue *[2+nparam];
317 subs_1d = new evalue *[2+nparam];
318 subs_0d[0] = NULL;
319 subs_0d[1] = NULL;
320 subs_1d[0] = NULL;
321 subs_1d[1] = NULL;
322 for (int i = 0; i < nparam; ++i) {
323 subs_0d[2+i] = evalue_var(i);
324 subs_1d[2+i] = evalue_var(1+i);
327 ~summator_2d() {
328 for (int i = 0; i < nparam; ++i) {
329 evalue_free(subs_0d[2+i]);
330 evalue_free(subs_1d[2+i]);
332 delete [] subs_0d;
333 delete [] subs_1d;
335 evalue *summate_over_pdomain(Polyhedron *P,
336 Param_Domain *PD,
337 struct barvinok_options *options);
338 void handle_facet(Param_Polyhedron *PP, Param_Domain *FD, Value *normal);
339 void integrate(Polyhedron *P, Param_Domain *PD);
340 virtual void handle(const signed_cone& sc, barvinok_options *options);
343 /* Replaces poly by its derivative along variable var */
344 static void evalue_derive(evalue *poly, int var)
346 if (value_notzero_p(poly->d)) {
347 value_set_si(poly->x.n, 0);
348 value_set_si(poly->d, 1);
349 return;
351 assert(poly->x.p->type == polynomial);
352 if (poly->x.p->pos-1 > var) {
353 free_evalue_refs(poly);
354 value_init(poly->d);
355 evalue_set_si(poly, 0, 1);
356 return;
359 if (poly->x.p->pos-1 < var) {
360 for (int i = 0; i < poly->x.p->size; ++i)
361 evalue_derive(&poly->x.p->arr[i], var);
362 reduce_evalue(poly);
363 return;
366 assert(poly->x.p->size > 1);
367 enode *p = poly->x.p;
368 free_evalue_refs(&p->arr[0]);
369 if (p->size == 2) {
370 value_clear(poly->d);
371 *poly = p->arr[1];
372 free(p);
373 return;
376 p->size--;
377 Value factor;
378 value_init(factor);
379 for (int i = 0; i < p->size; ++i) {
380 value_set_si(factor, i+1);
381 p->arr[i] = p->arr[i+1];
382 evalue_mul(&p->arr[i], factor);
384 value_clear(factor);
387 /* Check whether e is constant with respect to variable var */
388 static int evalue_is_constant(evalue *e, int var)
390 int i;
392 if (value_notzero_p(e->d))
393 return 1;
394 if (e->x.p->type == polynomial && e->x.p->pos-1 == var)
395 return 0;
396 assert(e->x.p->type == polynomial ||
397 e->x.p->type == fractional ||
398 e->x.p->type == flooring);
399 for (i = 0; i < e->x.p->size; ++i)
400 if (!evalue_is_constant(&e->x.p->arr[i], var))
401 return 0;
402 return 1;
405 /* Replaces poly by its anti-derivative with constant 0 along variable var */
406 static void evalue_anti_derive(evalue *poly, int var)
408 enode *p;
410 if (value_zero_p(poly->d) &&
411 poly->x.p->type == polynomial && poly->x.p->pos-1 < var) {
412 for (int i = 0; i < poly->x.p->size; ++i)
413 evalue_anti_derive(&poly->x.p->arr[i], var);
414 reduce_evalue(poly);
415 return;
418 if (evalue_is_constant(poly, var)) {
419 p = new_enode(polynomial, 2, 1+var);
420 evalue_set_si(&p->arr[0], 0, 1);
421 value_clear(p->arr[1].d);
422 p->arr[1] = *poly;
423 value_init(poly->d);
424 poly->x.p = p;
425 return;
427 assert(poly->x.p->type == polynomial);
429 p = new_enode(polynomial, poly->x.p->size+1, 1+var);
430 evalue_set_si(&p->arr[0], 0, 1);
431 for (int i = 0; i < poly->x.p->size; ++i) {
432 value_clear(p->arr[1+i].d);
433 p->arr[1+i] = poly->x.p->arr[i];
434 value_set_si(poly->d, 1+i);
435 evalue_div(&p->arr[1+i], poly->d);
437 free(poly->x.p);
438 poly->x.p = p;
439 value_set_si(poly->d, 0);
442 /* Computes offsets in the basis given by the rays of the cone
443 * to the integer point in the fundamental parallelepiped of
444 * the cone.
445 * The resulting evalues contain only the parameters as variables.
447 evalue **offsets_to_integer_point(Matrix *Rays, Matrix *vertex)
449 unsigned nparam = vertex->NbColumns-2;
450 evalue **t = new evalue *[2];
452 if (value_one_p(vertex->p[0][nparam+1])) {
453 t[0] = evalue_zero();
454 t[1] = evalue_zero();
455 } else {
456 Matrix *R2 = Matrix_Copy(Rays);
457 Matrix_Transposition(R2);
458 Matrix *inv = Matrix_Alloc(Rays->NbColumns, Rays->NbRows);
459 int ok = Matrix_Inverse(R2, inv);
460 assert(ok);
461 Matrix_Free(R2);
463 /* We want the fractional part of the negative relative coordinates */
464 Vector_Oppose(inv->p[0], inv->p[0], inv->NbColumns);
465 Vector_Oppose(inv->p[1], inv->p[1], inv->NbColumns);
467 Matrix *neg_rel = Matrix_Alloc(2, nparam+1);
468 vertex->NbColumns--;
469 Matrix_Product(inv, vertex, neg_rel);
470 Matrix_Free(inv);
471 vertex->NbColumns++;
472 t[0] = fractional_part(neg_rel->p[0], vertex->p[0][nparam+1],
473 nparam, NULL);
474 t[1] = fractional_part(neg_rel->p[1], vertex->p[0][nparam+1],
475 nparam, NULL);
476 Matrix_Free(neg_rel);
479 return t;
483 * Called from decompose_at_vertex.
485 * Handles a cone in the signed decomposition of the supporting
486 * cone of a vertex. The cone is assumed to be unimodular.
488 void summator_2d::handle(const signed_cone& sc, barvinok_options *options)
490 evalue **t;
491 unsigned degree = total_degree(polynomial, 2);
493 subs_0d[0] = affine2evalue(V->Vertex->p[0],
494 V->Vertex->p[0][nparam+1], nparam);
495 subs_0d[1] = affine2evalue(V->Vertex->p[1],
496 V->Vertex->p[1][nparam+1], nparam);
498 assert(sc.det == 1);
500 assert(V->Vertex->NbRows > 0);
501 Param_Vertex_Common_Denominator(V);
503 Matrix *Rays = zz2matrix(sc.rays);
505 t = offsets_to_integer_point(Rays, V->Vertex);
507 Vector *c = Vector_Alloc(3);
508 Inner_Product(Rays->p[0], Rays->p[1], 2, &c->p[0]);
509 Inner_Product(Rays->p[0], Rays->p[0], 2, &c->p[1]);
510 Inner_Product(Rays->p[1], Rays->p[1], 2, &c->p[2]);
512 mu_2d mu(degree, t[0], t[1], c->p[0], c->p[1], c->p[2]);
513 Vector_Free(c);
515 struct power power_r00(Rays->p[0][0], degree);
516 struct power power_r01(Rays->p[0][1], degree);
517 struct power power_r10(Rays->p[1][0], degree);
518 struct power power_r11(Rays->p[1][1], degree);
520 Value factor, tmp1, tmp2;
521 value_init(factor);
522 value_init(tmp1);
523 value_init(tmp2);
524 evalue *res = evalue_zero();
525 evalue *dx1 = evalue_dup(polynomial);
526 for (int i = 0; !EVALUE_IS_ZERO(*dx1); ++i) {
527 evalue *dx2 = evalue_dup(dx1);
528 for (int j = 0; !EVALUE_IS_ZERO(*dx2); ++j) {
529 evalue *cij = evalue_zero();
530 for (int n1 = 0; n1 <= i+j; ++n1) {
531 int n2 = i+j-n1;
532 value_set_si(factor, 0);
533 for (int k = max(0, i-n2); k <= i && k <= n1; ++k) {
534 int l = i-k;
535 value_multiply(tmp1, *power_r00[k], *power_r01[n1-k]);
536 value_multiply(tmp1, tmp1, *power_r10[l]);
537 value_multiply(tmp1, tmp1, *power_r11[n2-l]);
538 value_multiply(tmp1, tmp1, *binomial(n1, k));
539 value_multiply(tmp1, tmp1, *binomial(n2, l));
540 value_addto(factor, factor, tmp1);
542 if (value_zero_p(factor))
543 continue;
545 evalue *c = evalue_dup(mu.coefficient(n1, n2));
546 evalue_mul(c, factor);
547 eadd(c, cij);
548 evalue_free(c);
550 evalue *d = evalue_dup(dx2);
551 evalue_substitute(d, subs_0d);
552 emul(cij, d);
553 evalue_free(cij);
554 eadd(d, res);
555 evalue_free(d);
556 evalue_derive(dx2, 1);
558 evalue_free(dx2);
559 evalue_derive(dx1, 0);
561 evalue_free(dx1);
562 for (int i = 0; i < 2; ++i) {
563 evalue_free(subs_0d[i]);
564 subs_0d[i] = NULL;
565 evalue_free(t[i]);
567 delete [] t;
568 value_clear(factor);
569 value_clear(tmp1);
570 value_clear(tmp2);
571 Matrix_Free(Rays);
573 if (sc.sign < 0)
574 evalue_negate(res);
575 eadd(res, sum);
576 evalue_free(res);
579 evalue *summator_2d::summate_over_pdomain(Polyhedron *P,
580 Param_Domain *PD,
581 struct barvinok_options *options)
583 int j;
584 Param_Vertices *V;
586 assert(PP->V->Vertex->NbRows == 2);
588 FORALL_PVertex_in_ParamPolyhedron(V, PD, PP) // _i internal counter
589 decompose_at_vertex(V, _i, options);
590 END_FORALL_PVertex_in_ParamPolyhedron;
592 Vector *normal = Vector_Alloc(2);
593 for (j = P->NbEq; j < P->NbConstraints; ++j) {
594 Param_Domain *FD;
595 Vector_Copy(P->Constraint[j]+1, normal->p, 2);
596 if (value_zero_p(normal->p[0]) && value_zero_p(normal->p[1]))
597 continue;
599 FD = Param_Polyhedron_Facet(PP, PD, P, j);
600 Vector_Normalize(normal->p, 2);
601 handle_facet(PP, FD, normal->p);
602 Param_Domain_Free(FD);
604 Vector_Free(normal);
606 integrate(P, PD);
608 return sum;
611 void summator_2d::handle_facet(Param_Polyhedron *PP, Param_Domain *FD,
612 Value *normal)
614 Param_Vertices *V;
615 int nbV = 0;
616 Param_Vertices *vertex[2];
617 Value det;
618 unsigned degree = total_degree(polynomial, 2);
620 FORALL_PVertex_in_ParamPolyhedron(V, FD, PP)
621 vertex[nbV++] = V;
622 END_FORALL_PVertex_in_ParamPolyhedron;
624 assert(nbV == 2);
626 /* We can take either vertex[0] or vertex[1];
627 * the result should be the same
629 Param_Vertex_Common_Denominator(vertex[0]);
631 /* The extra variable in front is the coordinate along the facet. */
632 Vector *coef_normal = Vector_Alloc(1 + nparam + 2);
633 Vector_Combine(vertex[0]->Vertex->p[0], vertex[0]->Vertex->p[1],
634 coef_normal->p+1, normal[0], normal[1], nparam+1);
635 value_assign(coef_normal->p[1+nparam+1], vertex[0]->Vertex->p[0][nparam+1]);
636 Vector_Normalize(coef_normal->p, coef_normal->Size);
638 Vector *base = Vector_Alloc(2);
639 value_oppose(base->p[0], normal[1]);
640 value_assign(base->p[1], normal[0]);
642 value_init(det);
643 Inner_Product(normal, normal, 2, &det);
645 Vector *s = Vector_Alloc(1+nparam+2);
646 value_multiply(s->p[1+nparam+1], coef_normal->p[1+nparam+1], det);
647 value_multiply(s->p[0], base->p[0], s->p[1+nparam+1]);
648 Vector_Scale(coef_normal->p+1, s->p+1, normal[0], nparam+1);
649 subs_1d[0] = affine2evalue(s->p, s->p[1+nparam+1], 1+nparam);
650 value_multiply(s->p[0], base->p[1], s->p[1+nparam+1]);
651 Vector_Scale(coef_normal->p+1, s->p+1, normal[1], nparam+1);
652 subs_1d[1] = affine2evalue(s->p, s->p[1+nparam+1], 1+nparam);
653 Vector_Free(s);
655 evalue *t;
656 if (value_one_p(coef_normal->p[coef_normal->Size-1]))
657 t = evalue_zero();
658 else {
659 Vector_Oppose(coef_normal->p+1, coef_normal->p+1, nparam+1);
660 t = fractional_part(coef_normal->p,
661 coef_normal->p[coef_normal->Size-1],
662 1+nparam, NULL);
664 Vector_Free(coef_normal);
666 mu_1d mu(degree, t);
668 struct power power_normal0(normal[0], degree);
669 struct power power_normal1(normal[1], degree);
670 struct power power_det(det, degree);
672 Value(factor);
673 value_init(factor);
674 evalue *res = evalue_zero();
675 evalue *dx1 = evalue_dup(polynomial);
676 for (int i = 0; !EVALUE_IS_ZERO(*dx1); ++i) {
677 evalue *dx2 = evalue_dup(dx1);
678 for (int j = 0; !EVALUE_IS_ZERO(*dx2); ++j) {
679 value_multiply(factor, *power_normal0[i], *power_normal1[j]);
680 if (value_notzero_p(factor)) {
681 value_multiply(factor, factor, *binomial(i+j, i));
683 evalue *c = evalue_dup(mu.coefficient(i+j));
684 evalue_mul_div(c, factor, *power_det[i+j]);
686 evalue *d = evalue_dup(dx2);
687 evalue_substitute(d, subs_1d);
688 emul(c, d);
689 evalue_free(c);
690 eadd(d, res);
691 evalue_free(d);
693 evalue_derive(dx2, 1);
695 evalue_free(dx2);
696 evalue_derive(dx1, 0);
698 evalue_free(dx1);
699 evalue_free(t);
700 for (int i = 0; i < 2; ++i) {
701 evalue_free(subs_1d[i]);
702 subs_1d[i] = NULL;
704 value_clear(factor);
705 evalue_anti_derive(res, 0);
707 Matrix *z = Matrix_Alloc(2, nparam+2);
708 Vector *fixed_z = Vector_Alloc(2);
709 for (int i = 0; i < 2; ++i) {
710 Vector_Combine(vertex[i]->Vertex->p[0], vertex[i]->Vertex->p[1],
711 z->p[i], base->p[0], base->p[1], nparam+1);
712 value_multiply(z->p[i][nparam+1],
713 det, vertex[i]->Vertex->p[0][nparam+1]);
714 Inner_Product(z->p[i], inner, nparam+1, &fixed_z->p[i]);
716 Vector_Free(base);
717 /* Put on a common denominator */
718 value_multiply(fixed_z->p[0], fixed_z->p[0], z->p[1][nparam+1]);
719 value_multiply(fixed_z->p[1], fixed_z->p[1], z->p[0][nparam+1]);
720 /* Make sure z->p[0] is smaller than z->p[1] (for an internal
721 * point of the chamber and hence for all parameter values in
722 * the chamber), to ensure we integrate in the right direction.
724 if (value_lt(fixed_z->p[1], fixed_z->p[0]))
725 Vector_Exchange(z->p[0], z->p[1], nparam+2);
726 Vector_Free(fixed_z);
727 value_clear(det);
729 subs_0d[1] = affine2evalue(z->p[1], z->p[1][nparam+1], nparam);
730 evalue *up = evalue_dup(res);
731 evalue_substitute(up, subs_0d+1);
732 evalue_free(subs_0d[1]);
734 subs_0d[1] = affine2evalue(z->p[0], z->p[0][nparam+1], nparam);
735 evalue_substitute(res, subs_0d+1);
736 evalue_negate(res);
737 eadd(up, res);
738 evalue_free(subs_0d[1]);
739 subs_0d[1] = NULL;
740 evalue_free(up);
742 Matrix_Free(z);
744 eadd(res, sum);
745 evalue_free(res);
748 /* Integrate the polynomial over the whole polygon using
749 * the Green-Stokes theorem.
751 void summator_2d::integrate(Polyhedron *P, Param_Domain *PD)
753 Value tmp;
754 evalue *res = evalue_zero();
756 evalue *I = evalue_dup(polynomial);
757 evalue_anti_derive(I, 0);
759 value_init(tmp);
760 Vector *normal = Vector_Alloc(2);
761 Vector *dir = Vector_Alloc(2);
762 Matrix *v0v1 = Matrix_Alloc(2, nparam+2);
763 Vector *f_v0v1 = Vector_Alloc(2);
764 Vector *s = Vector_Alloc(1+nparam+2);
765 for (int j = P->NbEq; j < P->NbConstraints; ++j) {
766 Param_Domain *FD;
767 int nbV = 0;
768 Param_Vertices *vertex[2];
769 Vector_Copy(P->Constraint[j]+1, normal->p, 2);
771 if (value_zero_p(normal->p[0]))
772 continue;
774 Vector_Normalize(normal->p, 2);
775 value_assign(dir->p[0], normal->p[1]);
776 value_oppose(dir->p[1], normal->p[0]);
778 FD = Param_Polyhedron_Facet(PP, PD, P, j);
780 FORALL_PVertex_in_ParamPolyhedron(V, FD, PP)
781 vertex[nbV++] = V;
782 END_FORALL_PVertex_in_ParamPolyhedron;
784 assert(nbV == 2);
786 Param_Vertex_Common_Denominator(vertex[0]);
787 Param_Vertex_Common_Denominator(vertex[1]);
789 value_oppose(tmp, vertex[1]->Vertex->p[0][nparam+1]);
790 for (int i = 0; i < 2; ++i)
791 Vector_Combine(vertex[1]->Vertex->p[i],
792 vertex[0]->Vertex->p[i],
793 v0v1->p[i],
794 vertex[0]->Vertex->p[0][nparam+1], tmp, nparam+1);
795 value_multiply(v0v1->p[0][nparam+1],
796 vertex[0]->Vertex->p[0][nparam+1],
797 vertex[1]->Vertex->p[0][nparam+1]);
798 value_assign(v0v1->p[1][nparam+1], v0v1->p[0][nparam+1]);
800 /* Order vertices to ensure we integrate in the right
801 * direction, i.e., with the polytope "on the left".
803 for (int i = 0; i < 2; ++i)
804 Inner_Product(v0v1->p[i], inner, nparam+1, &f_v0v1->p[i]);
806 Inner_Product(dir->p, f_v0v1->p, 2, &tmp);
807 if (value_neg_p(tmp)) {
808 Param_Vertices *PV = vertex[0];
809 vertex[0] = vertex[1];
810 vertex[1] = PV;
811 for (int i = 0; i < 2; ++i)
812 Vector_Oppose(v0v1->p[i], v0v1->p[i], nparam+1);
814 value_oppose(tmp, normal->p[0]);
815 if (value_neg_p(tmp)) {
816 value_oppose(tmp, tmp);
817 Vector_Oppose(v0v1->p[1], v0v1->p[1], nparam+1);
819 value_multiply(tmp, tmp, v0v1->p[1][nparam+1]);
820 evalue *top = affine2evalue(v0v1->p[1], tmp, nparam);
822 value_multiply(s->p[0], normal->p[1], vertex[0]->Vertex->p[0][nparam+1]);
823 Vector_Copy(vertex[0]->Vertex->p[0], s->p+1, nparam+2);
824 subs_1d[0] = affine2evalue(s->p, s->p[1+nparam+1], 1+nparam);
825 value_multiply(s->p[0], normal->p[0], vertex[0]->Vertex->p[0][nparam+1]);
826 value_oppose(s->p[0], s->p[0]);
827 Vector_Copy(vertex[0]->Vertex->p[1], s->p+1, nparam+2);
828 subs_1d[1] = affine2evalue(s->p, s->p[1+nparam+1], 1+nparam);
830 evalue *d = evalue_dup(I);
831 evalue_substitute(d, subs_1d);
832 evalue_anti_derive(d, 0);
834 evalue_free(subs_1d[0]);
835 evalue_free(subs_1d[1]);
837 subs_0d[1] = top;
838 evalue_substitute(d, subs_0d+1);
839 evalue_mul(d, dir->p[1]);
840 evalue_free(subs_0d[1]);
842 eadd(d, res);
843 evalue_free(d);
845 Param_Domain_Free(FD);
847 Vector_Free(s);
848 Vector_Free(f_v0v1);
849 Matrix_Free(v0v1);
850 Vector_Free(normal);
851 Vector_Free(dir);
852 value_clear(tmp);
853 evalue_free(I);
855 eadd(res, sum);
856 evalue_free(res);
859 evalue *summate_over_1d_pdomain(evalue *e,
860 Param_Polyhedron *PP, Param_Domain *PD,
861 Polyhedron *P, Value *inner,
862 struct barvinok_options *options)
864 Param_Vertices *V;
865 int nbV = 0;
866 Param_Vertices *vertex[2];
867 unsigned nparam = PP->V->Vertex->NbColumns-2;
868 evalue *subs_0d[1+nparam];
869 evalue *a[2];
870 evalue *t[2];
871 unsigned degree = total_degree(e, 1);
873 for (int i = 0; i < nparam; ++i)
874 subs_0d[1+i] = evalue_var(i);
876 FORALL_PVertex_in_ParamPolyhedron(V, PD, PP)
877 vertex[nbV++] = V;
878 END_FORALL_PVertex_in_ParamPolyhedron;
879 assert(nbV == 2);
881 Vector *fixed = Vector_Alloc(2);
882 for (int i = 0; i < 2; ++i) {
883 Inner_Product(vertex[i]->Vertex->p[0], inner, nparam+1, &fixed->p[i]);
884 value_multiply(fixed->p[i],
885 fixed->p[i], vertex[1-i]->Vertex->p[0][nparam+1]);
887 if (value_lt(fixed->p[1], fixed->p[0])) {
888 V = vertex[0];
889 vertex[0] = vertex[1];
890 vertex[1] = V;
892 Vector_Free(fixed);
894 Vector *coef = Vector_Alloc(nparam+1);
895 for (int i = 0; i < 2; ++i)
896 a[i] = affine2evalue(vertex[i]->Vertex->p[0],
897 vertex[i]->Vertex->p[0][nparam+1], nparam);
898 if (value_one_p(vertex[0]->Vertex->p[0][nparam+1]))
899 t[0] = evalue_zero();
900 else {
901 Vector_Oppose(vertex[0]->Vertex->p[0], coef->p, nparam+1);
902 t[0] = fractional_part(coef->p, vertex[0]->Vertex->p[0][nparam+1],
903 nparam, NULL);
905 if (value_one_p(vertex[1]->Vertex->p[0][nparam+1]))
906 t[1] = evalue_zero();
907 else {
908 Vector_Copy(vertex[1]->Vertex->p[0], coef->p, nparam+1);
909 t[1] = fractional_part(coef->p, vertex[1]->Vertex->p[0][nparam+1],
910 nparam, NULL);
912 Vector_Free(coef);
914 evalue *I = evalue_dup(e);
915 evalue_anti_derive(I, 0);
916 evalue *up = evalue_dup(I);
917 subs_0d[0] = a[1];
918 evalue_substitute(up, subs_0d);
920 subs_0d[0] = a[0];
921 evalue_substitute(I, subs_0d);
922 evalue_negate(I);
923 eadd(up, I);
924 evalue_free(up);
926 evalue *res = I;
928 mu_1d mu0(degree, t[0]);
929 mu_1d mu1(degree, t[1]);
931 evalue *dx = evalue_dup(e);
932 for (int n = 0; !EVALUE_IS_ZERO(*dx); ++n) {
933 evalue *d;
935 d = evalue_dup(dx);
936 subs_0d[0] = a[0];
937 evalue_substitute(d, subs_0d);
938 emul(mu0.coefficient(n), d);
939 eadd(d, res);
940 evalue_free(d);
942 d = evalue_dup(dx);
943 subs_0d[0] = a[1];
944 evalue_substitute(d, subs_0d);
945 emul(mu1.coefficient(n), d);
946 if (n % 2)
947 evalue_negate(d);
948 eadd(d, res);
949 evalue_free(d);
951 evalue_derive(dx, 0);
953 evalue_free(dx);
955 for (int i = 0; i < nparam; ++i)
956 evalue_free(subs_0d[1+i]);
958 for (int i = 0; i < 2; ++i) {
959 evalue_free(a[i]);
960 evalue_free(t[i]);
963 return res;
966 static evalue *summate_over_domain(evalue *e, int nvar, Polyhedron *D,
967 struct barvinok_options *options)
969 Polyhedron *U;
970 Param_Polyhedron *PP;
971 evalue *res;
972 int nd = -1;
973 unsigned MaxRays;
974 struct evalue_section *s;
975 Param_Domain *PD;
977 MaxRays = options->MaxRays;
978 POL_UNSET(options->MaxRays, POL_INTEGER);
980 U = Universe_Polyhedron(D->Dimension - nvar);
981 PP = Polyhedron2Param_Polyhedron(D, U, options);
983 for (nd = 0, PD = PP->D; PD; ++nd, PD = PD->next);
984 s = ALLOCN(struct evalue_section, nd);
986 Polyhedron *TC = true_context(D, U, MaxRays);
987 FORALL_REDUCED_DOMAIN(PP, TC, nd, options, i, PD, rVD)
988 Polyhedron *CA, *F;
990 CA = align_context(PD->Domain, D->Dimension, options->MaxRays);
991 F = DomainIntersection(D, CA, options->MaxRays);
992 Domain_Free(CA);
994 Vector *inner = inner_point(rVD);
995 s[i].D = rVD;
997 if (nvar == 1) {
998 s[i].E = summate_over_1d_pdomain(e, PP, PD, F, inner->p+1, options);
999 } else if (nvar == 2) {
1000 summator_2d s2d(e, PP, inner->p+1, rVD->Dimension);
1002 s[i].E = s2d.summate_over_pdomain(F, PD, options);
1005 Polyhedron_Free(F);
1006 Vector_Free(inner);
1007 END_FORALL_REDUCED_DOMAIN
1008 Polyhedron_Free(TC);
1009 Polyhedron_Free(U);
1010 Param_Polyhedron_Free(PP);
1012 options->MaxRays = MaxRays;
1014 res = evalue_from_section_array(s, nd);
1015 free(s);
1017 return res;
1020 evalue *euler_summate(evalue *e, int nvar, struct barvinok_options *options)
1022 int i;
1023 evalue *res;
1025 assert(nvar <= 2);
1027 assert(nvar >= 0);
1028 if (nvar == 0 || EVALUE_IS_ZERO(*e))
1029 return evalue_dup(e);
1031 assert(value_zero_p(e->d));
1032 assert(e->x.p->type == partition);
1034 res = evalue_zero();
1036 for (i = 0; i < e->x.p->size/2; ++i) {
1037 evalue *t;
1038 t = summate_over_domain(&e->x.p->arr[2*i+1], nvar,
1039 EVALUE_DOMAIN(e->x.p->arr[2*i]), options);
1040 eadd(t, res);
1041 evalue_free(t);
1044 return res;