2 * Copyright 2008-2009 Katholieke Universiteit Leuven
4 * Use of this software is governed by the GNU LGPLv2.1 license
6 * Written by Sven Verdoolaege, K.U.Leuven, Departement
7 * Computerwetenschappen, Celestijnenlaan 200A, B-3001 Leuven, Belgium
10 #include <isl_ctx_private.h>
13 #include <isl_mat_private.h>
14 #include "isl_map_private.h"
15 #include <isl_dim_private.h>
17 isl_ctx *isl_mat_get_ctx(__isl_keep isl_mat *mat)
19 return mat ? mat->ctx : NULL;
22 struct isl_mat *isl_mat_alloc(struct isl_ctx *ctx,
23 unsigned n_row, unsigned n_col)
28 mat = isl_alloc_type(ctx, struct isl_mat);
33 mat->block = isl_blk_alloc(ctx, n_row * n_col);
34 if (isl_blk_is_error(mat->block))
36 mat->row = isl_alloc_array(ctx, isl_int *, n_row);
40 for (i = 0; i < n_row; ++i)
41 mat->row[i] = mat->block.data + i * n_col;
53 isl_blk_free(ctx, mat->block);
58 struct isl_mat *isl_mat_extend(struct isl_mat *mat,
59 unsigned n_row, unsigned n_col)
68 if (mat->max_col >= n_col && mat->n_row >= n_row) {
69 if (mat->n_col < n_col)
74 if (mat->max_col < n_col) {
75 struct isl_mat *new_mat;
77 if (n_row < mat->n_row)
79 new_mat = isl_mat_alloc(mat->ctx, n_row, n_col);
82 for (i = 0; i < mat->n_row; ++i)
83 isl_seq_cpy(new_mat->row[i], mat->row[i], mat->n_col);
88 mat = isl_mat_cow(mat);
92 old = mat->block.data;
93 mat->block = isl_blk_extend(mat->ctx, mat->block, n_row * mat->max_col);
94 if (isl_blk_is_error(mat->block))
96 row = isl_realloc_array(mat->ctx, mat->row, isl_int *, n_row);
101 for (i = 0; i < mat->n_row; ++i)
102 mat->row[i] = mat->block.data + (mat->row[i] - old);
103 for (i = mat->n_row; i < n_row; ++i)
104 mat->row[i] = mat->block.data + i * mat->max_col;
106 if (mat->n_col < n_col)
115 __isl_give isl_mat *isl_mat_sub_alloc6(isl_ctx *ctx, isl_int **row,
116 unsigned first_row, unsigned n_row, unsigned first_col, unsigned n_col)
121 mat = isl_alloc_type(ctx, struct isl_mat);
124 mat->row = isl_alloc_array(ctx, isl_int *, n_row);
127 for (i = 0; i < n_row; ++i)
128 mat->row[i] = row[first_row+i] + first_col;
134 mat->block = isl_blk_empty();
135 mat->flags = ISL_MAT_BORROWED;
142 __isl_give isl_mat *isl_mat_sub_alloc(__isl_keep isl_mat *mat,
143 unsigned first_row, unsigned n_row, unsigned first_col, unsigned n_col)
147 return isl_mat_sub_alloc6(mat->ctx, mat->row, first_row, n_row,
151 void isl_mat_sub_copy(struct isl_ctx *ctx, isl_int **dst, isl_int **src,
152 unsigned n_row, unsigned dst_col, unsigned src_col, unsigned n_col)
156 for (i = 0; i < n_row; ++i)
157 isl_seq_cpy(dst[i]+dst_col, src[i]+src_col, n_col);
160 void isl_mat_sub_neg(struct isl_ctx *ctx, isl_int **dst, isl_int **src,
161 unsigned n_row, unsigned dst_col, unsigned src_col, unsigned n_col)
165 for (i = 0; i < n_row; ++i)
166 isl_seq_neg(dst[i]+dst_col, src[i]+src_col, n_col);
169 struct isl_mat *isl_mat_copy(struct isl_mat *mat)
178 struct isl_mat *isl_mat_dup(struct isl_mat *mat)
181 struct isl_mat *mat2;
185 mat2 = isl_mat_alloc(mat->ctx, mat->n_row, mat->n_col);
188 for (i = 0; i < mat->n_row; ++i)
189 isl_seq_cpy(mat2->row[i], mat->row[i], mat->n_col);
193 struct isl_mat *isl_mat_cow(struct isl_mat *mat)
195 struct isl_mat *mat2;
199 if (mat->ref == 1 && !ISL_F_ISSET(mat, ISL_MAT_BORROWED))
202 mat2 = isl_mat_dup(mat);
207 void isl_mat_free(struct isl_mat *mat)
215 if (!ISL_F_ISSET(mat, ISL_MAT_BORROWED))
216 isl_blk_free(mat->ctx, mat->block);
217 isl_ctx_deref(mat->ctx);
222 int isl_mat_rows(__isl_keep isl_mat *mat)
224 return mat ? mat->n_row : -1;
227 int isl_mat_cols(__isl_keep isl_mat *mat)
229 return mat ? mat->n_col : -1;
232 int isl_mat_get_element(__isl_keep isl_mat *mat, int row, int col, isl_int *v)
236 if (row < 0 || row >= mat->n_row)
237 isl_die(mat->ctx, isl_error_invalid, "row out of range",
239 if (col < 0 || col >= mat->n_col)
240 isl_die(mat->ctx, isl_error_invalid, "column out of range",
242 isl_int_set(*v, mat->row[row][col]);
246 __isl_give isl_mat *isl_mat_set_element(__isl_take isl_mat *mat,
247 int row, int col, isl_int v)
249 mat = isl_mat_cow(mat);
252 if (row < 0 || row >= mat->n_row)
253 isl_die(mat->ctx, isl_error_invalid, "row out of range",
255 if (col < 0 || col >= mat->n_col)
256 isl_die(mat->ctx, isl_error_invalid, "column out of range",
258 isl_int_set(mat->row[row][col], v);
265 __isl_give isl_mat *isl_mat_set_element_si(__isl_take isl_mat *mat,
266 int row, int col, int v)
268 mat = isl_mat_cow(mat);
271 if (row < 0 || row >= mat->n_row)
272 isl_die(mat->ctx, isl_error_invalid, "row out of range",
274 if (col < 0 || col >= mat->n_col)
275 isl_die(mat->ctx, isl_error_invalid, "column out of range",
277 isl_int_set_si(mat->row[row][col], v);
284 struct isl_mat *isl_mat_identity(struct isl_ctx *ctx, unsigned n_row)
289 mat = isl_mat_alloc(ctx, n_row, n_row);
292 for (i = 0; i < n_row; ++i) {
293 isl_seq_clr(mat->row[i], i);
294 isl_int_set_si(mat->row[i][i], 1);
295 isl_seq_clr(mat->row[i]+i+1, n_row-(i+1));
301 struct isl_vec *isl_mat_vec_product(struct isl_mat *mat, struct isl_vec *vec)
304 struct isl_vec *prod;
309 isl_assert(mat->ctx, mat->n_col == vec->size, goto error);
311 prod = isl_vec_alloc(mat->ctx, mat->n_row);
315 for (i = 0; i < prod->size; ++i)
316 isl_seq_inner_product(mat->row[i], vec->el, vec->size,
317 &prod->block.data[i]);
327 __isl_give isl_vec *isl_mat_vec_inverse_product(__isl_take isl_mat *mat,
328 __isl_take isl_vec *vec)
330 struct isl_mat *vec_mat;
335 vec_mat = isl_mat_alloc(vec->ctx, vec->size, 1);
338 for (i = 0; i < vec->size; ++i)
339 isl_int_set(vec_mat->row[i][0], vec->el[i]);
340 vec_mat = isl_mat_inverse_product(mat, vec_mat);
344 vec = isl_vec_alloc(vec_mat->ctx, vec_mat->n_row);
346 for (i = 0; i < vec->size; ++i)
347 isl_int_set(vec->el[i], vec_mat->row[i][0]);
348 isl_mat_free(vec_mat);
356 struct isl_vec *isl_vec_mat_product(struct isl_vec *vec, struct isl_mat *mat)
359 struct isl_vec *prod;
364 isl_assert(mat->ctx, mat->n_row == vec->size, goto error);
366 prod = isl_vec_alloc(mat->ctx, mat->n_col);
370 for (i = 0; i < prod->size; ++i) {
371 isl_int_set_si(prod->el[i], 0);
372 for (j = 0; j < vec->size; ++j)
373 isl_int_addmul(prod->el[i], vec->el[j], mat->row[j][i]);
384 struct isl_mat *isl_mat_aff_direct_sum(struct isl_mat *left,
385 struct isl_mat *right)
393 isl_assert(left->ctx, left->n_row == right->n_row, goto error);
394 isl_assert(left->ctx, left->n_row >= 1, goto error);
395 isl_assert(left->ctx, left->n_col >= 1, goto error);
396 isl_assert(left->ctx, right->n_col >= 1, goto error);
397 isl_assert(left->ctx,
398 isl_seq_first_non_zero(left->row[0]+1, left->n_col-1) == -1,
400 isl_assert(left->ctx,
401 isl_seq_first_non_zero(right->row[0]+1, right->n_col-1) == -1,
404 sum = isl_mat_alloc(left->ctx, left->n_row, left->n_col + right->n_col - 1);
407 isl_int_lcm(sum->row[0][0], left->row[0][0], right->row[0][0]);
408 isl_int_divexact(left->row[0][0], sum->row[0][0], left->row[0][0]);
409 isl_int_divexact(right->row[0][0], sum->row[0][0], right->row[0][0]);
411 isl_seq_clr(sum->row[0]+1, sum->n_col-1);
412 for (i = 1; i < sum->n_row; ++i) {
413 isl_int_mul(sum->row[i][0], left->row[0][0], left->row[i][0]);
414 isl_int_addmul(sum->row[i][0],
415 right->row[0][0], right->row[i][0]);
416 isl_seq_scale(sum->row[i]+1, left->row[i]+1, left->row[0][0],
418 isl_seq_scale(sum->row[i]+left->n_col,
419 right->row[i]+1, right->row[0][0],
423 isl_int_divexact(left->row[0][0], sum->row[0][0], left->row[0][0]);
424 isl_int_divexact(right->row[0][0], sum->row[0][0], right->row[0][0]);
434 static void exchange(struct isl_mat *M, struct isl_mat **U,
435 struct isl_mat **Q, unsigned row, unsigned i, unsigned j)
438 for (r = row; r < M->n_row; ++r)
439 isl_int_swap(M->row[r][i], M->row[r][j]);
441 for (r = 0; r < (*U)->n_row; ++r)
442 isl_int_swap((*U)->row[r][i], (*U)->row[r][j]);
445 isl_mat_swap_rows(*Q, i, j);
448 static void subtract(struct isl_mat *M, struct isl_mat **U,
449 struct isl_mat **Q, unsigned row, unsigned i, unsigned j, isl_int m)
452 for (r = row; r < M->n_row; ++r)
453 isl_int_submul(M->row[r][j], m, M->row[r][i]);
455 for (r = 0; r < (*U)->n_row; ++r)
456 isl_int_submul((*U)->row[r][j], m, (*U)->row[r][i]);
459 for (r = 0; r < (*Q)->n_col; ++r)
460 isl_int_addmul((*Q)->row[i][r], m, (*Q)->row[j][r]);
464 static void oppose(struct isl_mat *M, struct isl_mat **U,
465 struct isl_mat **Q, unsigned row, unsigned col)
468 for (r = row; r < M->n_row; ++r)
469 isl_int_neg(M->row[r][col], M->row[r][col]);
471 for (r = 0; r < (*U)->n_row; ++r)
472 isl_int_neg((*U)->row[r][col], (*U)->row[r][col]);
475 isl_seq_neg((*Q)->row[col], (*Q)->row[col], (*Q)->n_col);
478 /* Given matrix M, compute
483 * with U and Q unimodular matrices and H a matrix in column echelon form
484 * such that on each echelon row the entries in the non-echelon column
485 * are non-negative (if neg == 0) or non-positive (if neg == 1)
486 * and stricly smaller (in absolute value) than the entries in the echelon
488 * If U or Q are NULL, then these matrices are not computed.
490 struct isl_mat *isl_mat_left_hermite(struct isl_mat *M, int neg,
491 struct isl_mat **U, struct isl_mat **Q)
506 *U = isl_mat_identity(M->ctx, M->n_col);
511 *Q = isl_mat_identity(M->ctx, M->n_col);
518 for (row = 0; row < M->n_row; ++row) {
520 first = isl_seq_abs_min_non_zero(M->row[row]+col, M->n_col-col);
525 exchange(M, U, Q, row, first, col);
526 if (isl_int_is_neg(M->row[row][col]))
527 oppose(M, U, Q, row, col);
529 while ((off = isl_seq_first_non_zero(M->row[row]+first,
530 M->n_col-first)) != -1) {
532 isl_int_fdiv_q(c, M->row[row][first], M->row[row][col]);
533 subtract(M, U, Q, row, col, first, c);
534 if (!isl_int_is_zero(M->row[row][first]))
535 exchange(M, U, Q, row, first, col);
539 for (i = 0; i < col; ++i) {
540 if (isl_int_is_zero(M->row[row][i]))
543 isl_int_cdiv_q(c, M->row[row][i], M->row[row][col]);
545 isl_int_fdiv_q(c, M->row[row][i], M->row[row][col]);
546 if (isl_int_is_zero(c))
548 subtract(M, U, Q, row, col, i, c);
568 struct isl_mat *isl_mat_right_kernel(struct isl_mat *mat)
571 struct isl_mat *U = NULL;
574 mat = isl_mat_left_hermite(mat, 0, &U, NULL);
578 for (i = 0, rank = 0; rank < mat->n_col; ++rank) {
579 while (i < mat->n_row && isl_int_is_zero(mat->row[i][rank]))
584 K = isl_mat_alloc(U->ctx, U->n_row, U->n_col - rank);
587 isl_mat_sub_copy(K->ctx, K->row, U->row, U->n_row, 0, rank, U->n_col-rank);
597 struct isl_mat *isl_mat_lin_to_aff(struct isl_mat *mat)
600 struct isl_mat *mat2;
604 mat2 = isl_mat_alloc(mat->ctx, 1+mat->n_row, 1+mat->n_col);
607 isl_int_set_si(mat2->row[0][0], 1);
608 isl_seq_clr(mat2->row[0]+1, mat->n_col);
609 for (i = 0; i < mat->n_row; ++i) {
610 isl_int_set_si(mat2->row[1+i][0], 0);
611 isl_seq_cpy(mat2->row[1+i]+1, mat->row[i], mat->n_col);
620 /* Given two matrices M1 and M2, return the block matrix
625 __isl_give isl_mat *isl_mat_diagonal(__isl_take isl_mat *mat1,
626 __isl_take isl_mat *mat2)
634 mat = isl_mat_alloc(mat1->ctx, mat1->n_row + mat2->n_row,
635 mat1->n_col + mat2->n_col);
638 for (i = 0; i < mat1->n_row; ++i) {
639 isl_seq_cpy(mat->row[i], mat1->row[i], mat1->n_col);
640 isl_seq_clr(mat->row[i] + mat1->n_col, mat2->n_col);
642 for (i = 0; i < mat2->n_row; ++i) {
643 isl_seq_clr(mat->row[mat1->n_row + i], mat1->n_col);
644 isl_seq_cpy(mat->row[mat1->n_row + i] + mat1->n_col,
645 mat2->row[i], mat2->n_col);
656 static int row_first_non_zero(isl_int **row, unsigned n_row, unsigned col)
660 for (i = 0; i < n_row; ++i)
661 if (!isl_int_is_zero(row[i][col]))
666 static int row_abs_min_non_zero(isl_int **row, unsigned n_row, unsigned col)
668 int i, min = row_first_non_zero(row, n_row, col);
671 for (i = min + 1; i < n_row; ++i) {
672 if (isl_int_is_zero(row[i][col]))
674 if (isl_int_abs_lt(row[i][col], row[min][col]))
680 static void inv_exchange(struct isl_mat *left, struct isl_mat *right,
681 unsigned i, unsigned j)
683 left = isl_mat_swap_rows(left, i, j);
684 right = isl_mat_swap_rows(right, i, j);
687 static void inv_oppose(
688 struct isl_mat *left, struct isl_mat *right, unsigned row)
690 isl_seq_neg(left->row[row]+row, left->row[row]+row, left->n_col-row);
691 isl_seq_neg(right->row[row], right->row[row], right->n_col);
694 static void inv_subtract(struct isl_mat *left, struct isl_mat *right,
695 unsigned row, unsigned i, isl_int m)
698 isl_seq_combine(left->row[i]+row,
699 left->ctx->one, left->row[i]+row,
700 m, left->row[row]+row,
702 isl_seq_combine(right->row[i], right->ctx->one, right->row[i],
703 m, right->row[row], right->n_col);
706 /* Compute inv(left)*right
708 struct isl_mat *isl_mat_inverse_product(struct isl_mat *left,
709 struct isl_mat *right)
717 isl_assert(left->ctx, left->n_row == left->n_col, goto error);
718 isl_assert(left->ctx, left->n_row == right->n_row, goto error);
720 if (left->n_row == 0) {
725 left = isl_mat_cow(left);
726 right = isl_mat_cow(right);
732 for (row = 0; row < left->n_row; ++row) {
733 int pivot, first, i, off;
734 pivot = row_abs_min_non_zero(left->row+row, left->n_row-row, row);
738 isl_assert(left->ctx, pivot >= 0, goto error);
742 inv_exchange(left, right, pivot, row);
743 if (isl_int_is_neg(left->row[row][row]))
744 inv_oppose(left, right, row);
746 while ((off = row_first_non_zero(left->row+first,
747 left->n_row-first, row)) != -1) {
749 isl_int_fdiv_q(a, left->row[first][row],
750 left->row[row][row]);
751 inv_subtract(left, right, row, first, a);
752 if (!isl_int_is_zero(left->row[first][row]))
753 inv_exchange(left, right, row, first);
757 for (i = 0; i < row; ++i) {
758 if (isl_int_is_zero(left->row[i][row]))
760 isl_int_gcd(a, left->row[row][row], left->row[i][row]);
761 isl_int_divexact(b, left->row[i][row], a);
762 isl_int_divexact(a, left->row[row][row], a);
764 isl_seq_combine(left->row[i] + i,
766 b, left->row[row] + i,
768 isl_seq_combine(right->row[i], a, right->row[i],
769 b, right->row[row], right->n_col);
774 isl_int_set(a, left->row[0][0]);
775 for (row = 1; row < left->n_row; ++row)
776 isl_int_lcm(a, a, left->row[row][row]);
777 if (isl_int_is_zero(a)){
779 isl_assert(left->ctx, 0, goto error);
781 for (row = 0; row < left->n_row; ++row) {
782 isl_int_divexact(left->row[row][row], a, left->row[row][row]);
783 if (isl_int_is_one(left->row[row][row]))
785 isl_seq_scale(right->row[row], right->row[row],
786 left->row[row][row], right->n_col);
798 void isl_mat_col_scale(struct isl_mat *mat, unsigned col, isl_int m)
802 for (i = 0; i < mat->n_row; ++i)
803 isl_int_mul(mat->row[i][col], mat->row[i][col], m);
806 void isl_mat_col_combine(struct isl_mat *mat, unsigned dst,
807 isl_int m1, unsigned src1, isl_int m2, unsigned src2)
813 for (i = 0; i < mat->n_row; ++i) {
814 isl_int_mul(tmp, m1, mat->row[i][src1]);
815 isl_int_addmul(tmp, m2, mat->row[i][src2]);
816 isl_int_set(mat->row[i][dst], tmp);
821 struct isl_mat *isl_mat_right_inverse(struct isl_mat *mat)
827 mat = isl_mat_cow(mat);
831 inv = isl_mat_identity(mat->ctx, mat->n_col);
832 inv = isl_mat_cow(inv);
838 for (row = 0; row < mat->n_row; ++row) {
839 int pivot, first, i, off;
840 pivot = isl_seq_abs_min_non_zero(mat->row[row]+row, mat->n_col-row);
844 isl_assert(mat->ctx, pivot >= 0, goto error);
848 exchange(mat, &inv, NULL, row, pivot, row);
849 if (isl_int_is_neg(mat->row[row][row]))
850 oppose(mat, &inv, NULL, row, row);
852 while ((off = isl_seq_first_non_zero(mat->row[row]+first,
853 mat->n_col-first)) != -1) {
855 isl_int_fdiv_q(a, mat->row[row][first],
857 subtract(mat, &inv, NULL, row, row, first, a);
858 if (!isl_int_is_zero(mat->row[row][first]))
859 exchange(mat, &inv, NULL, row, row, first);
863 for (i = 0; i < row; ++i) {
864 if (isl_int_is_zero(mat->row[row][i]))
866 isl_int_gcd(a, mat->row[row][row], mat->row[row][i]);
867 isl_int_divexact(b, mat->row[row][i], a);
868 isl_int_divexact(a, mat->row[row][row], a);
870 isl_mat_col_combine(mat, i, a, i, b, row);
871 isl_mat_col_combine(inv, i, a, i, b, row);
876 isl_int_set(a, mat->row[0][0]);
877 for (row = 1; row < mat->n_row; ++row)
878 isl_int_lcm(a, a, mat->row[row][row]);
879 if (isl_int_is_zero(a)){
883 for (row = 0; row < mat->n_row; ++row) {
884 isl_int_divexact(mat->row[row][row], a, mat->row[row][row]);
885 if (isl_int_is_one(mat->row[row][row]))
887 isl_mat_col_scale(inv, row, mat->row[row][row]);
900 struct isl_mat *isl_mat_transpose(struct isl_mat *mat)
902 struct isl_mat *transpose = NULL;
905 if (mat->n_col == mat->n_row) {
906 mat = isl_mat_cow(mat);
909 for (i = 0; i < mat->n_row; ++i)
910 for (j = i + 1; j < mat->n_col; ++j)
911 isl_int_swap(mat->row[i][j], mat->row[j][i]);
914 transpose = isl_mat_alloc(mat->ctx, mat->n_col, mat->n_row);
917 for (i = 0; i < mat->n_row; ++i)
918 for (j = 0; j < mat->n_col; ++j)
919 isl_int_set(transpose->row[j][i], mat->row[i][j]);
927 struct isl_mat *isl_mat_swap_cols(struct isl_mat *mat, unsigned i, unsigned j)
931 mat = isl_mat_cow(mat);
934 isl_assert(mat->ctx, i < mat->n_col, goto error);
935 isl_assert(mat->ctx, j < mat->n_col, goto error);
937 for (r = 0; r < mat->n_row; ++r)
938 isl_int_swap(mat->row[r][i], mat->row[r][j]);
945 struct isl_mat *isl_mat_swap_rows(struct isl_mat *mat, unsigned i, unsigned j)
951 mat = isl_mat_cow(mat);
955 mat->row[i] = mat->row[j];
960 struct isl_mat *isl_mat_product(struct isl_mat *left, struct isl_mat *right)
963 struct isl_mat *prod;
967 isl_assert(left->ctx, left->n_col == right->n_row, goto error);
968 prod = isl_mat_alloc(left->ctx, left->n_row, right->n_col);
971 if (left->n_col == 0) {
972 for (i = 0; i < prod->n_row; ++i)
973 isl_seq_clr(prod->row[i], prod->n_col);
978 for (i = 0; i < prod->n_row; ++i) {
979 for (j = 0; j < prod->n_col; ++j) {
980 isl_int_mul(prod->row[i][j],
981 left->row[i][0], right->row[0][j]);
982 for (k = 1; k < left->n_col; ++k)
983 isl_int_addmul(prod->row[i][j],
984 left->row[i][k], right->row[k][j]);
996 /* Replace the variables x in the rows q by x' given by x = M x',
997 * with M the matrix mat.
999 * If the number of new variables is greater than the original
1000 * number of variables, then the rows q have already been
1001 * preextended. If the new number is smaller, then the coefficients
1002 * of the divs, which are not changed, need to be shifted down.
1003 * The row q may be the equalities, the inequalities or the
1004 * div expressions. In the latter case, has_div is true and
1005 * we need to take into account the extra denominator column.
1007 static int preimage(struct isl_ctx *ctx, isl_int **q, unsigned n,
1008 unsigned n_div, int has_div, struct isl_mat *mat)
1014 if (mat->n_col >= mat->n_row)
1017 e = mat->n_row - mat->n_col;
1019 for (i = 0; i < n; ++i)
1020 isl_int_mul(q[i][0], q[i][0], mat->row[0][0]);
1021 t = isl_mat_sub_alloc6(mat->ctx, q, 0, n, has_div, mat->n_row);
1022 t = isl_mat_product(t, mat);
1025 for (i = 0; i < n; ++i) {
1026 isl_seq_swp_or_cpy(q[i] + has_div, t->row[i], t->n_col);
1027 isl_seq_cpy(q[i] + has_div + t->n_col,
1028 q[i] + has_div + t->n_col + e, n_div);
1029 isl_seq_clr(q[i] + has_div + t->n_col + n_div, e);
1035 /* Replace the variables x in bset by x' given by x = M x', with
1038 * If there are fewer variables x' then there are x, then we perform
1039 * the transformation in place, which that, in principle,
1040 * this frees up some extra variables as the number
1041 * of columns remains constant, but we would have to extend
1042 * the div array too as the number of rows in this array is assumed
1043 * to be equal to extra.
1045 struct isl_basic_set *isl_basic_set_preimage(struct isl_basic_set *bset,
1046 struct isl_mat *mat)
1048 struct isl_ctx *ctx;
1054 bset = isl_basic_set_cow(bset);
1058 isl_assert(ctx, bset->dim->nparam == 0, goto error);
1059 isl_assert(ctx, 1+bset->dim->n_out == mat->n_row, goto error);
1060 isl_assert(ctx, mat->n_col > 0, goto error);
1062 if (mat->n_col > mat->n_row) {
1063 bset = isl_basic_set_extend(bset, 0, mat->n_col-1, 0, 0, 0);
1066 } else if (mat->n_col < mat->n_row) {
1067 bset->dim = isl_dim_cow(bset->dim);
1070 bset->dim->n_out -= mat->n_row - mat->n_col;
1073 if (preimage(ctx, bset->eq, bset->n_eq, bset->n_div, 0,
1074 isl_mat_copy(mat)) < 0)
1077 if (preimage(ctx, bset->ineq, bset->n_ineq, bset->n_div, 0,
1078 isl_mat_copy(mat)) < 0)
1081 if (preimage(ctx, bset->div, bset->n_div, bset->n_div, 1, mat) < 0)
1084 ISL_F_CLR(bset, ISL_BASIC_SET_NO_IMPLICIT);
1085 ISL_F_CLR(bset, ISL_BASIC_SET_NO_REDUNDANT);
1086 ISL_F_CLR(bset, ISL_BASIC_SET_NORMALIZED);
1087 ISL_F_CLR(bset, ISL_BASIC_SET_NORMALIZED_DIVS);
1088 ISL_F_CLR(bset, ISL_BASIC_SET_ALL_EQUALITIES);
1090 bset = isl_basic_set_simplify(bset);
1091 bset = isl_basic_set_finalize(bset);
1097 isl_basic_set_free(bset);
1101 struct isl_set *isl_set_preimage(struct isl_set *set, struct isl_mat *mat)
1103 struct isl_ctx *ctx;
1106 set = isl_set_cow(set);
1111 for (i = 0; i < set->n; ++i) {
1112 set->p[i] = isl_basic_set_preimage(set->p[i],
1117 if (mat->n_col != mat->n_row) {
1118 set->dim = isl_dim_cow(set->dim);
1121 set->dim->n_out += mat->n_col;
1122 set->dim->n_out -= mat->n_row;
1125 ISL_F_CLR(set, ISL_SET_NORMALIZED);
1133 /* Replace the variables x starting at pos in the rows q
1134 * by x' with x = M x' with M the matrix mat.
1135 * That is, replace the corresponding coefficients c by c M.
1137 static int transform(isl_ctx *ctx, isl_int **q, unsigned n,
1138 unsigned pos, __isl_take isl_mat *mat)
1143 t = isl_mat_sub_alloc6(ctx, q, 0, n, pos, mat->n_row);
1144 t = isl_mat_product(t, mat);
1147 for (i = 0; i < n; ++i)
1148 isl_seq_swp_or_cpy(q[i] + pos, t->row[i], t->n_col);
1153 /* Replace the variables x of type "type" starting at "first" in "bset"
1154 * by x' with x = M x' with M the matrix trans.
1155 * That is, replace the corresponding coefficients c by c M.
1157 * The transformation matrix should be a square matrix.
1159 __isl_give isl_basic_set *isl_basic_set_transform_dims(
1160 __isl_take isl_basic_set *bset, enum isl_dim_type type, unsigned first,
1161 __isl_take isl_mat *trans)
1166 bset = isl_basic_set_cow(bset);
1167 if (!bset || !trans)
1170 ctx = isl_basic_set_get_ctx(bset);
1171 if (trans->n_row != trans->n_col)
1172 isl_die(trans->ctx, isl_error_invalid,
1173 "expecting square transformation matrix", goto error);
1174 if (first + trans->n_row > isl_basic_set_dim(bset, type))
1175 isl_die(trans->ctx, isl_error_invalid,
1176 "oversized transformation matrix", goto error);
1178 pos = isl_basic_set_offset(bset, type) + first;
1180 if (transform(ctx, bset->eq, bset->n_eq, pos, isl_mat_copy(trans)) < 0)
1182 if (transform(ctx, bset->ineq, bset->n_ineq, pos,
1183 isl_mat_copy(trans)) < 0)
1185 if (transform(ctx, bset->div, bset->n_div, 1 + pos,
1186 isl_mat_copy(trans)) < 0)
1189 ISL_F_CLR(bset, ISL_BASIC_SET_NORMALIZED);
1190 ISL_F_CLR(bset, ISL_BASIC_SET_NORMALIZED_DIVS);
1192 isl_mat_free(trans);
1195 isl_mat_free(trans);
1196 isl_basic_set_free(bset);
1200 void isl_mat_print_internal(__isl_keep isl_mat *mat, FILE *out, int indent)
1205 fprintf(out, "%*snull mat\n", indent, "");
1209 if (mat->n_row == 0)
1210 fprintf(out, "%*s[]\n", indent, "");
1212 for (i = 0; i < mat->n_row; ++i) {
1214 fprintf(out, "%*s[[", indent, "");
1216 fprintf(out, "%*s[", indent+1, "");
1217 for (j = 0; j < mat->n_col; ++j) {
1220 isl_int_print(out, mat->row[i][j], 0);
1222 if (i == mat->n_row-1)
1223 fprintf(out, "]]\n");
1225 fprintf(out, "]\n");
1229 void isl_mat_dump(__isl_keep isl_mat *mat)
1231 isl_mat_print_internal(mat, stderr, 0);
1234 struct isl_mat *isl_mat_drop_cols(struct isl_mat *mat, unsigned col, unsigned n)
1238 mat = isl_mat_cow(mat);
1242 if (col != mat->n_col-n) {
1243 for (r = 0; r < mat->n_row; ++r)
1244 isl_seq_cpy(mat->row[r]+col, mat->row[r]+col+n,
1245 mat->n_col - col - n);
1251 struct isl_mat *isl_mat_drop_rows(struct isl_mat *mat, unsigned row, unsigned n)
1255 mat = isl_mat_cow(mat);
1259 for (r = row; r+n < mat->n_row; ++r)
1260 mat->row[r] = mat->row[r+n];
1266 __isl_give isl_mat *isl_mat_insert_cols(__isl_take isl_mat *mat,
1267 unsigned col, unsigned n)
1276 ext = isl_mat_alloc(mat->ctx, mat->n_row, mat->n_col + n);
1280 isl_mat_sub_copy(mat->ctx, ext->row, mat->row, mat->n_row, 0, 0, col);
1281 isl_mat_sub_copy(mat->ctx, ext->row, mat->row, mat->n_row,
1282 col + n, col, mat->n_col - col);
1291 __isl_give isl_mat *isl_mat_insert_zero_cols(__isl_take isl_mat *mat,
1292 unsigned first, unsigned n)
1298 mat = isl_mat_insert_cols(mat, first, n);
1302 for (i = 0; i < mat->n_row; ++i)
1303 isl_seq_clr(mat->row[i] + first, n);
1308 __isl_give isl_mat *isl_mat_add_zero_cols(__isl_take isl_mat *mat, unsigned n)
1313 return isl_mat_insert_zero_cols(mat, mat->n_col, n);
1316 __isl_give isl_mat *isl_mat_insert_rows(__isl_take isl_mat *mat,
1317 unsigned row, unsigned n)
1326 ext = isl_mat_alloc(mat->ctx, mat->n_row + n, mat->n_col);
1330 isl_mat_sub_copy(mat->ctx, ext->row, mat->row, row, 0, 0, mat->n_col);
1331 isl_mat_sub_copy(mat->ctx, ext->row + row + n, mat->row + row,
1332 mat->n_row - row, 0, 0, mat->n_col);
1341 __isl_give isl_mat *isl_mat_add_rows(__isl_take isl_mat *mat, unsigned n)
1346 return isl_mat_insert_rows(mat, mat->n_row, n);
1349 __isl_give isl_mat *isl_mat_insert_zero_rows(__isl_take isl_mat *mat,
1350 unsigned row, unsigned n)
1354 mat = isl_mat_insert_rows(mat, row, n);
1358 for (i = 0; i < n; ++i)
1359 isl_seq_clr(mat->row[row + i], mat->n_col);
1364 __isl_give isl_mat *isl_mat_add_zero_rows(__isl_take isl_mat *mat, unsigned n)
1369 return isl_mat_insert_zero_rows(mat, mat->n_row, n);
1372 void isl_mat_col_submul(struct isl_mat *mat,
1373 int dst_col, isl_int f, int src_col)
1377 for (i = 0; i < mat->n_row; ++i)
1378 isl_int_submul(mat->row[i][dst_col], f, mat->row[i][src_col]);
1381 void isl_mat_col_add(__isl_keep isl_mat *mat, int dst_col, int src_col)
1388 for (i = 0; i < mat->n_row; ++i)
1389 isl_int_add(mat->row[i][dst_col],
1390 mat->row[i][dst_col], mat->row[i][src_col]);
1393 void isl_mat_col_mul(struct isl_mat *mat, int dst_col, isl_int f, int src_col)
1397 for (i = 0; i < mat->n_row; ++i)
1398 isl_int_mul(mat->row[i][dst_col], f, mat->row[i][src_col]);
1401 struct isl_mat *isl_mat_unimodular_complete(struct isl_mat *M, int row)
1404 struct isl_mat *H = NULL, *Q = NULL;
1409 isl_assert(M->ctx, M->n_row == M->n_col, goto error);
1411 H = isl_mat_left_hermite(isl_mat_copy(M), 0, NULL, &Q);
1412 M->n_row = M->n_col;
1415 for (r = 0; r < row; ++r)
1416 isl_assert(M->ctx, isl_int_is_one(H->row[r][r]), goto error);
1417 for (r = row; r < M->n_row; ++r)
1418 isl_seq_cpy(M->row[r], Q->row[r], M->n_col);
1429 __isl_give isl_mat *isl_mat_concat(__isl_take isl_mat *top,
1430 __isl_take isl_mat *bot)
1432 struct isl_mat *mat;
1437 isl_assert(top->ctx, top->n_col == bot->n_col, goto error);
1438 if (top->n_row == 0) {
1442 if (bot->n_row == 0) {
1447 mat = isl_mat_alloc(top->ctx, top->n_row + bot->n_row, top->n_col);
1450 isl_mat_sub_copy(mat->ctx, mat->row, top->row, top->n_row,
1452 isl_mat_sub_copy(mat->ctx, mat->row + top->n_row, bot->row, bot->n_row,
1463 int isl_mat_is_equal(__isl_keep isl_mat *mat1, __isl_keep isl_mat *mat2)
1470 if (mat1->n_row != mat2->n_row)
1473 if (mat1->n_col != mat2->n_col)
1476 for (i = 0; i < mat1->n_row; ++i)
1477 if (!isl_seq_eq(mat1->row[i], mat2->row[i], mat1->n_col))
1483 __isl_give isl_mat *isl_mat_from_row_vec(__isl_take isl_vec *vec)
1485 struct isl_mat *mat;
1489 mat = isl_mat_alloc(vec->ctx, 1, vec->size);
1493 isl_seq_cpy(mat->row[0], vec->el, vec->size);
1502 __isl_give isl_mat *isl_mat_vec_concat(__isl_take isl_mat *top,
1503 __isl_take isl_vec *bot)
1505 return isl_mat_concat(top, isl_mat_from_row_vec(bot));
1508 __isl_give isl_mat *isl_mat_move_cols(__isl_take isl_mat *mat,
1509 unsigned dst_col, unsigned src_col, unsigned n)
1515 if (n == 0 || dst_col == src_col)
1518 res = isl_mat_alloc(mat->ctx, mat->n_row, mat->n_col);
1522 if (dst_col < src_col) {
1523 isl_mat_sub_copy(res->ctx, res->row, mat->row, mat->n_row,
1525 isl_mat_sub_copy(res->ctx, res->row, mat->row, mat->n_row,
1526 dst_col, src_col, n);
1527 isl_mat_sub_copy(res->ctx, res->row, mat->row, mat->n_row,
1528 dst_col + n, dst_col, src_col - dst_col);
1529 isl_mat_sub_copy(res->ctx, res->row, mat->row, mat->n_row,
1530 src_col + n, src_col + n,
1531 res->n_col - src_col - n);
1533 isl_mat_sub_copy(res->ctx, res->row, mat->row, mat->n_row,
1535 isl_mat_sub_copy(res->ctx, res->row, mat->row, mat->n_row,
1536 src_col, src_col + n, dst_col - src_col);
1537 isl_mat_sub_copy(res->ctx, res->row, mat->row, mat->n_row,
1538 dst_col, src_col, n);
1539 isl_mat_sub_copy(res->ctx, res->row, mat->row, mat->n_row,
1540 dst_col + n, dst_col + n,
1541 res->n_col - dst_col - n);
1551 void isl_mat_gcd(__isl_keep isl_mat *mat, isl_int *gcd)
1556 isl_int_set_si(*gcd, 0);
1561 for (i = 0; i < mat->n_row; ++i) {
1562 isl_seq_gcd(mat->row[i], mat->n_col, &g);
1563 isl_int_gcd(*gcd, *gcd, g);
1568 __isl_give isl_mat *isl_mat_scale_down(__isl_take isl_mat *mat, isl_int m)
1572 if (isl_int_is_one(m))
1575 mat = isl_mat_cow(mat);
1579 for (i = 0; i < mat->n_row; ++i)
1580 isl_seq_scale_down(mat->row[i], mat->row[i], m, mat->n_col);
1585 __isl_give isl_mat *isl_mat_normalize(__isl_take isl_mat *mat)
1593 isl_mat_gcd(mat, &gcd);
1594 mat = isl_mat_scale_down(mat, gcd);
1600 __isl_give isl_mat *isl_mat_normalize_row(__isl_take isl_mat *mat, int row)
1602 mat = isl_mat_cow(mat);
1606 isl_seq_normalize(mat->ctx, mat->row[row], mat->n_col);
1611 /* Number of initial non-zero columns.
1613 int isl_mat_initial_non_zero_cols(__isl_keep isl_mat *mat)
1620 for (i = 0; i < mat->n_col; ++i)
1621 if (row_first_non_zero(mat->row, mat->n_row, i) < 0)