3 #include "isl_map_private.h"
4 #include "isl_equalities.h"
6 /* Given a set of modulo constraints
10 * this function computes a particular solution y_0
12 * The input is given as a matrix B = [ c A ] and a vector d.
14 * The output is matrix containing the solution y_0 or
15 * a zero-column matrix if the constraints admit no integer solution.
17 * The given set of constrains is equivalent to
21 * with D = diag d and x a fresh set of variables.
22 * Reducing both c and A modulo d does not change the
23 * value of y in the solution and may lead to smaller coefficients.
24 * Let M = [ D A ] and [ H 0 ] = M U, the Hermite normal form of M.
30 * [ H 0 ] U^{-1} [ y ] = - c
33 * [ B ] = U^{-1} [ y ]
37 * so B may be chosen arbitrarily, e.g., B = 0, and then
40 * U^{-1} [ y ] = [ 0 ]
48 * If any of the coordinates of this y are non-integer
49 * then the constraints admit no integer solution and
50 * a zero-column matrix is returned.
52 static struct isl_mat *particular_solution(struct isl_ctx *ctx,
53 struct isl_mat *B, struct isl_vec *d)
56 struct isl_mat *M = NULL;
57 struct isl_mat *C = NULL;
58 struct isl_mat *U = NULL;
59 struct isl_mat *H = NULL;
60 struct isl_mat *cst = NULL;
61 struct isl_mat *T = NULL;
63 M = isl_mat_alloc(ctx, B->n_row, B->n_row + B->n_col - 1);
64 C = isl_mat_alloc(ctx, 1 + B->n_row, 1);
67 isl_int_set_si(C->row[0][0], 1);
68 for (i = 0; i < B->n_row; ++i) {
69 isl_seq_clr(M->row[i], B->n_row);
70 isl_int_set(M->row[i][i], d->block.data[i]);
71 isl_int_neg(C->row[1 + i][0], B->row[i][0]);
72 isl_int_fdiv_r(C->row[1+i][0], C->row[1+i][0], M->row[i][i]);
73 for (j = 0; j < B->n_col - 1; ++j)
74 isl_int_fdiv_r(M->row[i][B->n_row + j],
75 B->row[i][1 + j], M->row[i][i]);
77 M = isl_mat_left_hermite(ctx, M, 0, &U, NULL);
80 H = isl_mat_sub_alloc(ctx, M->row, 0, B->n_row, 0, B->n_row);
81 H = isl_mat_lin_to_aff(ctx, H);
82 C = isl_mat_inverse_product(ctx, H, C);
85 for (i = 0; i < B->n_row; ++i) {
86 if (!isl_int_is_divisible_by(C->row[1+i][0], C->row[0][0]))
88 isl_int_divexact(C->row[1+i][0], C->row[1+i][0], C->row[0][0]);
91 cst = isl_mat_alloc(ctx, B->n_row, 0);
93 cst = isl_mat_sub_alloc(ctx, C->row, 1, B->n_row, 0, 1);
94 T = isl_mat_sub_alloc(ctx, U->row, B->n_row, B->n_col - 1, 0, B->n_row);
95 cst = isl_mat_product(ctx, T, cst);
101 isl_mat_free(ctx, M);
102 isl_mat_free(ctx, C);
103 isl_mat_free(ctx, U);
107 /* Given a set of modulo constraints
111 * this function returns an affine transformation T,
115 * that bijectively maps the integer vectors y' to integer
116 * vectors y that satisfy the modulo constraints.
118 * The implementation closely follows the description in Section 2.5.3
119 * of B. Meister, "Stating and Manipulating Periodicity in the Polytope
120 * Model. Applications to Program Analysis and Optimization".
122 * The input is given as a matrix B = [ c A ] and a vector d.
123 * Each element of the vector d corresponds to a row in B.
124 * The output is a lower triangular matrix.
125 * If no integer vector y satisfies the given constraints then
126 * a matrix with zero columns is returned.
128 * We first compute a particular solution y_0 to the given set of
129 * modulo constraints in particular_solution. If no such solution
130 * exists, then we return a zero-columned transformation matrix.
131 * Otherwise, we compute the generic solution to
135 * Let K be the right kernel of A, then any y = K y'' is a solution.
136 * Any multiple of a unit vector s_j e_j such that for each row i,
137 * we have that s_j a_{i,j} is a multiple of d_i, is also a generator
138 * for the set of solutions. The smallest such s_j can be obtained
139 * by taking the lcm of all the a_{i,j}/gcd(a_{i,j},d).
140 * That is, y = K y'' + S y''', with S = diag s is the general solution.
141 * To obtain a minimal representation we compute the Hermite normal
142 * form [ G 0 ] of [ S K ].
143 * The affine transformation matrix returned is then
148 * as any y = y_0 + G y' with y' integer is a solution to the original
149 * modulo constraints.
151 struct isl_mat *isl_mat_parameter_compression(struct isl_ctx *ctx,
152 struct isl_mat *B, struct isl_vec *d)
155 struct isl_mat *cst = NULL;
156 struct isl_mat *K = NULL;
157 struct isl_mat *M = NULL;
163 isl_assert(ctx, B->n_row == d->size, goto error);
164 cst = particular_solution(ctx, B, d);
167 if (cst->n_col == 0) {
168 T = isl_mat_alloc(ctx, B->n_col, 0);
169 isl_mat_free(ctx, cst);
170 isl_mat_free(ctx, B);
171 isl_vec_free(ctx, d);
174 T = isl_mat_sub_alloc(ctx, B->row, 0, B->n_row, 1, B->n_col - 1);
175 K = isl_mat_right_kernel(ctx, T);
178 M = isl_mat_alloc(ctx, B->n_col - 1, B->n_col - 1 + K->n_col);
181 isl_mat_sub_copy(ctx, M->row, K->row, M->n_row,
182 B->n_col - 1, 0, K->n_col);
183 for (i = 0; i < M->n_row; ++i) {
184 isl_seq_clr(M->row[i], B->n_col - 1);
185 isl_int_set_si(M->row[i][i], 1);
188 for (i = 0; i < B->n_row; ++i)
189 for (j = 0; j < B->n_col - 1; ++j) {
190 isl_int_gcd(v, B->row[i][1+j], d->block.data[i]);
191 isl_int_divexact(v, d->block.data[i], v);
192 isl_int_lcm(M->row[j][j], M->row[j][j], v);
195 M = isl_mat_left_hermite(ctx, M, 0, NULL, NULL);
198 T = isl_mat_alloc(ctx, 1 + B->n_col - 1, B->n_col);
201 isl_int_set_si(T->row[0][0], 1);
202 isl_seq_clr(T->row[0] + 1, T->n_col - 1);
203 isl_mat_sub_copy(ctx, T->row + 1, cst->row, cst->n_row, 0, 0, 1);
204 isl_mat_sub_copy(ctx, T->row + 1, M->row, M->n_row, 1, 0, T->n_col - 1);
205 isl_mat_free(ctx, K);
206 isl_mat_free(ctx, M);
207 isl_mat_free(ctx, cst);
208 isl_mat_free(ctx, B);
209 isl_vec_free(ctx, d);
212 isl_mat_free(ctx, K);
213 isl_mat_free(ctx, M);
214 isl_mat_free(ctx, cst);
215 isl_mat_free(ctx, B);
216 isl_vec_free(ctx, d);
220 /* Given a set of equalities
224 * this function computes unimodular transformation from a lower-dimensional
225 * space to the original space that bijectively maps the integer points x'
226 * in the lower-dimensional space to the integer points x in the original
227 * space that satisfy the equalities.
229 * The input is given as a matrix B = [ -c M ] and the out is a
230 * matrix that maps [1 x'] to [1 x].
231 * If T2 is not NULL, then *T2 is set to a matrix mapping [1 x] to [1 x'].
233 * First compute the (left) Hermite normal form of M,
235 * M [U1 U2] = M U = H = [H1 0]
237 * M = H Q = [H1 0] [Q1]
240 * with U, Q unimodular, Q = U^{-1} (and H lower triangular).
241 * Define the transformed variables as
243 * x = [U1 U2] [ x1' ] = [U1 U2] [Q1] x
246 * The equalities then become
248 * H1 x1' - c = 0 or x1' = H1^{-1} c = c'
250 * If any of the c' is non-integer, then the original set has no
251 * integer solutions (since the x' are a unimodular transformation
253 * Otherwise, the transformation is given by
255 * x = U1 H1^{-1} c + U2 x2'
257 * The inverse transformation is simply
261 struct isl_mat *isl_mat_variable_compression(struct isl_ctx *ctx,
262 struct isl_mat *B, struct isl_mat **T2)
265 struct isl_mat *H = NULL, *C = NULL, *H1, *U = NULL, *U1, *U2, *TC;
274 H = isl_mat_sub_alloc(ctx, B->row, 0, B->n_row, 1, dim);
275 H = isl_mat_left_hermite(ctx, H, 0, &U, T2);
276 if (!H || !U || (T2 && !*T2))
279 *T2 = isl_mat_drop_rows(ctx, *T2, 0, B->n_row);
280 *T2 = isl_mat_lin_to_aff(ctx, *T2);
284 C = isl_mat_alloc(ctx, 1+B->n_row, 1);
287 isl_int_set_si(C->row[0][0], 1);
288 isl_mat_sub_neg(ctx, C->row+1, B->row, B->n_row, 0, 0, 1);
289 H1 = isl_mat_sub_alloc(ctx, H->row, 0, H->n_row, 0, H->n_row);
290 H1 = isl_mat_lin_to_aff(ctx, H1);
291 TC = isl_mat_inverse_product(ctx, H1, C);
294 isl_mat_free(ctx, H);
295 if (!isl_int_is_one(TC->row[0][0])) {
296 for (i = 0; i < B->n_row; ++i) {
297 if (!isl_int_is_divisible_by(TC->row[1+i][0], TC->row[0][0])) {
298 isl_mat_free(ctx, B);
299 isl_mat_free(ctx, TC);
300 isl_mat_free(ctx, U);
302 isl_mat_free(ctx, *T2);
305 return isl_mat_alloc(ctx, 1 + dim, 0);
307 isl_seq_scale_down(TC->row[1+i], TC->row[1+i], TC->row[0][0], 1);
309 isl_int_set_si(TC->row[0][0], 1);
311 U1 = isl_mat_sub_alloc(ctx, U->row, 0, U->n_row, 0, B->n_row);
312 U1 = isl_mat_lin_to_aff(ctx, U1);
313 U2 = isl_mat_sub_alloc(ctx, U->row, 0, U->n_row,
314 B->n_row, U->n_row - B->n_row);
315 U2 = isl_mat_lin_to_aff(ctx, U2);
316 isl_mat_free(ctx, U);
317 TC = isl_mat_product(ctx, U1, TC);
318 TC = isl_mat_aff_direct_sum(ctx, TC, U2);
320 isl_mat_free(ctx, B);
324 isl_mat_free(ctx, B);
325 isl_mat_free(ctx, H);
326 isl_mat_free(ctx, U);
328 isl_mat_free(ctx, *T2);
334 /* Use the n equalities of bset to unimodularly transform the
335 * variables x such that n transformed variables x1' have a constant value
336 * and rewrite the constraints of bset in terms of the remaining
337 * transformed variables x2'. The matrix pointed to by T maps
338 * the new variables x2' back to the original variables x, while T2
339 * maps the original variables to the new variables.
341 static struct isl_basic_set *compress_variables(struct isl_ctx *ctx,
342 struct isl_basic_set *bset, struct isl_mat **T, struct isl_mat **T2)
344 struct isl_mat *B, *TC;
353 isl_assert(ctx, isl_basic_set_n_param(bset) == 0, goto error);
354 isl_assert(ctx, bset->n_div == 0, goto error);
355 dim = isl_basic_set_n_dim(bset);
356 isl_assert(ctx, bset->n_eq <= dim, goto error);
360 B = isl_mat_sub_alloc(ctx, bset->eq, 0, bset->n_eq, 0, 1 + dim);
361 TC = isl_mat_variable_compression(ctx, B, T2);
364 if (TC->n_col == 0) {
365 isl_mat_free(ctx, TC);
367 isl_mat_free(ctx, *T2);
370 return isl_basic_set_set_to_empty(bset);
373 bset = isl_basic_set_preimage(ctx, bset, T ? isl_mat_copy(ctx, TC) : TC);
378 isl_basic_set_free(bset);
382 struct isl_basic_set *isl_basic_set_remove_equalities(
383 struct isl_basic_set *bset, struct isl_mat **T, struct isl_mat **T2)
391 isl_assert(bset->ctx, isl_basic_set_n_param(bset) == 0, goto error);
392 bset = isl_basic_set_gauss(bset, NULL);
393 if (F_ISSET(bset, ISL_BASIC_SET_EMPTY))
395 bset = compress_variables(bset->ctx, bset, T, T2);
398 isl_basic_set_free(bset);
403 /* Check if dimension dim belongs to a residue class
404 * i_dim \equiv r mod m
405 * with m != 1 and if so return m in *modulo and r in *residue.
407 int isl_basic_set_dim_residue_class(struct isl_basic_set *bset,
408 int pos, isl_int *modulo, isl_int *residue)
411 struct isl_mat *H = NULL, *U = NULL, *C, *H1, *U1;
415 if (!bset || !modulo || !residue)
419 total = isl_basic_set_total_dim(bset);
420 nparam = isl_basic_set_n_param(bset);
421 H = isl_mat_sub_alloc(ctx, bset->eq, 0, bset->n_eq, 1, total);
422 H = isl_mat_left_hermite(ctx, H, 0, &U, NULL);
426 isl_seq_gcd(U->row[nparam + pos]+bset->n_eq,
427 total-bset->n_eq, modulo);
428 if (isl_int_is_zero(*modulo) || isl_int_is_one(*modulo)) {
429 isl_int_set_si(*residue, 0);
430 isl_mat_free(ctx, H);
431 isl_mat_free(ctx, U);
435 C = isl_mat_alloc(ctx, 1+bset->n_eq, 1);
438 isl_int_set_si(C->row[0][0], 1);
439 isl_mat_sub_neg(ctx, C->row+1, bset->eq, bset->n_eq, 0, 0, 1);
440 H1 = isl_mat_sub_alloc(ctx, H->row, 0, H->n_row, 0, H->n_row);
441 H1 = isl_mat_lin_to_aff(ctx, H1);
442 C = isl_mat_inverse_product(ctx, H1, C);
443 isl_mat_free(ctx, H);
444 U1 = isl_mat_sub_alloc(ctx, U->row, nparam+pos, 1, 0, bset->n_eq);
445 U1 = isl_mat_lin_to_aff(ctx, U1);
446 isl_mat_free(ctx, U);
447 C = isl_mat_product(ctx, U1, C);
450 if (!isl_int_is_divisible_by(C->row[1][0], C->row[0][0])) {
451 bset = isl_basic_set_copy(bset);
452 bset = isl_basic_set_set_to_empty(bset);
453 isl_basic_set_free(bset);
454 isl_int_set_si(*modulo, 0);
455 isl_int_set_si(*residue, 0);
458 isl_int_divexact(*residue, C->row[1][0], C->row[0][0]);
459 isl_int_fdiv_r(*residue, *residue, *modulo);
460 isl_mat_free(ctx, C);
463 isl_mat_free(ctx, H);
464 isl_mat_free(ctx, U);