isl_basic_set_opt: avoid invalid access on error path
[platform/upstream/isl.git] / isl_equalities.c
1 /*
2  * Copyright 2008-2009 Katholieke Universiteit Leuven
3  *
4  * Use of this software is governed by the MIT license
5  *
6  * Written by Sven Verdoolaege, K.U.Leuven, Departement
7  * Computerwetenschappen, Celestijnenlaan 200A, B-3001 Leuven, Belgium
8  */
9
10 #include <isl_mat_private.h>
11 #include <isl/seq.h>
12 #include "isl_map_private.h"
13 #include "isl_equalities.h"
14
15 /* Given a set of modulo constraints
16  *
17  *              c + A y = 0 mod d
18  *
19  * this function computes a particular solution y_0
20  *
21  * The input is given as a matrix B = [ c A ] and a vector d.
22  *
23  * The output is matrix containing the solution y_0 or
24  * a zero-column matrix if the constraints admit no integer solution.
25  *
26  * The given set of constrains is equivalent to
27  *
28  *              c + A y = -D x
29  *
30  * with D = diag d and x a fresh set of variables.
31  * Reducing both c and A modulo d does not change the
32  * value of y in the solution and may lead to smaller coefficients.
33  * Let M = [ D A ] and [ H 0 ] = M U, the Hermite normal form of M.
34  * Then
35  *                [ x ]
36  *              M [ y ] = - c
37  * and so
38  *                             [ x ]
39  *              [ H 0 ] U^{-1} [ y ] = - c
40  * Let
41  *              [ A ]          [ x ]
42  *              [ B ] = U^{-1} [ y ]
43  * then
44  *              H A + 0 B = -c
45  *
46  * so B may be chosen arbitrarily, e.g., B = 0, and then
47  *
48  *                     [ x ] = [ -c ]
49  *              U^{-1} [ y ] = [  0 ]
50  * or
51  *              [ x ]     [ -c ]
52  *              [ y ] = U [  0 ]
53  * specifically,
54  *
55  *              y = U_{2,1} (-c)
56  *
57  * If any of the coordinates of this y are non-integer
58  * then the constraints admit no integer solution and
59  * a zero-column matrix is returned.
60  */
61 static struct isl_mat *particular_solution(struct isl_mat *B, struct isl_vec *d)
62 {
63         int i, j;
64         struct isl_mat *M = NULL;
65         struct isl_mat *C = NULL;
66         struct isl_mat *U = NULL;
67         struct isl_mat *H = NULL;
68         struct isl_mat *cst = NULL;
69         struct isl_mat *T = NULL;
70
71         M = isl_mat_alloc(B->ctx, B->n_row, B->n_row + B->n_col - 1);
72         C = isl_mat_alloc(B->ctx, 1 + B->n_row, 1);
73         if (!M || !C)
74                 goto error;
75         isl_int_set_si(C->row[0][0], 1);
76         for (i = 0; i < B->n_row; ++i) {
77                 isl_seq_clr(M->row[i], B->n_row);
78                 isl_int_set(M->row[i][i], d->block.data[i]);
79                 isl_int_neg(C->row[1 + i][0], B->row[i][0]);
80                 isl_int_fdiv_r(C->row[1+i][0], C->row[1+i][0], M->row[i][i]);
81                 for (j = 0; j < B->n_col - 1; ++j)
82                         isl_int_fdiv_r(M->row[i][B->n_row + j],
83                                         B->row[i][1 + j], M->row[i][i]);
84         }
85         M = isl_mat_left_hermite(M, 0, &U, NULL);
86         if (!M || !U)
87                 goto error;
88         H = isl_mat_sub_alloc(M, 0, B->n_row, 0, B->n_row);
89         H = isl_mat_lin_to_aff(H);
90         C = isl_mat_inverse_product(H, C);
91         if (!C)
92                 goto error;
93         for (i = 0; i < B->n_row; ++i) {
94                 if (!isl_int_is_divisible_by(C->row[1+i][0], C->row[0][0]))
95                         break;
96                 isl_int_divexact(C->row[1+i][0], C->row[1+i][0], C->row[0][0]);
97         }
98         if (i < B->n_row)
99                 cst = isl_mat_alloc(B->ctx, B->n_row, 0);
100         else
101                 cst = isl_mat_sub_alloc(C, 1, B->n_row, 0, 1);
102         T = isl_mat_sub_alloc(U, B->n_row, B->n_col - 1, 0, B->n_row);
103         cst = isl_mat_product(T, cst);
104         isl_mat_free(M);
105         isl_mat_free(C);
106         isl_mat_free(U);
107         return cst;
108 error:
109         isl_mat_free(M);
110         isl_mat_free(C);
111         isl_mat_free(U);
112         return NULL;
113 }
114
115 /* Compute and return the matrix
116  *
117  *              U_1^{-1} diag(d_1, 1, ..., 1)
118  *
119  * with U_1 the unimodular completion of the first (and only) row of B.
120  * The columns of this matrix generate the lattice that satisfies
121  * the single (linear) modulo constraint.
122  */
123 static struct isl_mat *parameter_compression_1(
124                         struct isl_mat *B, struct isl_vec *d)
125 {
126         struct isl_mat *U;
127
128         U = isl_mat_alloc(B->ctx, B->n_col - 1, B->n_col - 1);
129         if (!U)
130                 return NULL;
131         isl_seq_cpy(U->row[0], B->row[0] + 1, B->n_col - 1);
132         U = isl_mat_unimodular_complete(U, 1);
133         U = isl_mat_right_inverse(U);
134         if (!U)
135                 return NULL;
136         isl_mat_col_mul(U, 0, d->block.data[0], 0);
137         U = isl_mat_lin_to_aff(U);
138         return U;
139 }
140
141 /* Compute a common lattice of solutions to the linear modulo
142  * constraints specified by B and d.
143  * See also the documentation of isl_mat_parameter_compression.
144  * We put the matrix
145  * 
146  *              A = [ L_1^{-T} L_2^{-T} ... L_k^{-T} ]
147  *
148  * on a common denominator.  This denominator D is the lcm of modulos d.
149  * Since L_i = U_i^{-1} diag(d_i, 1, ... 1), we have
150  * L_i^{-T} = U_i^T diag(d_i, 1, ... 1)^{-T} = U_i^T diag(1/d_i, 1, ..., 1).
151  * Putting this on the common denominator, we have
152  * D * L_i^{-T} = U_i^T diag(D/d_i, D, ..., D).
153  */
154 static struct isl_mat *parameter_compression_multi(
155                         struct isl_mat *B, struct isl_vec *d)
156 {
157         int i, j, k;
158         isl_int D;
159         struct isl_mat *A = NULL, *U = NULL;
160         struct isl_mat *T;
161         unsigned size;
162
163         isl_int_init(D);
164
165         isl_vec_lcm(d, &D);
166
167         size = B->n_col - 1;
168         A = isl_mat_alloc(B->ctx, size, B->n_row * size);
169         U = isl_mat_alloc(B->ctx, size, size);
170         if (!U || !A)
171                 goto error;
172         for (i = 0; i < B->n_row; ++i) {
173                 isl_seq_cpy(U->row[0], B->row[i] + 1, size);
174                 U = isl_mat_unimodular_complete(U, 1);
175                 if (!U)
176                         goto error;
177                 isl_int_divexact(D, D, d->block.data[i]);
178                 for (k = 0; k < U->n_col; ++k)
179                         isl_int_mul(A->row[k][i*size+0], D, U->row[0][k]);
180                 isl_int_mul(D, D, d->block.data[i]);
181                 for (j = 1; j < U->n_row; ++j)
182                         for (k = 0; k < U->n_col; ++k)
183                                 isl_int_mul(A->row[k][i*size+j],
184                                                 D, U->row[j][k]);
185         }
186         A = isl_mat_left_hermite(A, 0, NULL, NULL);
187         T = isl_mat_sub_alloc(A, 0, A->n_row, 0, A->n_row);
188         T = isl_mat_lin_to_aff(T);
189         if (!T)
190                 goto error;
191         isl_int_set(T->row[0][0], D);
192         T = isl_mat_right_inverse(T);
193         if (!T)
194                 goto error;
195         isl_assert(T->ctx, isl_int_is_one(T->row[0][0]), goto error);
196         T = isl_mat_transpose(T);
197         isl_mat_free(A);
198         isl_mat_free(U);
199
200         isl_int_clear(D);
201         return T;
202 error:
203         isl_mat_free(A);
204         isl_mat_free(U);
205         isl_int_clear(D);
206         return NULL;
207 }
208
209 /* Given a set of modulo constraints
210  *
211  *              c + A y = 0 mod d
212  *
213  * this function returns an affine transformation T,
214  *
215  *              y = T y'
216  *
217  * that bijectively maps the integer vectors y' to integer
218  * vectors y that satisfy the modulo constraints.
219  *
220  * This function is inspired by Section 2.5.3
221  * of B. Meister, "Stating and Manipulating Periodicity in the Polytope
222  * Model.  Applications to Program Analysis and Optimization".
223  * However, the implementation only follows the algorithm of that
224  * section for computing a particular solution and not for computing
225  * a general homogeneous solution.  The latter is incomplete and
226  * may remove some valid solutions.
227  * Instead, we use an adaptation of the algorithm in Section 7 of
228  * B. Meister, S. Verdoolaege, "Polynomial Approximations in the Polytope
229  * Model: Bringing the Power of Quasi-Polynomials to the Masses".
230  *
231  * The input is given as a matrix B = [ c A ] and a vector d.
232  * Each element of the vector d corresponds to a row in B.
233  * The output is a lower triangular matrix.
234  * If no integer vector y satisfies the given constraints then
235  * a matrix with zero columns is returned.
236  *
237  * We first compute a particular solution y_0 to the given set of
238  * modulo constraints in particular_solution.  If no such solution
239  * exists, then we return a zero-columned transformation matrix.
240  * Otherwise, we compute the generic solution to
241  *
242  *              A y = 0 mod d
243  *
244  * That is we want to compute G such that
245  *
246  *              y = G y''
247  *
248  * with y'' integer, describes the set of solutions.
249  *
250  * We first remove the common factors of each row.
251  * In particular if gcd(A_i,d_i) != 1, then we divide the whole
252  * row i (including d_i) by this common factor.  If afterwards gcd(A_i) != 1,
253  * then we divide this row of A by the common factor, unless gcd(A_i) = 0.
254  * In the later case, we simply drop the row (in both A and d).
255  *
256  * If there are no rows left in A, then G is the identity matrix. Otherwise,
257  * for each row i, we now determine the lattice of integer vectors
258  * that satisfies this row.  Let U_i be the unimodular extension of the
259  * row A_i.  This unimodular extension exists because gcd(A_i) = 1.
260  * The first component of
261  *
262  *              y' = U_i y
263  *
264  * needs to be a multiple of d_i.  Let y' = diag(d_i, 1, ..., 1) y''.
265  * Then,
266  *
267  *              y = U_i^{-1} diag(d_i, 1, ..., 1) y''
268  *
269  * for arbitrary integer vectors y''.  That is, y belongs to the lattice
270  * generated by the columns of L_i = U_i^{-1} diag(d_i, 1, ..., 1).
271  * If there is only one row, then G = L_1.
272  *
273  * If there is more than one row left, we need to compute the intersection
274  * of the lattices.  That is, we need to compute an L such that
275  *
276  *              L = L_i L_i'    for all i
277  *
278  * with L_i' some integer matrices.  Let A be constructed as follows
279  *
280  *              A = [ L_1^{-T} L_2^{-T} ... L_k^{-T} ]
281  *
282  * and computed the Hermite Normal Form of A = [ H 0 ] U
283  * Then,
284  *
285  *              L_i^{-T} = H U_{1,i}
286  *
287  * or
288  *
289  *              H^{-T} = L_i U_{1,i}^T
290  *
291  * In other words G = L = H^{-T}.
292  * To ensure that G is lower triangular, we compute and use its Hermite
293  * normal form.
294  *
295  * The affine transformation matrix returned is then
296  *
297  *              [  1   0  ]
298  *              [ y_0  G  ]
299  *
300  * as any y = y_0 + G y' with y' integer is a solution to the original
301  * modulo constraints.
302  */
303 struct isl_mat *isl_mat_parameter_compression(
304                         struct isl_mat *B, struct isl_vec *d)
305 {
306         int i;
307         struct isl_mat *cst = NULL;
308         struct isl_mat *T = NULL;
309         isl_int D;
310
311         if (!B || !d)
312                 goto error;
313         isl_assert(B->ctx, B->n_row == d->size, goto error);
314         cst = particular_solution(B, d);
315         if (!cst)
316                 goto error;
317         if (cst->n_col == 0) {
318                 T = isl_mat_alloc(B->ctx, B->n_col, 0);
319                 isl_mat_free(cst);
320                 isl_mat_free(B);
321                 isl_vec_free(d);
322                 return T;
323         }
324         isl_int_init(D);
325         /* Replace a*g*row = 0 mod g*m by row = 0 mod m */
326         for (i = 0; i < B->n_row; ++i) {
327                 isl_seq_gcd(B->row[i] + 1, B->n_col - 1, &D);
328                 if (isl_int_is_one(D))
329                         continue;
330                 if (isl_int_is_zero(D)) {
331                         B = isl_mat_drop_rows(B, i, 1);
332                         d = isl_vec_cow(d);
333                         if (!B || !d)
334                                 goto error2;
335                         isl_seq_cpy(d->block.data+i, d->block.data+i+1,
336                                                         d->size - (i+1));
337                         d->size--;
338                         i--;
339                         continue;
340                 }
341                 B = isl_mat_cow(B);
342                 if (!B)
343                         goto error2;
344                 isl_seq_scale_down(B->row[i] + 1, B->row[i] + 1, D, B->n_col-1);
345                 isl_int_gcd(D, D, d->block.data[i]);
346                 d = isl_vec_cow(d);
347                 if (!d)
348                         goto error2;
349                 isl_int_divexact(d->block.data[i], d->block.data[i], D);
350         }
351         isl_int_clear(D);
352         if (B->n_row == 0)
353                 T = isl_mat_identity(B->ctx, B->n_col);
354         else if (B->n_row == 1)
355                 T = parameter_compression_1(B, d);
356         else
357                 T = parameter_compression_multi(B, d);
358         T = isl_mat_left_hermite(T, 0, NULL, NULL);
359         if (!T)
360                 goto error;
361         isl_mat_sub_copy(T->ctx, T->row + 1, cst->row, cst->n_row, 0, 0, 1);
362         isl_mat_free(cst);
363         isl_mat_free(B);
364         isl_vec_free(d);
365         return T;
366 error2:
367         isl_int_clear(D);
368 error:
369         isl_mat_free(cst);
370         isl_mat_free(B);
371         isl_vec_free(d);
372         return NULL;
373 }
374
375 /* Given a set of equalities
376  *
377  *              M x - c = 0
378  *
379  * this function computes a unimodular transformation from a lower-dimensional
380  * space to the original space that bijectively maps the integer points x'
381  * in the lower-dimensional space to the integer points x in the original
382  * space that satisfy the equalities.
383  *
384  * The input is given as a matrix B = [ -c M ] and the output is a
385  * matrix that maps [1 x'] to [1 x].
386  * If T2 is not NULL, then *T2 is set to a matrix mapping [1 x] to [1 x'].
387  *
388  * First compute the (left) Hermite normal form of M,
389  *
390  *              M [U1 U2] = M U = H = [H1 0]
391  * or
392  *                            M = H Q = [H1 0] [Q1]
393  *                                             [Q2]
394  *
395  * with U, Q unimodular, Q = U^{-1} (and H lower triangular).
396  * Define the transformed variables as
397  *
398  *              x = [U1 U2] [ x1' ] = [U1 U2] [Q1] x
399  *                          [ x2' ]           [Q2]
400  *
401  * The equalities then become
402  *
403  *              H1 x1' - c = 0   or   x1' = H1^{-1} c = c'
404  *
405  * If any of the c' is non-integer, then the original set has no
406  * integer solutions (since the x' are a unimodular transformation
407  * of the x) and a zero-column matrix is returned.
408  * Otherwise, the transformation is given by
409  *
410  *              x = U1 H1^{-1} c + U2 x2'
411  *
412  * The inverse transformation is simply
413  *
414  *              x2' = Q2 x
415  */
416 struct isl_mat *isl_mat_variable_compression(struct isl_mat *B,
417         struct isl_mat **T2)
418 {
419         int i;
420         struct isl_mat *H = NULL, *C = NULL, *H1, *U = NULL, *U1, *U2, *TC;
421         unsigned dim;
422
423         if (T2)
424                 *T2 = NULL;
425         if (!B)
426                 goto error;
427
428         dim = B->n_col - 1;
429         H = isl_mat_sub_alloc(B, 0, B->n_row, 1, dim);
430         H = isl_mat_left_hermite(H, 0, &U, T2);
431         if (!H || !U || (T2 && !*T2))
432                 goto error;
433         if (T2) {
434                 *T2 = isl_mat_drop_rows(*T2, 0, B->n_row);
435                 *T2 = isl_mat_lin_to_aff(*T2);
436                 if (!*T2)
437                         goto error;
438         }
439         C = isl_mat_alloc(B->ctx, 1+B->n_row, 1);
440         if (!C)
441                 goto error;
442         isl_int_set_si(C->row[0][0], 1);
443         isl_mat_sub_neg(C->ctx, C->row+1, B->row, B->n_row, 0, 0, 1);
444         H1 = isl_mat_sub_alloc(H, 0, H->n_row, 0, H->n_row);
445         H1 = isl_mat_lin_to_aff(H1);
446         TC = isl_mat_inverse_product(H1, C);
447         if (!TC)
448                 goto error;
449         isl_mat_free(H);
450         if (!isl_int_is_one(TC->row[0][0])) {
451                 for (i = 0; i < B->n_row; ++i) {
452                         if (!isl_int_is_divisible_by(TC->row[1+i][0], TC->row[0][0])) {
453                                 struct isl_ctx *ctx = B->ctx;
454                                 isl_mat_free(B);
455                                 isl_mat_free(TC);
456                                 isl_mat_free(U);
457                                 if (T2) {
458                                         isl_mat_free(*T2);
459                                         *T2 = NULL;
460                                 }
461                                 return isl_mat_alloc(ctx, 1 + dim, 0);
462                         }
463                         isl_seq_scale_down(TC->row[1+i], TC->row[1+i], TC->row[0][0], 1);
464                 }
465                 isl_int_set_si(TC->row[0][0], 1);
466         }
467         U1 = isl_mat_sub_alloc(U, 0, U->n_row, 0, B->n_row);
468         U1 = isl_mat_lin_to_aff(U1);
469         U2 = isl_mat_sub_alloc(U, 0, U->n_row, B->n_row, U->n_row - B->n_row);
470         U2 = isl_mat_lin_to_aff(U2);
471         isl_mat_free(U);
472         TC = isl_mat_product(U1, TC);
473         TC = isl_mat_aff_direct_sum(TC, U2);
474
475         isl_mat_free(B);
476
477         return TC;
478 error:
479         isl_mat_free(B);
480         isl_mat_free(H);
481         isl_mat_free(U);
482         if (T2) {
483                 isl_mat_free(*T2);
484                 *T2 = NULL;
485         }
486         return NULL;
487 }
488
489 /* Use the n equalities of bset to unimodularly transform the
490  * variables x such that n transformed variables x1' have a constant value
491  * and rewrite the constraints of bset in terms of the remaining
492  * transformed variables x2'.  The matrix pointed to by T maps
493  * the new variables x2' back to the original variables x, while T2
494  * maps the original variables to the new variables.
495  */
496 static struct isl_basic_set *compress_variables(
497         struct isl_basic_set *bset, struct isl_mat **T, struct isl_mat **T2)
498 {
499         struct isl_mat *B, *TC;
500         unsigned dim;
501
502         if (T)
503                 *T = NULL;
504         if (T2)
505                 *T2 = NULL;
506         if (!bset)
507                 goto error;
508         isl_assert(bset->ctx, isl_basic_set_n_param(bset) == 0, goto error);
509         isl_assert(bset->ctx, bset->n_div == 0, goto error);
510         dim = isl_basic_set_n_dim(bset);
511         isl_assert(bset->ctx, bset->n_eq <= dim, goto error);
512         if (bset->n_eq == 0)
513                 return bset;
514
515         B = isl_mat_sub_alloc6(bset->ctx, bset->eq, 0, bset->n_eq, 0, 1 + dim);
516         TC = isl_mat_variable_compression(B, T2);
517         if (!TC)
518                 goto error;
519         if (TC->n_col == 0) {
520                 isl_mat_free(TC);
521                 if (T2) {
522                         isl_mat_free(*T2);
523                         *T2 = NULL;
524                 }
525                 return isl_basic_set_set_to_empty(bset);
526         }
527
528         bset = isl_basic_set_preimage(bset, T ? isl_mat_copy(TC) : TC);
529         if (T)
530                 *T = TC;
531         return bset;
532 error:
533         isl_basic_set_free(bset);
534         return NULL;
535 }
536
537 struct isl_basic_set *isl_basic_set_remove_equalities(
538         struct isl_basic_set *bset, struct isl_mat **T, struct isl_mat **T2)
539 {
540         if (T)
541                 *T = NULL;
542         if (T2)
543                 *T2 = NULL;
544         if (!bset)
545                 return NULL;
546         isl_assert(bset->ctx, isl_basic_set_n_param(bset) == 0, goto error);
547         bset = isl_basic_set_gauss(bset, NULL);
548         if (ISL_F_ISSET(bset, ISL_BASIC_SET_EMPTY))
549                 return bset;
550         bset = compress_variables(bset, T, T2);
551         return bset;
552 error:
553         isl_basic_set_free(bset);
554         *T = NULL;
555         return NULL;
556 }
557
558 /* Check if dimension dim belongs to a residue class
559  *              i_dim \equiv r mod m
560  * with m != 1 and if so return m in *modulo and r in *residue.
561  * As a special case, when i_dim has a fixed value v, then
562  * *modulo is set to 0 and *residue to v.
563  *
564  * If i_dim does not belong to such a residue class, then *modulo
565  * is set to 1 and *residue is set to 0.
566  */
567 int isl_basic_set_dim_residue_class(struct isl_basic_set *bset,
568         int pos, isl_int *modulo, isl_int *residue)
569 {
570         struct isl_ctx *ctx;
571         struct isl_mat *H = NULL, *U = NULL, *C, *H1, *U1;
572         unsigned total;
573         unsigned nparam;
574
575         if (!bset || !modulo || !residue)
576                 return -1;
577
578         if (isl_basic_set_plain_dim_is_fixed(bset, pos, residue)) {
579                 isl_int_set_si(*modulo, 0);
580                 return 0;
581         }
582
583         ctx = bset->ctx;
584         total = isl_basic_set_total_dim(bset);
585         nparam = isl_basic_set_n_param(bset);
586         H = isl_mat_sub_alloc6(bset->ctx, bset->eq, 0, bset->n_eq, 1, total);
587         H = isl_mat_left_hermite(H, 0, &U, NULL);
588         if (!H)
589                 return -1;
590
591         isl_seq_gcd(U->row[nparam + pos]+bset->n_eq,
592                         total-bset->n_eq, modulo);
593         if (isl_int_is_zero(*modulo))
594                 isl_int_set_si(*modulo, 1);
595         if (isl_int_is_one(*modulo)) {
596                 isl_int_set_si(*residue, 0);
597                 isl_mat_free(H);
598                 isl_mat_free(U);
599                 return 0;
600         }
601
602         C = isl_mat_alloc(bset->ctx, 1+bset->n_eq, 1);
603         if (!C)
604                 goto error;
605         isl_int_set_si(C->row[0][0], 1);
606         isl_mat_sub_neg(C->ctx, C->row+1, bset->eq, bset->n_eq, 0, 0, 1);
607         H1 = isl_mat_sub_alloc(H, 0, H->n_row, 0, H->n_row);
608         H1 = isl_mat_lin_to_aff(H1);
609         C = isl_mat_inverse_product(H1, C);
610         isl_mat_free(H);
611         U1 = isl_mat_sub_alloc(U, nparam+pos, 1, 0, bset->n_eq);
612         U1 = isl_mat_lin_to_aff(U1);
613         isl_mat_free(U);
614         C = isl_mat_product(U1, C);
615         if (!C)
616                 goto error;
617         if (!isl_int_is_divisible_by(C->row[1][0], C->row[0][0])) {
618                 bset = isl_basic_set_copy(bset);
619                 bset = isl_basic_set_set_to_empty(bset);
620                 isl_basic_set_free(bset);
621                 isl_int_set_si(*modulo, 1);
622                 isl_int_set_si(*residue, 0);
623                 return 0;
624         }
625         isl_int_divexact(*residue, C->row[1][0], C->row[0][0]);
626         isl_int_fdiv_r(*residue, *residue, *modulo);
627         isl_mat_free(C);
628         return 0;
629 error:
630         isl_mat_free(H);
631         isl_mat_free(U);
632         return -1;
633 }
634
635 /* Check if dimension dim belongs to a residue class
636  *              i_dim \equiv r mod m
637  * with m != 1 and if so return m in *modulo and r in *residue.
638  * As a special case, when i_dim has a fixed value v, then
639  * *modulo is set to 0 and *residue to v.
640  *
641  * If i_dim does not belong to such a residue class, then *modulo
642  * is set to 1 and *residue is set to 0.
643  */
644 int isl_set_dim_residue_class(struct isl_set *set,
645         int pos, isl_int *modulo, isl_int *residue)
646 {
647         isl_int m;
648         isl_int r;
649         int i;
650
651         if (!set || !modulo || !residue)
652                 return -1;
653
654         if (set->n == 0) {
655                 isl_int_set_si(*modulo, 0);
656                 isl_int_set_si(*residue, 0);
657                 return 0;
658         }
659
660         if (isl_basic_set_dim_residue_class(set->p[0], pos, modulo, residue)<0)
661                 return -1;
662
663         if (set->n == 1)
664                 return 0;
665
666         if (isl_int_is_one(*modulo))
667                 return 0;
668
669         isl_int_init(m);
670         isl_int_init(r);
671
672         for (i = 1; i < set->n; ++i) {
673                 if (isl_basic_set_dim_residue_class(set->p[i], pos, &m, &r) < 0)
674                         goto error;
675                 isl_int_gcd(*modulo, *modulo, m);
676                 isl_int_sub(m, *residue, r);
677                 isl_int_gcd(*modulo, *modulo, m);
678                 if (!isl_int_is_zero(*modulo))
679                         isl_int_fdiv_r(*residue, *residue, *modulo);
680                 if (isl_int_is_one(*modulo))
681                         break;
682         }
683
684         isl_int_clear(m);
685         isl_int_clear(r);
686
687         return 0;
688 error:
689         isl_int_clear(m);
690         isl_int_clear(r);
691         return -1;
692 }