4 #include "isl_map_private.h"
6 struct isl_mat *isl_mat_alloc(struct isl_ctx *ctx,
7 unsigned n_row, unsigned n_col)
12 mat = isl_alloc_type(ctx, struct isl_mat);
17 mat->block = isl_blk_alloc(ctx, n_row * n_col);
18 if (isl_blk_is_error(mat->block))
20 mat->row = isl_alloc_array(ctx, isl_int *, n_row);
24 for (i = 0; i < n_row; ++i)
25 mat->row[i] = mat->block.data + i * n_col;
37 isl_blk_free(ctx, mat->block);
42 struct isl_mat *isl_mat_extend(struct isl_mat *mat,
43 unsigned n_row, unsigned n_col)
51 if (mat->max_col >= n_col && mat->n_row >= n_row) {
52 if (mat->n_col < n_col)
57 if (mat->max_col < n_col) {
58 struct isl_mat *new_mat;
60 new_mat = isl_mat_alloc(mat->ctx, n_row, n_col);
63 for (i = 0; i < mat->n_row; ++i)
64 isl_seq_cpy(new_mat->row[i], mat->row[i], mat->n_col);
69 mat = isl_mat_cow(mat);
73 assert(mat->ref == 1);
74 old = mat->block.data;
75 mat->block = isl_blk_extend(mat->ctx, mat->block, n_row * mat->max_col);
76 if (isl_blk_is_error(mat->block))
78 mat->row = isl_realloc_array(mat->ctx, mat->row, isl_int *, n_row);
82 for (i = 0; i < mat->n_row; ++i)
83 mat->row[i] = mat->block.data + (mat->row[i] - old);
84 for (i = mat->n_row; i < n_row; ++i)
85 mat->row[i] = mat->block.data + i * mat->max_col;
87 if (mat->n_col < n_col)
96 struct isl_mat *isl_mat_sub_alloc(struct isl_ctx *ctx, isl_int **row,
97 unsigned first_row, unsigned n_row, unsigned first_col, unsigned n_col)
102 mat = isl_alloc_type(ctx, struct isl_mat);
105 mat->row = isl_alloc_array(ctx, isl_int *, n_row);
108 for (i = 0; i < n_row; ++i)
109 mat->row[i] = row[first_row+i] + first_col;
115 mat->block = isl_blk_empty();
116 mat->flags = ISL_MAT_BORROWED;
123 void isl_mat_sub_copy(struct isl_ctx *ctx, isl_int **dst, isl_int **src,
124 unsigned n_row, unsigned dst_col, unsigned src_col, unsigned n_col)
128 for (i = 0; i < n_row; ++i)
129 isl_seq_cpy(dst[i]+dst_col, src[i]+src_col, n_col);
132 void isl_mat_sub_neg(struct isl_ctx *ctx, isl_int **dst, isl_int **src,
133 unsigned n_row, unsigned dst_col, unsigned src_col, unsigned n_col)
137 for (i = 0; i < n_row; ++i)
138 isl_seq_neg(dst[i]+dst_col, src[i]+src_col, n_col);
141 struct isl_mat *isl_mat_copy(struct isl_mat *mat)
150 struct isl_mat *isl_mat_dup(struct isl_mat *mat)
153 struct isl_mat *mat2;
157 mat2 = isl_mat_alloc(mat->ctx, mat->n_row, mat->n_col);
160 for (i = 0; i < mat->n_row; ++i)
161 isl_seq_cpy(mat2->row[i], mat->row[i], mat->n_col);
165 struct isl_mat *isl_mat_cow(struct isl_mat *mat)
167 struct isl_mat *mat2;
171 if (mat->ref == 1 && !ISL_F_ISSET(mat, ISL_MAT_BORROWED))
174 mat2 = isl_mat_dup(mat);
179 void isl_mat_free(struct isl_mat *mat)
187 if (!ISL_F_ISSET(mat, ISL_MAT_BORROWED))
188 isl_blk_free(mat->ctx, mat->block);
189 isl_ctx_deref(mat->ctx);
194 struct isl_mat *isl_mat_identity(struct isl_ctx *ctx, unsigned n_row)
199 mat = isl_mat_alloc(ctx, n_row, n_row);
202 for (i = 0; i < n_row; ++i) {
203 isl_seq_clr(mat->row[i], i);
204 isl_int_set_si(mat->row[i][i], 1);
205 isl_seq_clr(mat->row[i]+i+1, n_row-(i+1));
211 struct isl_vec *isl_mat_vec_product(struct isl_mat *mat, struct isl_vec *vec)
214 struct isl_vec *prod;
219 isl_assert(ctx, mat->n_col == vec->size, goto error);
221 prod = isl_vec_alloc(mat->ctx, mat->n_row);
225 for (i = 0; i < prod->size; ++i)
226 isl_seq_inner_product(mat->row[i], vec->el, vec->size,
227 &prod->block.data[i]);
237 struct isl_mat *isl_mat_aff_direct_sum(struct isl_mat *left,
238 struct isl_mat *right)
246 isl_assert(ctx, left->n_row == right->n_row, goto error);
247 isl_assert(ctx, left->n_row >= 1, goto error);
248 isl_assert(ctx, left->n_col >= 1, goto error);
249 isl_assert(ctx, right->n_col >= 1, goto error);
251 isl_seq_first_non_zero(left->row[0]+1, left->n_col-1) == -1,
254 isl_seq_first_non_zero(right->row[0]+1, right->n_col-1) == -1,
257 sum = isl_mat_alloc(left->ctx, left->n_row, left->n_col + right->n_col - 1);
260 isl_int_lcm(sum->row[0][0], left->row[0][0], right->row[0][0]);
261 isl_int_divexact(left->row[0][0], sum->row[0][0], left->row[0][0]);
262 isl_int_divexact(right->row[0][0], sum->row[0][0], right->row[0][0]);
264 isl_seq_clr(sum->row[0]+1, sum->n_col-1);
265 for (i = 1; i < sum->n_row; ++i) {
266 isl_int_mul(sum->row[i][0], left->row[0][0], left->row[i][0]);
267 isl_int_addmul(sum->row[i][0],
268 right->row[0][0], right->row[i][0]);
269 isl_seq_scale(sum->row[i]+1, left->row[i]+1, left->row[0][0],
271 isl_seq_scale(sum->row[i]+left->n_col,
272 right->row[i]+1, right->row[0][0],
276 isl_int_divexact(left->row[0][0], sum->row[0][0], left->row[0][0]);
277 isl_int_divexact(right->row[0][0], sum->row[0][0], right->row[0][0]);
287 static void exchange(struct isl_mat *M, struct isl_mat **U,
288 struct isl_mat **Q, unsigned row, unsigned i, unsigned j)
291 for (r = row; r < M->n_row; ++r)
292 isl_int_swap(M->row[r][i], M->row[r][j]);
294 for (r = 0; r < (*U)->n_row; ++r)
295 isl_int_swap((*U)->row[r][i], (*U)->row[r][j]);
298 isl_mat_swap_rows(*Q, i, j);
301 static void subtract(struct isl_mat *M, struct isl_mat **U,
302 struct isl_mat **Q, unsigned row, unsigned i, unsigned j, isl_int m)
305 for (r = row; r < M->n_row; ++r)
306 isl_int_submul(M->row[r][j], m, M->row[r][i]);
308 for (r = 0; r < (*U)->n_row; ++r)
309 isl_int_submul((*U)->row[r][j], m, (*U)->row[r][i]);
312 for (r = 0; r < (*Q)->n_col; ++r)
313 isl_int_addmul((*Q)->row[i][r], m, (*Q)->row[j][r]);
317 static void oppose(struct isl_mat *M, struct isl_mat **U,
318 struct isl_mat **Q, unsigned row, unsigned col)
321 for (r = row; r < M->n_row; ++r)
322 isl_int_neg(M->row[r][col], M->row[r][col]);
324 for (r = 0; r < (*U)->n_row; ++r)
325 isl_int_neg((*U)->row[r][col], (*U)->row[r][col]);
328 isl_seq_neg((*Q)->row[col], (*Q)->row[col], (*Q)->n_col);
331 /* Given matrix M, compute
336 * with U and Q unimodular matrices and H a matrix in column echelon form
337 * such that on each echelon row the entries in the non-echelon column
338 * are non-negative (if neg == 0) or non-positive (if neg == 1)
339 * and stricly smaller (in absolute value) than the entries in the echelon
341 * If U or Q are NULL, then these matrices are not computed.
343 struct isl_mat *isl_mat_left_hermite(struct isl_mat *M, int neg,
344 struct isl_mat **U, struct isl_mat **Q)
359 *U = isl_mat_identity(M->ctx, M->n_col);
364 *Q = isl_mat_identity(M->ctx, M->n_col);
371 for (row = 0; row < M->n_row; ++row) {
373 first = isl_seq_abs_min_non_zero(M->row[row]+col, M->n_col-col);
378 exchange(M, U, Q, row, first, col);
379 if (isl_int_is_neg(M->row[row][col]))
380 oppose(M, U, Q, row, col);
382 while ((off = isl_seq_first_non_zero(M->row[row]+first,
383 M->n_col-first)) != -1) {
385 isl_int_fdiv_q(c, M->row[row][first], M->row[row][col]);
386 subtract(M, U, Q, row, col, first, c);
387 if (!isl_int_is_zero(M->row[row][first]))
388 exchange(M, U, Q, row, first, col);
392 for (i = 0; i < col; ++i) {
393 if (isl_int_is_zero(M->row[row][i]))
396 isl_int_cdiv_q(c, M->row[row][i], M->row[row][col]);
398 isl_int_fdiv_q(c, M->row[row][i], M->row[row][col]);
399 if (isl_int_is_zero(c))
401 subtract(M, U, Q, row, col, i, c);
420 struct isl_mat *isl_mat_right_kernel(struct isl_mat *mat)
423 struct isl_mat *U = NULL;
426 mat = isl_mat_left_hermite(mat, 0, &U, NULL);
430 for (i = 0, rank = 0; rank < mat->n_col; ++rank) {
431 while (i < mat->n_row && isl_int_is_zero(mat->row[i][rank]))
436 K = isl_mat_alloc(U->ctx, U->n_row, U->n_col - rank);
439 isl_mat_sub_copy(K->ctx, K->row, U->row, U->n_row, 0, rank, U->n_col-rank);
449 struct isl_mat *isl_mat_lin_to_aff(struct isl_mat *mat)
452 struct isl_mat *mat2;
456 mat2 = isl_mat_alloc(mat->ctx, 1+mat->n_row, 1+mat->n_col);
459 isl_int_set_si(mat2->row[0][0], 1);
460 isl_seq_clr(mat2->row[0]+1, mat->n_col);
461 for (i = 0; i < mat->n_row; ++i) {
462 isl_int_set_si(mat2->row[1+i][0], 0);
463 isl_seq_cpy(mat2->row[1+i]+1, mat->row[i], mat->n_col);
469 static int row_first_non_zero(isl_int **row, unsigned n_row, unsigned col)
473 for (i = 0; i < n_row; ++i)
474 if (!isl_int_is_zero(row[i][col]))
479 static int row_abs_min_non_zero(isl_int **row, unsigned n_row, unsigned col)
481 int i, min = row_first_non_zero(row, n_row, col);
484 for (i = min + 1; i < n_row; ++i) {
485 if (isl_int_is_zero(row[i][col]))
487 if (isl_int_abs_lt(row[i][col], row[min][col]))
493 static void inv_exchange(struct isl_mat *left, struct isl_mat *right,
494 unsigned i, unsigned j)
496 left = isl_mat_swap_rows(left, i, j);
497 right = isl_mat_swap_rows(right, i, j);
500 static void inv_oppose(
501 struct isl_mat *left, struct isl_mat *right, unsigned row)
503 isl_seq_neg(left->row[row]+row, left->row[row]+row, left->n_col-row);
504 isl_seq_neg(right->row[row], right->row[row], right->n_col);
507 static void inv_subtract(struct isl_mat *left, struct isl_mat *right,
508 unsigned row, unsigned i, isl_int m)
511 isl_seq_combine(left->row[i]+row,
512 left->ctx->one, left->row[i]+row,
513 m, left->row[row]+row,
515 isl_seq_combine(right->row[i], right->ctx->one, right->row[i],
516 m, right->row[row], right->n_col);
519 /* Compute inv(left)*right
521 struct isl_mat *isl_mat_inverse_product(struct isl_mat *left,
522 struct isl_mat *right)
530 isl_assert(left->ctx, left->n_row == left->n_col, goto error);
531 isl_assert(left->ctx, left->n_row == right->n_row, goto error);
533 if (left->n_row == 0) {
538 left = isl_mat_cow(left);
539 right = isl_mat_cow(right);
545 for (row = 0; row < left->n_row; ++row) {
546 int pivot, first, i, off;
547 pivot = row_abs_min_non_zero(left->row+row, left->n_row-row, row);
551 isl_assert(ctx, pivot >= 0, goto error);
555 inv_exchange(left, right, pivot, row);
556 if (isl_int_is_neg(left->row[row][row]))
557 inv_oppose(left, right, row);
559 while ((off = row_first_non_zero(left->row+first,
560 left->n_row-first, row)) != -1) {
562 isl_int_fdiv_q(a, left->row[first][row],
563 left->row[row][row]);
564 inv_subtract(left, right, row, first, a);
565 if (!isl_int_is_zero(left->row[first][row]))
566 inv_exchange(left, right, row, first);
570 for (i = 0; i < row; ++i) {
571 if (isl_int_is_zero(left->row[i][row]))
573 isl_int_gcd(a, left->row[row][row], left->row[i][row]);
574 isl_int_divexact(b, left->row[i][row], a);
575 isl_int_divexact(a, left->row[row][row], a);
577 isl_seq_combine(left->row[i]+row,
579 b, left->row[row]+row,
581 isl_seq_combine(right->row[i], a, right->row[i],
582 b, right->row[row], right->n_col);
587 isl_int_set(a, left->row[0][0]);
588 for (row = 1; row < left->n_row; ++row)
589 isl_int_lcm(a, a, left->row[row][row]);
590 if (isl_int_is_zero(a)){
592 isl_assert(ctx, 0, goto error);
594 for (row = 0; row < left->n_row; ++row) {
595 isl_int_divexact(left->row[row][row], a, left->row[row][row]);
596 if (isl_int_is_one(left->row[row][row]))
598 isl_seq_scale(right->row[row], right->row[row],
599 left->row[row][row], right->n_col);
611 void isl_mat_col_scale(struct isl_mat *mat, unsigned col, isl_int m)
615 for (i = 0; i < mat->n_row; ++i)
616 isl_int_mul(mat->row[i][col], mat->row[i][col], m);
619 void isl_mat_col_combine(struct isl_mat *mat, unsigned dst,
620 isl_int m1, unsigned src1, isl_int m2, unsigned src2)
626 for (i = 0; i < mat->n_row; ++i) {
627 isl_int_mul(tmp, m1, mat->row[i][src1]);
628 isl_int_addmul(tmp, m2, mat->row[i][src2]);
629 isl_int_set(mat->row[i][dst], tmp);
634 struct isl_mat *isl_mat_right_inverse(struct isl_mat *mat)
640 mat = isl_mat_cow(mat);
644 inv = isl_mat_identity(mat->ctx, mat->n_col);
645 inv = isl_mat_cow(inv);
651 for (row = 0; row < mat->n_row; ++row) {
652 int pivot, first, i, off;
653 pivot = isl_seq_abs_min_non_zero(mat->row[row]+row, mat->n_col-row);
661 exchange(mat, &inv, NULL, row, pivot, row);
662 if (isl_int_is_neg(mat->row[row][row]))
663 oppose(mat, &inv, NULL, row, row);
665 while ((off = isl_seq_first_non_zero(mat->row[row]+first,
666 mat->n_col-first)) != -1) {
668 isl_int_fdiv_q(a, mat->row[row][first],
670 subtract(mat, &inv, NULL, row, row, first, a);
671 if (!isl_int_is_zero(mat->row[row][first]))
672 exchange(mat, &inv, NULL, row, row, first);
676 for (i = 0; i < row; ++i) {
677 if (isl_int_is_zero(mat->row[row][i]))
679 isl_int_gcd(a, mat->row[row][row], mat->row[row][i]);
680 isl_int_divexact(b, mat->row[row][i], a);
681 isl_int_divexact(a, mat->row[row][row], a);
683 isl_mat_col_combine(mat, i, a, i, b, row);
684 isl_mat_col_combine(inv, i, a, i, b, row);
689 isl_int_set(a, mat->row[0][0]);
690 for (row = 1; row < mat->n_row; ++row)
691 isl_int_lcm(a, a, mat->row[row][row]);
692 if (isl_int_is_zero(a)){
696 for (row = 0; row < mat->n_row; ++row) {
697 isl_int_divexact(mat->row[row][row], a, mat->row[row][row]);
698 if (isl_int_is_one(mat->row[row][row]))
700 isl_mat_col_scale(inv, row, mat->row[row][row]);
712 struct isl_mat *isl_mat_transpose(struct isl_mat *mat)
714 struct isl_mat *transpose = NULL;
717 if (mat->n_col == mat->n_row) {
718 mat = isl_mat_cow(mat);
721 for (i = 0; i < mat->n_row; ++i)
722 for (j = i + 1; j < mat->n_col; ++j)
723 isl_int_swap(mat->row[i][j], mat->row[j][i]);
726 transpose = isl_mat_alloc(mat->ctx, mat->n_col, mat->n_row);
729 for (i = 0; i < mat->n_row; ++i)
730 for (j = 0; j < mat->n_col; ++j)
731 isl_int_set(transpose->row[j][i], mat->row[i][j]);
739 struct isl_mat *isl_mat_swap_cols(struct isl_mat *mat, unsigned i, unsigned j)
743 mat = isl_mat_cow(mat);
746 isl_assert(ctx, i < mat->n_col, goto error);
747 isl_assert(ctx, j < mat->n_col, goto error);
749 for (r = 0; r < mat->n_row; ++r)
750 isl_int_swap(mat->row[r][i], mat->row[r][j]);
757 struct isl_mat *isl_mat_swap_rows(struct isl_mat *mat, unsigned i, unsigned j)
763 mat = isl_mat_cow(mat);
767 mat->row[i] = mat->row[j];
772 struct isl_mat *isl_mat_product(struct isl_mat *left, struct isl_mat *right)
775 struct isl_mat *prod;
779 isl_assert(ctx, left->n_col == right->n_row, goto error);
780 prod = isl_mat_alloc(left->ctx, left->n_row, right->n_col);
783 if (left->n_col == 0) {
784 for (i = 0; i < prod->n_row; ++i)
785 isl_seq_clr(prod->row[i], prod->n_col);
788 for (i = 0; i < prod->n_row; ++i) {
789 for (j = 0; j < prod->n_col; ++j) {
790 isl_int_mul(prod->row[i][j],
791 left->row[i][0], right->row[0][j]);
792 for (k = 1; k < left->n_col; ++k)
793 isl_int_addmul(prod->row[i][j],
794 left->row[i][k], right->row[k][j]);
806 /* Replace the variables x in the rows q by x' given by x = M x',
807 * with M the matrix mat.
809 * If the number of new variables is greater than the original
810 * number of variables, then the rows q have already been
811 * preextended. If the new number is smaller, then the coefficients
812 * of the divs, which are not changed, need to be shifted down.
813 * The row q may be the equalities, the inequalities or the
814 * div expressions. In the latter case, has_div is true and
815 * we need to take into account the extra denominator column.
817 static int preimage(struct isl_ctx *ctx, isl_int **q, unsigned n,
818 unsigned n_div, int has_div, struct isl_mat *mat)
824 if (mat->n_col >= mat->n_row)
827 e = mat->n_row - mat->n_col;
829 for (i = 0; i < n; ++i)
830 isl_int_mul(q[i][0], q[i][0], mat->row[0][0]);
831 t = isl_mat_sub_alloc(mat->ctx, q, 0, n, has_div, mat->n_row);
832 t = isl_mat_product(t, mat);
835 for (i = 0; i < n; ++i) {
836 isl_seq_swp_or_cpy(q[i] + has_div, t->row[i], t->n_col);
837 isl_seq_cpy(q[i] + has_div + t->n_col,
838 q[i] + has_div + t->n_col + e, n_div);
839 isl_seq_clr(q[i] + has_div + t->n_col + n_div, e);
845 /* Replace the variables x in bset by x' given by x = M x', with
848 * If there are fewer variables x' then there are x, then we perform
849 * the transformation in place, which that, in principle,
850 * this frees up some extra variables as the number
851 * of columns remains constant, but we would have to extend
852 * the div array too as the number of rows in this array is assumed
853 * to be equal to extra.
855 struct isl_basic_set *isl_basic_set_preimage(struct isl_basic_set *bset,
864 bset = isl_basic_set_cow(bset);
868 isl_assert(ctx, bset->dim->nparam == 0, goto error);
869 isl_assert(ctx, 1+bset->dim->n_out == mat->n_row, goto error);
871 if (mat->n_col > mat->n_row)
872 bset = isl_basic_set_extend(bset, 0, mat->n_col-1, 0,
874 else if (mat->n_col < mat->n_row) {
875 bset->dim = isl_dim_cow(bset->dim);
878 bset->dim->n_out -= mat->n_row - mat->n_col;
881 if (preimage(ctx, bset->eq, bset->n_eq, bset->n_div, 0,
882 isl_mat_copy(mat)) < 0)
885 if (preimage(ctx, bset->ineq, bset->n_ineq, bset->n_div, 0,
886 isl_mat_copy(mat)) < 0)
889 if (preimage(ctx, bset->div, bset->n_div, bset->n_div, 1, mat) < 0)
892 ISL_F_CLR(bset, ISL_BASIC_SET_NO_IMPLICIT);
893 ISL_F_CLR(bset, ISL_BASIC_SET_NO_REDUNDANT);
894 ISL_F_CLR(bset, ISL_BASIC_SET_NORMALIZED);
895 ISL_F_CLR(bset, ISL_BASIC_SET_NORMALIZED_DIVS);
896 ISL_F_CLR(bset, ISL_BASIC_SET_ALL_EQUALITIES);
898 bset = isl_basic_set_simplify(bset);
899 bset = isl_basic_set_finalize(bset);
905 isl_basic_set_free(bset);
909 struct isl_set *isl_set_preimage(struct isl_set *set, struct isl_mat *mat)
914 set = isl_set_cow(set);
919 for (i = 0; i < set->n; ++i) {
920 set->p[i] = isl_basic_set_preimage(set->p[i],
925 if (mat->n_col != mat->n_row) {
926 set->dim = isl_dim_cow(set->dim);
929 set->dim->n_out += mat->n_col;
930 set->dim->n_out -= mat->n_row;
933 ISL_F_CLR(set, ISL_SET_NORMALIZED);
941 void isl_mat_dump(struct isl_mat *mat, FILE *out, int indent)
946 fprintf(out, "%*snull mat\n", indent, "");
951 fprintf(out, "%*s[]\n", indent, "");
953 for (i = 0; i < mat->n_row; ++i) {
955 fprintf(out, "%*s[[", indent, "");
957 fprintf(out, "%*s[", indent+1, "");
958 for (j = 0; j < mat->n_col; ++j) {
961 isl_int_print(out, mat->row[i][j], 0);
963 if (i == mat->n_row-1)
964 fprintf(out, "]]\n");
970 struct isl_mat *isl_mat_drop_cols(struct isl_mat *mat, unsigned col, unsigned n)
974 mat = isl_mat_cow(mat);
978 if (col != mat->n_col-n) {
979 for (r = 0; r < mat->n_row; ++r)
980 isl_seq_cpy(mat->row[r]+col, mat->row[r]+col+n,
981 mat->n_col - col - n);
987 struct isl_mat *isl_mat_drop_rows(struct isl_mat *mat, unsigned row, unsigned n)
991 mat = isl_mat_cow(mat);
995 for (r = row; r+n < mat->n_row; ++r)
996 mat->row[r] = mat->row[r+n];
1002 void isl_mat_col_submul(struct isl_mat *mat,
1003 int dst_col, isl_int f, int src_col)
1007 for (i = 0; i < mat->n_row; ++i)
1008 isl_int_submul(mat->row[i][dst_col], f, mat->row[i][src_col]);
1011 void isl_mat_col_mul(struct isl_mat *mat, int dst_col, isl_int f, int src_col)
1015 for (i = 0; i < mat->n_row; ++i)
1016 isl_int_mul(mat->row[i][dst_col], f, mat->row[i][src_col]);
1019 struct isl_mat *isl_mat_unimodular_complete(struct isl_mat *M, int row)
1022 struct isl_mat *H = NULL, *Q = NULL;
1024 isl_assert(ctx, M->n_row == M->n_col, goto error);
1026 H = isl_mat_left_hermite(isl_mat_copy(M), 0, NULL, &Q);
1027 M->n_row = M->n_col;
1030 for (r = 0; r < row; ++r)
1031 isl_assert(ctx, isl_int_is_one(H->row[r][r]), goto error);
1032 for (r = row; r < M->n_row; ++r)
1033 isl_seq_cpy(M->row[r], Q->row[r], M->n_col);