power.h: extract from euler.cc
[barvinok.git] / euler.cc
blob1f4f602a11578769152065e7b5797b9c6291c72e
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 "power.h"
16 #include "reduce_domain.h"
18 /* Compute total degree in first nvar variables */
19 static unsigned total_degree(const evalue *e, unsigned nvar)
21 int i;
22 unsigned max_degree;
24 if (value_notzero_p(e->d))
25 return 0;
26 assert(e->x.p->type == polynomial);
27 if (e->x.p->pos-1 >= nvar)
28 return 0;
30 max_degree = total_degree(&e->x.p->arr[0], nvar);
31 for (i = 1; i < e->x.p->size; ++i) {
32 unsigned degree = i + total_degree(&e->x.p->arr[i], nvar);
33 if (degree > max_degree)
34 max_degree = degree;
36 return max_degree;
39 struct fact {
40 Value *fact;
41 unsigned size;
42 unsigned n;
44 struct fact fact;
46 Value *factorial(unsigned n)
48 if (n < fact.n)
49 return &fact.fact[n];
51 if (n >= fact.size) {
52 int size = 3*(n + 5)/2;
54 fact.fact = (Value *)realloc(fact.fact, size*sizeof(Value));
55 fact.size = size;
57 for (int i = fact.n; i <= n; ++i) {
58 value_init(fact.fact[i]);
59 if (!i)
60 value_set_si(fact.fact[0], 1);
61 else
62 mpz_mul_ui(fact.fact[i], fact.fact[i-1], i);
64 fact.n = n+1;
65 return &fact.fact[n];
68 struct binom {
69 Vector **binom;
70 unsigned size;
71 unsigned n;
73 struct binom binom;
75 Value *binomial(unsigned n, unsigned k)
77 if (n < binom.n)
78 return &binom.binom[n]->p[k];
80 if (n >= binom.size) {
81 int size = 3*(n + 5)/2;
83 binom.binom = (Vector **)realloc(binom.binom, size*sizeof(Vector *));
84 binom.size = size;
86 for (int i = binom.n; i <= n; ++i) {
87 binom.binom[i] = Vector_Alloc(i+1);
88 if (!i)
89 value_set_si(binom.binom[0]->p[0], 1);
90 else {
91 value_set_si(binom.binom[i]->p[0], 1);
92 value_set_si(binom.binom[i]->p[i], 1);
93 for (int j = 1; j < i; ++j)
94 value_addto(binom.binom[i]->p[j],
95 binom.binom[i-1]->p[j-1], binom.binom[i-1]->p[j]);
98 binom.n = n+1;
99 return &binom.binom[n]->p[k];
103 * Computes the coefficients of
105 * \mu(-t + R_+)(\xi) = \sum_{n=0)^\infty -b(n+1, t)/(n+1)! \xi^n
107 * where b(n, t) is the Bernoulli polynomial of degree n in t
108 * and t(p) is an expression (a fractional part) of the parameters p
109 * such that 0 <= t(p) < 1 for all values of the parameters.
110 * The coefficients are computed on demand up to (and including)
111 * the maximal degree max_degree.
113 struct mu_1d {
114 unsigned max_degree;
115 evalue **coefficients;
116 const evalue *t;
118 mu_1d(unsigned max_degree, evalue *t) : max_degree(max_degree), t(t) {
119 coefficients = new evalue *[max_degree+1];
120 for (int i = 0; i < max_degree+1; ++i)
121 coefficients[i] = NULL;
123 ~mu_1d() {
124 for (int i = 0; i < max_degree+1; ++i)
125 if (coefficients[i])
126 evalue_free(coefficients[i]);
127 delete [] coefficients;
129 void compute_coefficient(unsigned n);
130 const evalue *coefficient(unsigned n) {
131 if (!coefficients[n])
132 compute_coefficient(n);
133 return coefficients[n];
137 void mu_1d::compute_coefficient(unsigned n)
139 struct poly_list *bernoulli = bernoulli_compute(n+1);
140 evalue *c = evalue_polynomial(bernoulli->poly[n+1], t);
141 evalue_negate(c);
142 evalue_div(c, *factorial(n+1));
144 coefficients[n] = c;
148 * Computes the coefficients of
150 * \mu(a)(y) = \sum_{n_1} \sum_{n_2} c_{n_1,n_2} y^{n_1} y^{n_2}
152 * with c_{n1,n2} given by
154 * b(n1+1,t1)/(n1+1)! b(n2+1,t2)/(n2+1)!
155 * - b(n1+n2+2,t2)/(n1+n2+2)! (-c1)^{n1+1}
156 * - b(n1+n2+2,t1)/(n1+n2+2)! (-c2)^{n2+1}
158 * where b(n, t) is the Bernoulli polynomial of degree n in t,
159 * t1(p) and t2(p) are expressions (fractional parts) of the parameters p
160 * such that 0 <= t1(p), t2(p) < 1 for all values of the parameters
161 * and c1 = cn/c1d and c2 = cn/c2d.
162 * The coefficients are computed on demand up to (and including)
163 * the maximal degree (n1,n2) = (max_degree,max_degree).
165 * bernoulli_t[i][j] contains b(j+1, t_i)/(j+1)!
167 struct mu_2d {
168 unsigned max_degree;
169 evalue ***coefficients;
170 /* bernoulli_t[i][n] stores b(n+1, t_i)/(n+1)! */
171 evalue **bernoulli_t[2];
172 /* stores the powers of -cn */
173 struct power *power_cn;
174 struct power *power_c1d;
175 struct power *power_c2d;
176 const evalue *t[2];
178 mu_2d(unsigned max_degree, evalue *t1, evalue *t2,
179 Value cn, Value c1d, Value c2d) : max_degree(max_degree) {
180 t[0] = t1;
181 t[1] = t2;
182 coefficients = new evalue **[max_degree+1];
183 for (int i = 0; i < max_degree+1; ++i) {
184 coefficients[i] = new evalue *[max_degree+1];
185 for (int j = 0; j < max_degree+1; ++j)
186 coefficients[i][j] = NULL;
188 for (int i = 0; i < 2; ++i) {
189 bernoulli_t[i] = new evalue *[max_degree+2];
190 for (int j = 0; j < max_degree+2; ++j)
191 bernoulli_t[i][j] = NULL;
193 value_oppose(cn, cn);
194 power_cn = new struct power(cn, max_degree+1);
195 value_oppose(cn, cn);
196 power_c1d = new struct power(c1d, max_degree+1);
197 power_c2d = new struct power(c2d, max_degree+1);
199 ~mu_2d() {
200 for (int i = 0; i < max_degree+1; ++i) {
201 for (int j = 0; j < max_degree+1; ++j)
202 if (coefficients[i][j])
203 evalue_free(coefficients[i][j]);
204 delete [] coefficients[i];
206 delete [] coefficients;
207 for (int i = 0; i < 2; ++i)
208 for (int j = 0; j < max_degree+2; ++j)
209 if (bernoulli_t[i][j])
210 evalue_free(bernoulli_t[i][j]);
211 for (int i = 0; i < 2; ++i)
212 delete [] bernoulli_t[i];
213 delete power_cn;
214 delete power_c1d;
215 delete power_c2d;
217 const evalue *get_bernoulli(int n, int i);
219 void compute_coefficient(unsigned n1, unsigned n2);
220 const evalue *coefficient(unsigned n1, unsigned n2) {
221 if (!coefficients[n1][n2])
222 compute_coefficient(n1, n2);
223 return coefficients[n1][n2];
228 * Returns b(n, t_i)/n!
230 const evalue *mu_2d::get_bernoulli(int n, int i)
232 if (!bernoulli_t[i][n-1]) {
233 struct poly_list *bernoulli = bernoulli_compute(n);
234 bernoulli_t[i][n-1] = evalue_polynomial(bernoulli->poly[n], t[i]);
235 evalue_div(bernoulli_t[i][n-1], *factorial(n));
237 return bernoulli_t[i][n-1];
240 void mu_2d::compute_coefficient(unsigned n1, unsigned n2)
242 evalue *c = evalue_dup(get_bernoulli(n1+1, 0));
243 emul(get_bernoulli(n2+1, 1), c);
245 if (value_notzero_p(*(*power_cn)[1])) {
246 evalue *t;
247 Value neg_power;
249 value_init(neg_power);
251 t = evalue_dup(get_bernoulli(n1+n2+2, 1));
252 value_multiply(neg_power,
253 *(*power_cn)[n1+1], *binomial(n1+n2+1, n1+1));
254 value_oppose(neg_power, neg_power);
255 evalue_mul_div(t, neg_power, *(*power_c1d)[n1+1]);
256 eadd(t, c);
257 evalue_free(t);
259 t = evalue_dup(get_bernoulli(n1+n2+2, 0));
260 value_multiply(neg_power,
261 *(*power_cn)[n2+1], *binomial(n1+n2+1, n2+1));
262 value_oppose(neg_power, neg_power);
263 evalue_mul_div(t, neg_power, *(*power_c2d)[n2+1]);
264 eadd(t, c);
265 evalue_free(t);
267 value_clear(neg_power);
270 coefficients[n1][n2] = c;
273 /* Later: avoid recomputation of mu of faces that appear in
274 * more than one domain.
276 struct summator_2d : public signed_cone_consumer, public vertex_decomposer {
277 const evalue *polynomial;
278 Value *inner;
279 unsigned nparam;
281 /* substitutions to use when result is 0-dimensional (only parameters) */
282 evalue **subs_0d;
283 /* substitutions to use when result is 1-dimensional */
284 evalue **subs_1d;
285 evalue *sum;
287 summator_2d(evalue *e, Param_Polyhedron *PP, Value *inner,
288 unsigned nparam) :
289 polynomial(e), vertex_decomposer(PP, *this),
290 inner(inner), nparam(nparam) {
291 sum = evalue_zero();
293 subs_0d = new evalue *[2+nparam];
294 subs_1d = new evalue *[2+nparam];
295 subs_0d[0] = NULL;
296 subs_0d[1] = NULL;
297 subs_1d[0] = NULL;
298 subs_1d[1] = NULL;
299 for (int i = 0; i < nparam; ++i) {
300 subs_0d[2+i] = evalue_var(i);
301 subs_1d[2+i] = evalue_var(1+i);
304 ~summator_2d() {
305 for (int i = 0; i < nparam; ++i) {
306 evalue_free(subs_0d[2+i]);
307 evalue_free(subs_1d[2+i]);
309 delete [] subs_0d;
310 delete [] subs_1d;
312 evalue *summate_over_pdomain(Param_Polyhedron *PP, unsigned *facets,
313 Param_Domain *PD,
314 struct barvinok_options *options);
315 void handle_facet(Param_Polyhedron *PP, Param_Domain *FD, Value *normal);
316 void integrate(Param_Polyhedron *PP, unsigned *facets, Param_Domain *PD);
317 virtual void handle(const signed_cone& sc, barvinok_options *options);
320 /* Replaces poly by its derivative along variable var */
321 static void evalue_derive(evalue *poly, int var)
323 if (value_notzero_p(poly->d)) {
324 value_set_si(poly->x.n, 0);
325 value_set_si(poly->d, 1);
326 return;
328 assert(poly->x.p->type == polynomial);
329 if (poly->x.p->pos-1 > var) {
330 free_evalue_refs(poly);
331 value_init(poly->d);
332 evalue_set_si(poly, 0, 1);
333 return;
336 if (poly->x.p->pos-1 < var) {
337 for (int i = 0; i < poly->x.p->size; ++i)
338 evalue_derive(&poly->x.p->arr[i], var);
339 reduce_evalue(poly);
340 return;
343 assert(poly->x.p->size >= 1);
344 enode *p = poly->x.p;
345 free_evalue_refs(&p->arr[0]);
346 if (p->size == 1) {
347 free(p);
348 evalue_set_si(poly, 0, 1);
349 return;
351 if (p->size == 2) {
352 value_clear(poly->d);
353 *poly = p->arr[1];
354 free(p);
355 return;
358 p->size--;
359 Value factor;
360 value_init(factor);
361 for (int i = 0; i < p->size; ++i) {
362 value_set_si(factor, i+1);
363 p->arr[i] = p->arr[i+1];
364 evalue_mul(&p->arr[i], factor);
366 value_clear(factor);
369 /* Check whether e is constant with respect to variable var */
370 static int evalue_is_constant(evalue *e, int var)
372 int i;
374 if (value_notzero_p(e->d))
375 return 1;
376 if (e->x.p->type == polynomial && e->x.p->pos-1 == var)
377 return 0;
378 assert(e->x.p->type == polynomial ||
379 e->x.p->type == fractional ||
380 e->x.p->type == flooring);
381 for (i = 0; i < e->x.p->size; ++i)
382 if (!evalue_is_constant(&e->x.p->arr[i], var))
383 return 0;
384 return 1;
387 /* Replaces poly by its anti-derivative with constant 0 along variable var */
388 static void evalue_anti_derive(evalue *poly, int var)
390 enode *p;
392 if (value_zero_p(poly->d) &&
393 poly->x.p->type == polynomial && poly->x.p->pos-1 < var) {
394 for (int i = 0; i < poly->x.p->size; ++i)
395 evalue_anti_derive(&poly->x.p->arr[i], var);
396 reduce_evalue(poly);
397 return;
400 if (evalue_is_constant(poly, var)) {
401 p = new_enode(polynomial, 2, 1+var);
402 evalue_set_si(&p->arr[0], 0, 1);
403 value_clear(p->arr[1].d);
404 p->arr[1] = *poly;
405 value_init(poly->d);
406 poly->x.p = p;
407 return;
409 assert(poly->x.p->type == polynomial);
411 p = new_enode(polynomial, poly->x.p->size+1, 1+var);
412 evalue_set_si(&p->arr[0], 0, 1);
413 for (int i = 0; i < poly->x.p->size; ++i) {
414 value_clear(p->arr[1+i].d);
415 p->arr[1+i] = poly->x.p->arr[i];
416 value_set_si(poly->d, 1+i);
417 evalue_div(&p->arr[1+i], poly->d);
419 free(poly->x.p);
420 poly->x.p = p;
421 value_set_si(poly->d, 0);
424 /* Computes offsets in the basis given by the rays of the cone
425 * to the integer point in the fundamental parallelepiped of
426 * the cone.
427 * The resulting evalues contain only the parameters as variables.
429 evalue **offsets_to_integer_point(Matrix *Rays, Matrix *vertex)
431 unsigned nparam = vertex->NbColumns-2;
432 evalue **t = new evalue *[2];
434 if (value_one_p(vertex->p[0][nparam+1])) {
435 t[0] = evalue_zero();
436 t[1] = evalue_zero();
437 } else {
438 Matrix *R2 = Matrix_Copy(Rays);
439 Matrix_Transposition(R2);
440 Matrix *inv = Matrix_Alloc(Rays->NbColumns, Rays->NbRows);
441 int ok = Matrix_Inverse(R2, inv);
442 assert(ok);
443 Matrix_Free(R2);
445 /* We want the fractional part of the negative relative coordinates */
446 Vector_Oppose(inv->p[0], inv->p[0], inv->NbColumns);
447 Vector_Oppose(inv->p[1], inv->p[1], inv->NbColumns);
449 Matrix *neg_rel = Matrix_Alloc(2, nparam+1);
450 vertex->NbColumns--;
451 Matrix_Product(inv, vertex, neg_rel);
452 Matrix_Free(inv);
453 vertex->NbColumns++;
454 t[0] = fractional_part(neg_rel->p[0], vertex->p[0][nparam+1],
455 nparam, NULL);
456 t[1] = fractional_part(neg_rel->p[1], vertex->p[0][nparam+1],
457 nparam, NULL);
458 Matrix_Free(neg_rel);
461 return t;
465 * Called from decompose_at_vertex.
467 * Handles a cone in the signed decomposition of the supporting
468 * cone of a vertex. The cone is assumed to be unimodular.
470 void summator_2d::handle(const signed_cone& sc, barvinok_options *options)
472 evalue **t;
473 unsigned degree = total_degree(polynomial, 2);
475 subs_0d[0] = affine2evalue(V->Vertex->p[0],
476 V->Vertex->p[0][nparam+1], nparam);
477 subs_0d[1] = affine2evalue(V->Vertex->p[1],
478 V->Vertex->p[1][nparam+1], nparam);
480 assert(sc.det == 1);
482 assert(V->Vertex->NbRows > 0);
483 Param_Vertex_Common_Denominator(V);
485 Matrix *Rays = zz2matrix(sc.rays);
487 t = offsets_to_integer_point(Rays, V->Vertex);
489 Vector *c = Vector_Alloc(3);
490 Inner_Product(Rays->p[0], Rays->p[1], 2, &c->p[0]);
491 Inner_Product(Rays->p[0], Rays->p[0], 2, &c->p[1]);
492 Inner_Product(Rays->p[1], Rays->p[1], 2, &c->p[2]);
494 mu_2d mu(degree, t[0], t[1], c->p[0], c->p[1], c->p[2]);
495 Vector_Free(c);
497 struct power power_r00(Rays->p[0][0], degree);
498 struct power power_r01(Rays->p[0][1], degree);
499 struct power power_r10(Rays->p[1][0], degree);
500 struct power power_r11(Rays->p[1][1], degree);
502 Value factor, tmp1, tmp2;
503 value_init(factor);
504 value_init(tmp1);
505 value_init(tmp2);
506 evalue *res = evalue_zero();
507 evalue *dx1 = evalue_dup(polynomial);
508 for (int i = 0; !EVALUE_IS_ZERO(*dx1); ++i) {
509 evalue *dx2 = evalue_dup(dx1);
510 for (int j = 0; !EVALUE_IS_ZERO(*dx2); ++j) {
511 evalue *cij = evalue_zero();
512 for (int n1 = 0; n1 <= i+j; ++n1) {
513 int n2 = i+j-n1;
514 value_set_si(factor, 0);
515 for (int k = max(0, i-n2); k <= i && k <= n1; ++k) {
516 int l = i-k;
517 value_multiply(tmp1, *power_r00[k], *power_r01[n1-k]);
518 value_multiply(tmp1, tmp1, *power_r10[l]);
519 value_multiply(tmp1, tmp1, *power_r11[n2-l]);
520 value_multiply(tmp1, tmp1, *binomial(n1, k));
521 value_multiply(tmp1, tmp1, *binomial(n2, l));
522 value_addto(factor, factor, tmp1);
524 if (value_zero_p(factor))
525 continue;
527 evalue *c = evalue_dup(mu.coefficient(n1, n2));
528 evalue_mul(c, factor);
529 eadd(c, cij);
530 evalue_free(c);
532 evalue *d = evalue_dup(dx2);
533 evalue_substitute(d, subs_0d);
534 emul(cij, d);
535 evalue_free(cij);
536 eadd(d, res);
537 evalue_free(d);
538 evalue_derive(dx2, 1);
540 evalue_free(dx2);
541 evalue_derive(dx1, 0);
543 evalue_free(dx1);
544 for (int i = 0; i < 2; ++i) {
545 evalue_free(subs_0d[i]);
546 subs_0d[i] = NULL;
547 evalue_free(t[i]);
549 delete [] t;
550 value_clear(factor);
551 value_clear(tmp1);
552 value_clear(tmp2);
553 Matrix_Free(Rays);
555 if (sc.sign < 0)
556 evalue_negate(res);
557 eadd(res, sum);
558 evalue_free(res);
561 evalue *summator_2d::summate_over_pdomain(Param_Polyhedron *PP,
562 unsigned *facets, Param_Domain *PD,
563 struct barvinok_options *options)
565 Param_Vertices *V;
566 int i, ix;
567 unsigned bx;
569 assert(PP->V->Vertex->NbRows == 2);
571 FORALL_PVertex_in_ParamPolyhedron(V, PD, PP) // _i internal counter
572 decompose_at_vertex(V, _i, options);
573 END_FORALL_PVertex_in_ParamPolyhedron;
575 Vector *normal = Vector_Alloc(2);
576 for (i = 0, ix = 0, bx = MSB; i < PP->Constraints->NbRows; ++i) {
577 Param_Domain *FD;
579 if (!(facets[ix] & bx)) {
580 NEXT(ix, bx);
581 continue;
584 Vector_Copy(PP->Constraints->p[i]+1, normal->p, 2);
585 if (value_zero_p(normal->p[0]) && value_zero_p(normal->p[1]))
586 continue;
588 FD = Param_Polyhedron_Facet(PP, PD, PP->Constraints->p[i]);
589 Vector_Normalize(normal->p, 2);
590 handle_facet(PP, FD, normal->p);
591 Param_Domain_Free(FD);
592 NEXT(ix, bx);
594 Vector_Free(normal);
596 integrate(PP, facets, PD);
598 return sum;
601 void summator_2d::handle_facet(Param_Polyhedron *PP, Param_Domain *FD,
602 Value *normal)
604 Param_Vertices *V;
605 int nbV = 0;
606 Param_Vertices *vertex[2];
607 Value det;
608 unsigned degree = total_degree(polynomial, 2);
610 FORALL_PVertex_in_ParamPolyhedron(V, FD, PP)
611 vertex[nbV++] = V;
612 END_FORALL_PVertex_in_ParamPolyhedron;
614 assert(nbV == 2);
616 /* We can take either vertex[0] or vertex[1];
617 * the result should be the same
619 Param_Vertex_Common_Denominator(vertex[0]);
621 /* The extra variable in front is the coordinate along the facet. */
622 Vector *coef_normal = Vector_Alloc(1 + nparam + 2);
623 Vector_Combine(vertex[0]->Vertex->p[0], vertex[0]->Vertex->p[1],
624 coef_normal->p+1, normal[0], normal[1], nparam+1);
625 value_assign(coef_normal->p[1+nparam+1], vertex[0]->Vertex->p[0][nparam+1]);
626 Vector_Normalize(coef_normal->p, coef_normal->Size);
628 Vector *base = Vector_Alloc(2);
629 value_oppose(base->p[0], normal[1]);
630 value_assign(base->p[1], normal[0]);
632 value_init(det);
633 Inner_Product(normal, normal, 2, &det);
635 Vector *s = Vector_Alloc(1+nparam+2);
636 value_multiply(s->p[1+nparam+1], coef_normal->p[1+nparam+1], det);
637 value_multiply(s->p[0], base->p[0], s->p[1+nparam+1]);
638 Vector_Scale(coef_normal->p+1, s->p+1, normal[0], nparam+1);
639 subs_1d[0] = affine2evalue(s->p, s->p[1+nparam+1], 1+nparam);
640 value_multiply(s->p[0], base->p[1], s->p[1+nparam+1]);
641 Vector_Scale(coef_normal->p+1, s->p+1, normal[1], nparam+1);
642 subs_1d[1] = affine2evalue(s->p, s->p[1+nparam+1], 1+nparam);
643 Vector_Free(s);
645 evalue *t;
646 if (value_one_p(coef_normal->p[coef_normal->Size-1]))
647 t = evalue_zero();
648 else {
649 Vector_Oppose(coef_normal->p+1, coef_normal->p+1, nparam+1);
650 t = fractional_part(coef_normal->p,
651 coef_normal->p[coef_normal->Size-1],
652 1+nparam, NULL);
654 Vector_Free(coef_normal);
656 mu_1d mu(degree, t);
658 struct power power_normal0(normal[0], degree);
659 struct power power_normal1(normal[1], degree);
660 struct power power_det(det, degree);
662 Value(factor);
663 value_init(factor);
664 evalue *res = evalue_zero();
665 evalue *dx1 = evalue_dup(polynomial);
666 for (int i = 0; !EVALUE_IS_ZERO(*dx1); ++i) {
667 evalue *dx2 = evalue_dup(dx1);
668 for (int j = 0; !EVALUE_IS_ZERO(*dx2); ++j) {
669 value_multiply(factor, *power_normal0[i], *power_normal1[j]);
670 if (value_notzero_p(factor)) {
671 value_multiply(factor, factor, *binomial(i+j, i));
673 evalue *c = evalue_dup(mu.coefficient(i+j));
674 evalue_mul_div(c, factor, *power_det[i+j]);
676 evalue *d = evalue_dup(dx2);
677 evalue_substitute(d, subs_1d);
678 emul(c, d);
679 evalue_free(c);
680 eadd(d, res);
681 evalue_free(d);
683 evalue_derive(dx2, 1);
685 evalue_free(dx2);
686 evalue_derive(dx1, 0);
688 evalue_free(dx1);
689 evalue_free(t);
690 for (int i = 0; i < 2; ++i) {
691 evalue_free(subs_1d[i]);
692 subs_1d[i] = NULL;
694 value_clear(factor);
695 evalue_anti_derive(res, 0);
697 Matrix *z = Matrix_Alloc(2, nparam+2);
698 Vector *fixed_z = Vector_Alloc(2);
699 for (int i = 0; i < 2; ++i) {
700 Vector_Combine(vertex[i]->Vertex->p[0], vertex[i]->Vertex->p[1],
701 z->p[i], base->p[0], base->p[1], nparam+1);
702 value_multiply(z->p[i][nparam+1],
703 det, vertex[i]->Vertex->p[0][nparam+1]);
704 Inner_Product(z->p[i], inner, nparam+1, &fixed_z->p[i]);
706 Vector_Free(base);
707 /* Put on a common denominator */
708 value_multiply(fixed_z->p[0], fixed_z->p[0], z->p[1][nparam+1]);
709 value_multiply(fixed_z->p[1], fixed_z->p[1], z->p[0][nparam+1]);
710 /* Make sure z->p[0] is smaller than z->p[1] (for an internal
711 * point of the chamber and hence for all parameter values in
712 * the chamber), to ensure we integrate in the right direction.
714 if (value_lt(fixed_z->p[1], fixed_z->p[0]))
715 Vector_Exchange(z->p[0], z->p[1], nparam+2);
716 Vector_Free(fixed_z);
717 value_clear(det);
719 subs_0d[1] = affine2evalue(z->p[1], z->p[1][nparam+1], nparam);
720 evalue *up = evalue_dup(res);
721 evalue_substitute(up, subs_0d+1);
722 evalue_free(subs_0d[1]);
724 subs_0d[1] = affine2evalue(z->p[0], z->p[0][nparam+1], nparam);
725 evalue_substitute(res, subs_0d+1);
726 evalue_negate(res);
727 eadd(up, res);
728 evalue_free(subs_0d[1]);
729 subs_0d[1] = NULL;
730 evalue_free(up);
732 Matrix_Free(z);
734 eadd(res, sum);
735 evalue_free(res);
738 /* Integrate the polynomial over the whole polygon using
739 * the Green-Stokes theorem.
741 void summator_2d::integrate(Param_Polyhedron *PP, unsigned *facets,
742 Param_Domain *PD)
744 Value tmp;
745 evalue *res = evalue_zero();
746 int i, ix;
747 unsigned bx;
749 evalue *I = evalue_dup(polynomial);
750 evalue_anti_derive(I, 0);
752 value_init(tmp);
753 Vector *normal = Vector_Alloc(2);
754 Vector *dir = Vector_Alloc(2);
755 Matrix *v0v1 = Matrix_Alloc(2, nparam+2);
756 Vector *f_v0v1 = Vector_Alloc(2);
757 Vector *s = Vector_Alloc(1+nparam+2);
758 for (i = 0, ix = 0, bx = MSB; i < PP->Constraints->NbRows; ++i) {
759 Param_Domain *FD;
760 int nbV = 0;
761 Param_Vertices *vertex[2];
763 if (!(facets[ix] & bx)) {
764 NEXT(ix, bx);
765 continue;
768 Vector_Copy(PP->Constraints->p[i]+1, normal->p, 2);
770 if (value_zero_p(normal->p[0]))
771 continue;
773 Vector_Normalize(normal->p, 2);
774 value_assign(dir->p[0], normal->p[1]);
775 value_oppose(dir->p[1], normal->p[0]);
777 FD = Param_Polyhedron_Facet(PP, PD, PP->Constraints->p[i]);
779 FORALL_PVertex_in_ParamPolyhedron(V, FD, PP)
780 vertex[nbV++] = V;
781 END_FORALL_PVertex_in_ParamPolyhedron;
783 assert(nbV == 2);
785 Param_Vertex_Common_Denominator(vertex[0]);
786 Param_Vertex_Common_Denominator(vertex[1]);
788 value_oppose(tmp, vertex[1]->Vertex->p[0][nparam+1]);
789 for (int i = 0; i < 2; ++i)
790 Vector_Combine(vertex[1]->Vertex->p[i],
791 vertex[0]->Vertex->p[i],
792 v0v1->p[i],
793 vertex[0]->Vertex->p[0][nparam+1], tmp, nparam+1);
794 value_multiply(v0v1->p[0][nparam+1],
795 vertex[0]->Vertex->p[0][nparam+1],
796 vertex[1]->Vertex->p[0][nparam+1]);
797 value_assign(v0v1->p[1][nparam+1], v0v1->p[0][nparam+1]);
799 /* Order vertices to ensure we integrate in the right
800 * direction, i.e., with the polytope "on the left".
802 for (int i = 0; i < 2; ++i)
803 Inner_Product(v0v1->p[i], inner, nparam+1, &f_v0v1->p[i]);
805 Inner_Product(dir->p, f_v0v1->p, 2, &tmp);
806 if (value_neg_p(tmp)) {
807 Param_Vertices *PV = vertex[0];
808 vertex[0] = vertex[1];
809 vertex[1] = PV;
810 for (int i = 0; i < 2; ++i)
811 Vector_Oppose(v0v1->p[i], v0v1->p[i], nparam+1);
813 value_oppose(tmp, normal->p[0]);
814 if (value_neg_p(tmp)) {
815 value_oppose(tmp, tmp);
816 Vector_Oppose(v0v1->p[1], v0v1->p[1], nparam+1);
818 value_multiply(tmp, tmp, v0v1->p[1][nparam+1]);
819 evalue *top = affine2evalue(v0v1->p[1], tmp, nparam);
821 value_multiply(s->p[0], normal->p[1], vertex[0]->Vertex->p[0][nparam+1]);
822 Vector_Copy(vertex[0]->Vertex->p[0], s->p+1, nparam+2);
823 subs_1d[0] = affine2evalue(s->p, s->p[1+nparam+1], 1+nparam);
824 value_multiply(s->p[0], normal->p[0], vertex[0]->Vertex->p[0][nparam+1]);
825 value_oppose(s->p[0], s->p[0]);
826 Vector_Copy(vertex[0]->Vertex->p[1], s->p+1, nparam+2);
827 subs_1d[1] = affine2evalue(s->p, s->p[1+nparam+1], 1+nparam);
829 evalue *d = evalue_dup(I);
830 evalue_substitute(d, subs_1d);
831 evalue_anti_derive(d, 0);
833 evalue_free(subs_1d[0]);
834 evalue_free(subs_1d[1]);
836 subs_0d[1] = top;
837 evalue_substitute(d, subs_0d+1);
838 evalue_mul(d, dir->p[1]);
839 evalue_free(subs_0d[1]);
841 eadd(d, res);
842 evalue_free(d);
844 Param_Domain_Free(FD);
845 NEXT(ix, bx);
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 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 #define INT_BITS (sizeof(unsigned) * 8)
968 static unsigned *active_constraints(Param_Polyhedron *PP, Param_Domain *D)
970 int len = (PP->Constraints->NbRows+INT_BITS-1)/INT_BITS;
971 unsigned *facets = (unsigned *)calloc(len, sizeof(unsigned));
972 Param_Vertices *V;
974 FORALL_PVertex_in_ParamPolyhedron(V, D, PP)
975 if (!V->Facets)
976 Param_Vertex_Set_Facets(PP, V);
977 for (int i = 0; i < len; ++i)
978 facets[i] |= V->Facets[i];
979 END_FORALL_PVertex_in_ParamPolyhedron;
981 return facets;
984 static evalue *summate_over_domain(evalue *e, int nvar, Polyhedron *D,
985 struct barvinok_options *options)
987 Polyhedron *U;
988 Param_Polyhedron *PP;
989 evalue *res;
990 int nd = -1;
991 unsigned MaxRays;
992 struct evalue_section *s;
993 Param_Domain *PD;
995 MaxRays = options->MaxRays;
996 POL_UNSET(options->MaxRays, POL_INTEGER);
998 U = Universe_Polyhedron(D->Dimension - nvar);
999 PP = Polyhedron2Param_Polyhedron(D, U, options);
1001 for (nd = 0, PD = PP->D; PD; ++nd, PD = PD->next);
1002 s = ALLOCN(struct evalue_section, nd);
1004 Polyhedron *TC = true_context(D, U, MaxRays);
1005 FORALL_REDUCED_DOMAIN(PP, TC, nd, options, i, PD, rVD)
1006 unsigned *facets;
1008 facets = active_constraints(PP, PD);
1010 Vector *inner = inner_point(rVD);
1011 s[i].D = rVD;
1013 if (nvar == 1) {
1014 s[i].E = summate_over_1d_pdomain(e, PP, PD, inner->p+1, options);
1015 } else if (nvar == 2) {
1016 summator_2d s2d(e, PP, inner->p+1, rVD->Dimension);
1018 s[i].E = s2d.summate_over_pdomain(PP, facets, PD, options);
1021 Vector_Free(inner);
1022 free(facets);
1023 END_FORALL_REDUCED_DOMAIN
1024 Polyhedron_Free(TC);
1025 Polyhedron_Free(U);
1026 Param_Polyhedron_Free(PP);
1028 options->MaxRays = MaxRays;
1030 res = evalue_from_section_array(s, nd);
1031 free(s);
1033 return res;
1036 evalue *euler_summate(evalue *e, unsigned nvar,
1037 struct barvinok_options *options)
1039 int i;
1040 evalue *res;
1042 assert(nvar <= 2);
1044 assert(nvar >= 0);
1045 if (nvar == 0 || EVALUE_IS_ZERO(*e))
1046 return evalue_dup(e);
1048 assert(value_zero_p(e->d));
1049 assert(e->x.p->type == partition);
1051 res = evalue_zero();
1053 for (i = 0; i < e->x.p->size/2; ++i) {
1054 evalue *t;
1055 t = summate_over_domain(&e->x.p->arr[2*i+1], nvar,
1056 EVALUE_DOMAIN(e->x.p->arr[2*i]), options);
1057 eadd(t, res);
1058 evalue_free(t);
1061 return res;