volume.c: drop redundant arguments to volume_simplex
[barvinok.git] / lattice_point.cc
blob397345bb55dab33f701b8068ce2f99bbcc0670b7
1 #include <assert.h>
2 #include <NTL/mat_ZZ.h>
3 #include <NTL/vec_ZZ.h>
4 #include <barvinok/barvinok.h>
5 #include <barvinok/evalue.h>
6 #include <barvinok/util.h>
7 #include "config.h"
8 #include "conversion.h"
9 #include "lattice_point.h"
11 using std::cerr;
12 using std::endl;
14 #define ALLOC(type) (type*)malloc(sizeof(type))
16 /* returns an evalue that corresponds to
18 * c/(*den) x_param
20 static evalue *term(int param, ZZ& c, Value *den = NULL)
22 evalue *EP = new evalue();
23 value_init(EP->d);
24 value_set_si(EP->d,0);
25 EP->x.p = new_enode(polynomial, 2, param + 1);
26 evalue_set_si(&EP->x.p->arr[0], 0, 1);
27 value_init(EP->x.p->arr[1].x.n);
28 if (den == NULL)
29 value_set_si(EP->x.p->arr[1].d, 1);
30 else
31 value_assign(EP->x.p->arr[1].d, *den);
32 zz2value(c, EP->x.p->arr[1].x.n);
33 return EP;
36 /* returns an evalue that corresponds to
38 * sum_i p[i] * x_i
40 evalue *multi_monom(vec_ZZ& p)
42 evalue *X = new evalue();
43 value_init(X->d);
44 value_init(X->x.n);
45 unsigned nparam = p.length()-1;
46 zz2value(p[nparam], X->x.n);
47 value_set_si(X->d, 1);
48 for (int i = 0; i < nparam; ++i) {
49 if (p[i] == 0)
50 continue;
51 evalue *T = term(i, p[i]);
52 eadd(T, X);
53 free_evalue_refs(T);
54 delete T;
56 return X;
60 * Check whether mapping polyhedron P on the affine combination
61 * num yields a range that has a fixed quotient on integer
62 * division by d
63 * If zero is true, then we are only interested in the quotient
64 * for the cases where the remainder is zero.
65 * Returns NULL if false and a newly allocated value if true.
67 static Value *fixed_quotient(Polyhedron *P, vec_ZZ& num, Value d, bool zero)
69 Value* ret = NULL;
70 int len = num.length();
71 Matrix *T = Matrix_Alloc(2, len);
72 zz2values(num, T->p[0]);
73 value_set_si(T->p[1][len-1], 1);
74 Polyhedron *I = Polyhedron_Image(P, T, P->NbConstraints);
75 Matrix_Free(T);
77 int i;
78 for (i = 0; i < I->NbRays; ++i)
79 if (value_zero_p(I->Ray[i][2])) {
80 Polyhedron_Free(I);
81 return NULL;
84 Value min, max;
85 value_init(min);
86 value_init(max);
87 int bounded = line_minmax(I, &min, &max);
88 assert(bounded);
90 if (zero)
91 mpz_cdiv_q(min, min, d);
92 else
93 mpz_fdiv_q(min, min, d);
94 mpz_fdiv_q(max, max, d);
96 if (value_eq(min, max)) {
97 ret = ALLOC(Value);
98 value_init(*ret);
99 value_assign(*ret, min);
101 value_clear(min);
102 value_clear(max);
103 return ret;
107 * Normalize linear expression coef modulo m
108 * Removes common factor and reduces coefficients
109 * Returns index of first non-zero coefficient or len
111 int normal_mod(Value *coef, int len, Value *m)
113 Value gcd;
114 value_init(gcd);
116 Vector_Gcd(coef, len, &gcd);
117 Gcd(gcd, *m, &gcd);
118 Vector_AntiScale(coef, coef, gcd, len);
120 value_division(*m, *m, gcd);
121 value_clear(gcd);
123 if (value_one_p(*m))
124 return len;
126 int j;
127 for (j = 0; j < len; ++j)
128 mpz_fdiv_r(coef[j], coef[j], *m);
129 for (j = 0; j < len; ++j)
130 if (value_notzero_p(coef[j]))
131 break;
133 return j;
136 static bool mod_needed(Polyhedron *PD, vec_ZZ& num, Value d, evalue *E)
138 Value *q = fixed_quotient(PD, num, d, false);
140 if (!q)
141 return true;
143 value_oppose(*q, *q);
144 evalue EV;
145 value_init(EV.d);
146 value_set_si(EV.d, 1);
147 value_init(EV.x.n);
148 value_multiply(EV.x.n, *q, d);
149 eadd(&EV, E);
150 free_evalue_refs(&EV);
151 value_clear(*q);
152 free(q);
153 return false;
156 /* modifies f argument ! */
157 static void ceil_mod(Value *coef, int len, Value d, ZZ& f, evalue *EP, Polyhedron *PD)
159 Value m;
160 value_init(m);
161 value_set_si(m, -1);
163 Vector_Scale(coef, coef, m, len);
165 value_assign(m, d);
166 int j = normal_mod(coef, len, &m);
168 if (j == len) {
169 value_clear(m);
170 return;
173 vec_ZZ num;
174 values2zz(coef, num, len);
176 ZZ g;
177 value2zz(m, g);
179 evalue tmp;
180 value_init(tmp.d);
181 evalue_set_si(&tmp, 0, 1);
183 int p = j;
184 if (g % 2 == 0)
185 while (j < len-1 && (num[j] == g/2 || num[j] == 0))
186 ++j;
187 if ((j < len-1 && num[j] > g/2) || (j == len-1 && num[j] >= (g+1)/2)) {
188 for (int k = j; k < len-1; ++k)
189 if (num[k] != 0)
190 num[k] = g - num[k];
191 num[len-1] = g - 1 - num[len-1];
192 value_assign(tmp.d, m);
193 ZZ t = f*(g-1);
194 zz2value(t, tmp.x.n);
195 eadd(&tmp, EP);
196 f = -f;
199 if (p >= len-1) {
200 ZZ t = num[len-1] * f;
201 zz2value(t, tmp.x.n);
202 value_assign(tmp.d, m);
203 eadd(&tmp, EP);
204 } else {
205 evalue *E = multi_monom(num);
206 evalue EV;
207 value_init(EV.d);
209 if (PD && !mod_needed(PD, num, m, E)) {
210 value_init(EV.x.n);
211 zz2value(f, EV.x.n);
212 value_assign(EV.d, m);
213 emul(&EV, E);
214 eadd(E, EP);
215 } else {
216 value_init(EV.x.n);
217 value_set_si(EV.x.n, 1);
218 value_assign(EV.d, m);
219 emul(&EV, E);
220 value_clear(EV.x.n);
221 value_set_si(EV.d, 0);
222 EV.x.p = new_enode(fractional, 3, -1);
223 evalue_copy(&EV.x.p->arr[0], E);
224 evalue_set_si(&EV.x.p->arr[1], 0, 1);
225 value_init(EV.x.p->arr[2].x.n);
226 zz2value(f, EV.x.p->arr[2].x.n);
227 value_set_si(EV.x.p->arr[2].d, 1);
229 eadd(&EV, EP);
232 free_evalue_refs(&EV);
233 free_evalue_refs(E);
234 delete E;
237 free_evalue_refs(&tmp);
239 out:
240 value_clear(m);
243 static void ceil(Value *coef, int len, Value d, ZZ& f,
244 evalue *EP, Polyhedron *PD, barvinok_options *options)
246 ceil_mod(coef, len, d, f, EP, PD);
247 if (options->lookup_table)
248 evalue_mod2table(EP, len-1);
251 evalue* bv_ceil3(Value *coef, int len, Value d, Polyhedron *P)
253 Vector *val = Vector_Alloc(len);
255 Value t;
256 value_init(t);
257 value_set_si(t, -1);
258 Vector_Scale(coef, val->p, t, len);
259 value_absolute(t, d);
261 vec_ZZ num;
262 values2zz(val->p, num, len);
263 evalue *EP = multi_monom(num);
265 evalue tmp;
266 value_init(tmp.d);
267 value_init(tmp.x.n);
268 value_set_si(tmp.x.n, 1);
269 value_assign(tmp.d, t);
271 emul(&tmp, EP);
273 ZZ one;
274 one = 1;
275 ceil_mod(val->p, len, t, one, EP, P);
276 value_clear(t);
278 /* copy EP to malloc'ed evalue */
279 evalue *E = ALLOC(evalue);
280 *E = *EP;
281 delete EP;
283 free_evalue_refs(&tmp);
284 Vector_Free(val);
286 return E;
289 void lattice_point(Value* values, const mat_ZZ& rays, vec_ZZ& vertex, int *closed)
291 unsigned dim = rays.NumRows();
292 if (value_one_p(values[dim]) && !closed)
293 values2zz(values, vertex, dim);
294 else {
295 Matrix* Rays = rays2matrix(rays);
296 Matrix *inv = Matrix_Alloc(Rays->NbRows, Rays->NbColumns);
297 int ok = Matrix_Inverse(Rays, inv);
298 assert(ok);
299 Matrix_Free(Rays);
300 Rays = rays2matrix(rays);
301 Vector *lambda = Vector_Alloc(dim+1);
302 Vector_Matrix_Product(values, inv, lambda->p);
303 Matrix_Free(inv);
304 for (int j = 0; j < dim; ++j)
305 if (!closed || closed[j])
306 mpz_cdiv_q(lambda->p[j], lambda->p[j], lambda->p[dim]);
307 else {
308 value_addto(lambda->p[j], lambda->p[j], lambda->p[dim]);
309 mpz_fdiv_q(lambda->p[j], lambda->p[j], lambda->p[dim]);
311 value_set_si(lambda->p[dim], 1);
312 Vector *A = Vector_Alloc(dim+1);
313 Vector_Matrix_Product(lambda->p, Rays, A->p);
314 Vector_Free(lambda);
315 Matrix_Free(Rays);
316 values2zz(A->p, vertex, dim);
317 Vector_Free(A);
321 /* Compute the lattice points in the vertex cone at "values" with rays "rays".
322 * The lattice points are returned in "vertex".
324 * Rays has the generators as rows and so does W.
325 * We first compute { m-v, u_i^* } with m = k W, where k runs through
326 * the cosets.
327 * We compute
328 * [k 1] [ d1*W 0 ] [ U' 0 ] = [k 1] T2
329 * [ -v d1 ] [ 0 d2 ]
330 * where d1 and d2 are the denominators of v and U^{-1}=U'/d2.
331 * Then lambda = { k } (componentwise)
332 * We compute x - floor(x) = {x} = { a/b } as fdiv_r(a,b)/b
333 * For open rays/facets, we need values in (0,1] rather than [0,1),
334 * so we compute {{x}} = x - ceil(x-1) = a/b - ceil((a-b)/b)
335 * = (a - b cdiv_q(a-b,b) - b + b)/b
336 * = (cdiv_r(a,b)+b)/b
337 * Finally, we compute v + lambda * U
338 * The denominator of lambda can be d1*d2, that of lambda2 = lambda*U
339 * can be at most d1, since it is integer if v = 0.
340 * The denominator of v + lambda2 is 1.
342 void lattice_point(Value* values, const mat_ZZ& rays, mat_ZZ& vertex,
343 unsigned long det, int *closed)
345 unsigned dim = rays.NumRows();
346 vertex.SetDims(det, dim);
347 if (det == 1) {
348 lattice_point(values, rays, vertex[0], closed);
349 return;
351 Matrix* Rays = rays2matrix2(rays);
352 Matrix *U, *W, *D;
353 Smith(Rays, &U, &W, &D);
354 Matrix_Free(Rays);
355 Matrix_Free(U);
357 /* Sanity check */
358 unsigned long det2 = 1;
359 for (int i = 0 ; i < D->NbRows; ++i)
360 det2 *= mpz_get_ui(D->p[i][i]);
361 assert(det == det2);
363 Matrix *T = Matrix_Alloc(W->NbRows+1, W->NbColumns+1);
364 for (int i = 0; i < W->NbRows; ++i)
365 Vector_Scale(W->p[i], T->p[i], values[dim], W->NbColumns);
366 Matrix_Free(W);
367 Value tmp;
368 value_init(tmp);
369 value_set_si(tmp, -1);
370 Vector_Scale(values, T->p[dim], tmp, dim);
371 value_clear(tmp);
372 value_assign(T->p[dim][dim], values[dim]);
374 Rays = rays2matrix(rays);
375 Matrix *inv = Matrix_Alloc(Rays->NbRows, Rays->NbColumns);
376 int ok = Matrix_Inverse(Rays, inv);
377 assert(ok);
378 Matrix_Free(Rays);
380 Matrix *T2 = Matrix_Alloc(dim+1, dim+1);
381 Matrix_Product(T, inv, T2);
382 Matrix_Free(T);
384 Rays = rays2matrix(rays);
386 Vector *k = Vector_Alloc(dim+1);
387 value_set_si(k->p[dim], 1);
388 Vector *lambda = Vector_Alloc(dim+1);
389 Vector *lambda2 = Vector_Alloc(dim+1);
390 for (unsigned long i = 0; i < det; ++i) {
391 unsigned long val = i;
392 for (int j = 0; j < dim; ++j) {
393 value_set_si(k->p[j], val % mpz_get_ui(D->p[j][j]));
394 val /= mpz_get_ui(D->p[j][j]);
396 Vector_Matrix_Product(k->p, T2, lambda->p);
397 for (int j = 0; j < dim; ++j)
398 if (!closed || closed[j])
399 mpz_fdiv_r(lambda->p[j], lambda->p[j], lambda->p[dim]);
400 else {
401 mpz_cdiv_r(lambda->p[j], lambda->p[j], lambda->p[dim]);
402 value_addto(lambda->p[j], lambda->p[j], lambda->p[dim]);
404 Vector_Matrix_Product(lambda->p, Rays, lambda2->p);
405 for (int j = 0; j < dim; ++j)
406 assert(mpz_divisible_p(lambda2->p[j], inv->p[dim][dim]));
407 Vector_AntiScale(lambda2->p, lambda2->p, inv->p[dim][dim], dim+1);
408 Vector_Add(lambda2->p, values, lambda2->p, dim);
409 for (int j = 0; j < dim; ++j)
410 assert(mpz_divisible_p(lambda2->p[j], values[dim]));
411 Vector_AntiScale(lambda2->p, lambda2->p, values[dim], dim+1);
412 values2zz(lambda2->p, vertex[i], dim);
414 Vector_Free(k);
415 Vector_Free(lambda);
416 Vector_Free(lambda2);
417 Matrix_Free(D);
418 Matrix_Free(Rays);
419 Matrix_Free(inv);
421 Matrix_Free(T2);
424 static void vertex_period(
425 const mat_ZZ& rays, vec_ZZ& lambda, Matrix *T,
426 Value lcm, int p, Vector *val,
427 evalue *E, evalue* ev,
428 ZZ& offset)
430 unsigned nparam = T->NbRows - 1;
431 unsigned dim = rays.NumRows();
432 Value tmp;
433 ZZ nump;
435 if (p == nparam) {
436 vec_ZZ vertex;
437 ZZ num, l;
438 Vector * values = Vector_Alloc(dim + 1);
439 Vector_Matrix_Product(val->p, T, values->p);
440 value_assign(values->p[dim], lcm);
441 lattice_point(values->p, rays, vertex, NULL);
442 num = vertex * lambda;
443 value2zz(lcm, l);
444 num *= l;
445 num += offset;
446 value_init(ev->x.n);
447 zz2value(num, ev->x.n);
448 value_assign(ev->d, lcm);
449 Vector_Free(values);
450 return;
453 value_init(tmp);
454 vec_ZZ vertex;
455 values2zz(T->p[p], vertex, dim);
456 nump = vertex * lambda;
457 if (First_Non_Zero(val->p, p) == -1) {
458 value_assign(tmp, lcm);
459 evalue *ET = term(p, nump, &tmp);
460 eadd(ET, E);
461 free_evalue_refs(ET);
462 delete ET;
465 value_assign(tmp, lcm);
466 if (First_Non_Zero(T->p[p], dim) != -1)
467 Vector_Gcd(T->p[p], dim, &tmp);
468 Gcd(tmp, lcm, &tmp);
469 if (value_lt(tmp, lcm)) {
470 ZZ count;
472 value_division(tmp, lcm, tmp);
473 value_set_si(ev->d, 0);
474 ev->x.p = new_enode(periodic, VALUE_TO_INT(tmp), p+1);
475 value2zz(tmp, count);
476 do {
477 value_decrement(tmp, tmp);
478 --count;
479 ZZ new_offset = offset - count * nump;
480 value_assign(val->p[p], tmp);
481 vertex_period(rays, lambda, T, lcm, p+1, val, E,
482 &ev->x.p->arr[VALUE_TO_INT(tmp)], new_offset);
483 } while (value_pos_p(tmp));
484 } else
485 vertex_period(rays, lambda, T, lcm, p+1, val, E, ev, offset);
486 value_clear(tmp);
489 /* Returns the power of (t+1) in the term of a rational generating function,
490 * i.e., the scalar product of the actual lattice point and lambda.
491 * The lattice point is the unique lattice point in the fundamental parallelepiped
492 * of the unimodual cone i shifted to the parametric vertex W/lcm.
494 * The rows of W refer to the coordinates of the vertex
495 * The first nparam columns are the coefficients of the parameters
496 * and the final column is the constant term.
497 * lcm is the common denominator of all coefficients.
499 * PD is the parameter domain, which, if != NULL, may be used to simply the
500 * resulting expression.
502 static evalue* lattice_point_fractional(const mat_ZZ& rays, vec_ZZ& lambda,
503 Matrix *W, Value lcm, Polyhedron *PD)
505 unsigned nparam = W->NbColumns - 1;
507 Matrix* Rays = rays2matrix2(rays);
508 Matrix *T = Transpose(Rays);
509 Matrix *T2 = Matrix_Copy(T);
510 Matrix *inv = Matrix_Alloc(T2->NbRows, T2->NbColumns);
511 int ok = Matrix_Inverse(T2, inv);
512 assert(ok);
513 Matrix_Free(Rays);
514 Matrix_Free(T2);
515 mat_ZZ vertex;
516 matrix2zz(W, vertex, W->NbRows, W->NbColumns);
518 vec_ZZ num;
519 num = lambda * vertex;
521 evalue *EP = multi_monom(num);
523 evalue tmp;
524 value_init(tmp.d);
525 value_init(tmp.x.n);
526 value_set_si(tmp.x.n, 1);
527 value_assign(tmp.d, lcm);
529 emul(&tmp, EP);
531 Matrix *L = Matrix_Alloc(inv->NbRows, W->NbColumns);
532 Matrix_Product(inv, W, L);
534 mat_ZZ RT;
535 matrix2zz(T, RT, T->NbRows, T->NbColumns);
536 Matrix_Free(T);
538 vec_ZZ p = lambda * RT;
540 for (int i = 0; i < L->NbRows; ++i) {
541 ceil_mod(L->p[i], nparam+1, lcm, p[i], EP, PD);
544 Matrix_Free(L);
546 Matrix_Free(inv);
547 free_evalue_refs(&tmp);
548 return EP;
551 static evalue* lattice_point_table(const mat_ZZ& rays, vec_ZZ& lambda, Matrix *W,
552 Value lcm, Polyhedron *PD)
554 Matrix *T = Transpose(W);
555 unsigned nparam = T->NbRows - 1;
557 evalue *EP = new evalue();
558 value_init(EP->d);
559 evalue_set_si(EP, 0, 1);
561 evalue ev;
562 Vector *val = Vector_Alloc(nparam+1);
563 value_set_si(val->p[nparam], 1);
564 ZZ offset(INIT_VAL, 0);
565 value_init(ev.d);
566 vertex_period(rays, lambda, T, lcm, 0, val, EP, &ev, offset);
567 Vector_Free(val);
568 eadd(&ev, EP);
569 free_evalue_refs(&ev);
571 Matrix_Free(T);
573 reduce_evalue(EP);
575 return EP;
578 evalue* lattice_point(const mat_ZZ& rays, vec_ZZ& lambda, Matrix *W,
579 Value lcm, Polyhedron *PD, barvinok_options *options)
581 if (options->lookup_table)
582 return lattice_point_table(rays, lambda, W, lcm, PD);
583 else
584 return lattice_point_fractional(rays, lambda, W, lcm, PD);
587 /* returns the unique lattice point in the fundamental parallelepiped
588 * of the unimodual cone C shifted to the parametric vertex V.
590 * The return values num and E_vertex are such that
591 * coordinate i of this lattice point is equal to
593 * num[i] + E_vertex[i]
595 void lattice_point(Param_Vertices *V, const mat_ZZ& rays, vec_ZZ& num,
596 evalue **E_vertex, barvinok_options *options)
598 unsigned nparam = V->Vertex->NbColumns - 2;
599 unsigned dim = rays.NumCols();
600 vec_ZZ vertex;
601 vertex.SetLength(nparam+1);
603 Value lcm, tmp;
604 value_init(lcm);
605 value_init(tmp);
606 value_set_si(lcm, 1);
608 for (int j = 0; j < V->Vertex->NbRows; ++j) {
609 value_lcm(lcm, V->Vertex->p[j][nparam+1], &lcm);
612 if (value_notone_p(lcm)) {
613 Matrix * mv = Matrix_Alloc(dim, nparam+1);
614 for (int j = 0 ; j < dim; ++j) {
615 value_division(tmp, lcm, V->Vertex->p[j][nparam+1]);
616 Vector_Scale(V->Vertex->p[j], mv->p[j], tmp, nparam+1);
619 Matrix* Rays = rays2matrix2(rays);
620 Matrix *T = Transpose(Rays);
621 Matrix *T2 = Matrix_Copy(T);
622 Matrix *inv = Matrix_Alloc(T2->NbRows, T2->NbColumns);
623 int ok = Matrix_Inverse(T2, inv);
624 assert(ok);
625 Matrix_Free(Rays);
626 Matrix_Free(T2);
627 Matrix *L = Matrix_Alloc(inv->NbRows, mv->NbColumns);
628 Matrix_Product(inv, mv, L);
629 Matrix_Free(inv);
631 evalue f;
632 value_init(f.d);
633 value_init(f.x.n);
635 ZZ one;
637 evalue *remainders[dim];
638 for (int i = 0; i < dim; ++i) {
639 remainders[i] = evalue_zero();
640 one = 1;
641 ceil(L->p[i], nparam+1, lcm, one, remainders[i], 0, options);
643 Matrix_Free(L);
646 for (int i = 0; i < V->Vertex->NbRows; ++i) {
647 values2zz(mv->p[i], vertex, nparam+1);
648 E_vertex[i] = multi_monom(vertex);
649 num[i] = 0;
651 value_set_si(f.x.n, 1);
652 value_assign(f.d, lcm);
654 emul(&f, E_vertex[i]);
656 for (int j = 0; j < dim; ++j) {
657 if (value_zero_p(T->p[i][j]))
658 continue;
659 evalue cp;
660 value_init(cp.d);
661 evalue_copy(&cp, remainders[j]);
662 if (value_notone_p(T->p[i][j])) {
663 value_set_si(f.d, 1);
664 value_assign(f.x.n, T->p[i][j]);
665 emul(&f, &cp);
667 eadd(&cp, E_vertex[i]);
668 free_evalue_refs(&cp);
671 for (int i = 0; i < dim; ++i) {
672 free_evalue_refs(remainders[i]);
673 free(remainders[i]);
676 free_evalue_refs(&f);
678 Matrix_Free(T);
679 Matrix_Free(mv);
680 value_clear(lcm);
681 value_clear(tmp);
682 return;
684 value_clear(lcm);
685 value_clear(tmp);
687 for (int i = 0; i < V->Vertex->NbRows; ++i) {
688 /* fixed value */
689 if (First_Non_Zero(V->Vertex->p[i], nparam) == -1) {
690 E_vertex[i] = 0;
691 value2zz(V->Vertex->p[i][nparam], num[i]);
692 } else {
693 values2zz(V->Vertex->p[i], vertex, nparam+1);
694 E_vertex[i] = multi_monom(vertex);
695 num[i] = 0;