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
12 #include <isl_mat_private.h>
13 #include "isl_map_private.h"
14 #include <isl_dim_private.h>
16 struct isl_mat *isl_mat_alloc(struct isl_ctx *ctx,
17 unsigned n_row, unsigned n_col)
22 mat = isl_alloc_type(ctx, struct isl_mat);
27 mat->block = isl_blk_alloc(ctx, n_row * n_col);
28 if (isl_blk_is_error(mat->block))
30 mat->row = isl_alloc_array(ctx, isl_int *, n_row);
34 for (i = 0; i < n_row; ++i)
35 mat->row[i] = mat->block.data + i * n_col;
47 isl_blk_free(ctx, mat->block);
52 struct isl_mat *isl_mat_extend(struct isl_mat *mat,
53 unsigned n_row, unsigned n_col)
61 if (mat->max_col >= n_col && mat->n_row >= n_row) {
62 if (mat->n_col < n_col)
67 if (mat->max_col < n_col) {
68 struct isl_mat *new_mat;
70 if (n_row < mat->n_row)
72 new_mat = isl_mat_alloc(mat->ctx, n_row, n_col);
75 for (i = 0; i < mat->n_row; ++i)
76 isl_seq_cpy(new_mat->row[i], mat->row[i], mat->n_col);
81 mat = isl_mat_cow(mat);
85 old = mat->block.data;
86 mat->block = isl_blk_extend(mat->ctx, mat->block, n_row * mat->max_col);
87 if (isl_blk_is_error(mat->block))
89 mat->row = isl_realloc_array(mat->ctx, mat->row, isl_int *, n_row);
93 for (i = 0; i < mat->n_row; ++i)
94 mat->row[i] = mat->block.data + (mat->row[i] - old);
95 for (i = mat->n_row; i < n_row; ++i)
96 mat->row[i] = mat->block.data + i * mat->max_col;
98 if (mat->n_col < n_col)
107 struct isl_mat *isl_mat_sub_alloc(struct isl_ctx *ctx, isl_int **row,
108 unsigned first_row, unsigned n_row, unsigned first_col, unsigned n_col)
113 mat = isl_alloc_type(ctx, struct isl_mat);
116 mat->row = isl_alloc_array(ctx, isl_int *, n_row);
119 for (i = 0; i < n_row; ++i)
120 mat->row[i] = row[first_row+i] + first_col;
126 mat->block = isl_blk_empty();
127 mat->flags = ISL_MAT_BORROWED;
134 void isl_mat_sub_copy(struct isl_ctx *ctx, isl_int **dst, isl_int **src,
135 unsigned n_row, unsigned dst_col, unsigned src_col, unsigned n_col)
139 for (i = 0; i < n_row; ++i)
140 isl_seq_cpy(dst[i]+dst_col, src[i]+src_col, n_col);
143 void isl_mat_sub_neg(struct isl_ctx *ctx, isl_int **dst, isl_int **src,
144 unsigned n_row, unsigned dst_col, unsigned src_col, unsigned n_col)
148 for (i = 0; i < n_row; ++i)
149 isl_seq_neg(dst[i]+dst_col, src[i]+src_col, n_col);
152 struct isl_mat *isl_mat_copy(struct isl_mat *mat)
161 struct isl_mat *isl_mat_dup(struct isl_mat *mat)
164 struct isl_mat *mat2;
168 mat2 = isl_mat_alloc(mat->ctx, mat->n_row, mat->n_col);
171 for (i = 0; i < mat->n_row; ++i)
172 isl_seq_cpy(mat2->row[i], mat->row[i], mat->n_col);
176 struct isl_mat *isl_mat_cow(struct isl_mat *mat)
178 struct isl_mat *mat2;
182 if (mat->ref == 1 && !ISL_F_ISSET(mat, ISL_MAT_BORROWED))
185 mat2 = isl_mat_dup(mat);
190 void isl_mat_free(struct isl_mat *mat)
198 if (!ISL_F_ISSET(mat, ISL_MAT_BORROWED))
199 isl_blk_free(mat->ctx, mat->block);
200 isl_ctx_deref(mat->ctx);
205 int isl_mat_rows(__isl_keep isl_mat *mat)
207 return mat ? mat->n_row : -1;
210 int isl_mat_cols(__isl_keep isl_mat *mat)
212 return mat ? mat->n_col : -1;
215 int isl_mat_get_element(__isl_keep isl_mat *mat, int row, int col, isl_int *v)
219 if (row < 0 || row >= mat->n_row)
220 isl_die(mat->ctx, isl_error_invalid, "row out of range",
222 if (col < 0 || col >= mat->n_col)
223 isl_die(mat->ctx, isl_error_invalid, "column out of range",
225 isl_int_set(*v, mat->row[row][col]);
229 __isl_give isl_mat *isl_mat_set_element(__isl_take isl_mat *mat,
230 int row, int col, isl_int v)
232 mat = isl_mat_cow(mat);
235 if (row < 0 || row >= mat->n_row)
236 isl_die(mat->ctx, isl_error_invalid, "row out of range",
238 if (col < 0 || col >= mat->n_col)
239 isl_die(mat->ctx, isl_error_invalid, "column out of range",
241 isl_int_set(mat->row[row][col], v);
248 struct isl_mat *isl_mat_identity(struct isl_ctx *ctx, unsigned n_row)
253 mat = isl_mat_alloc(ctx, n_row, n_row);
256 for (i = 0; i < n_row; ++i) {
257 isl_seq_clr(mat->row[i], i);
258 isl_int_set_si(mat->row[i][i], 1);
259 isl_seq_clr(mat->row[i]+i+1, n_row-(i+1));
265 struct isl_vec *isl_mat_vec_product(struct isl_mat *mat, struct isl_vec *vec)
268 struct isl_vec *prod;
273 isl_assert(mat->ctx, mat->n_col == vec->size, goto error);
275 prod = isl_vec_alloc(mat->ctx, mat->n_row);
279 for (i = 0; i < prod->size; ++i)
280 isl_seq_inner_product(mat->row[i], vec->el, vec->size,
281 &prod->block.data[i]);
291 __isl_give isl_vec *isl_mat_vec_inverse_product(__isl_take isl_mat *mat,
292 __isl_take isl_vec *vec)
294 struct isl_mat *vec_mat;
299 vec_mat = isl_mat_alloc(vec->ctx, vec->size, 1);
302 for (i = 0; i < vec->size; ++i)
303 isl_int_set(vec_mat->row[i][0], vec->el[i]);
304 vec_mat = isl_mat_inverse_product(mat, vec_mat);
308 vec = isl_vec_alloc(vec_mat->ctx, vec_mat->n_row);
310 for (i = 0; i < vec->size; ++i)
311 isl_int_set(vec->el[i], vec_mat->row[i][0]);
312 isl_mat_free(vec_mat);
320 struct isl_vec *isl_vec_mat_product(struct isl_vec *vec, struct isl_mat *mat)
323 struct isl_vec *prod;
328 isl_assert(mat->ctx, mat->n_row == vec->size, goto error);
330 prod = isl_vec_alloc(mat->ctx, mat->n_col);
334 for (i = 0; i < prod->size; ++i) {
335 isl_int_set_si(prod->el[i], 0);
336 for (j = 0; j < vec->size; ++j)
337 isl_int_addmul(prod->el[i], vec->el[j], mat->row[j][i]);
348 struct isl_mat *isl_mat_aff_direct_sum(struct isl_mat *left,
349 struct isl_mat *right)
357 isl_assert(left->ctx, left->n_row == right->n_row, goto error);
358 isl_assert(left->ctx, left->n_row >= 1, goto error);
359 isl_assert(left->ctx, left->n_col >= 1, goto error);
360 isl_assert(left->ctx, right->n_col >= 1, goto error);
361 isl_assert(left->ctx,
362 isl_seq_first_non_zero(left->row[0]+1, left->n_col-1) == -1,
364 isl_assert(left->ctx,
365 isl_seq_first_non_zero(right->row[0]+1, right->n_col-1) == -1,
368 sum = isl_mat_alloc(left->ctx, left->n_row, left->n_col + right->n_col - 1);
371 isl_int_lcm(sum->row[0][0], left->row[0][0], right->row[0][0]);
372 isl_int_divexact(left->row[0][0], sum->row[0][0], left->row[0][0]);
373 isl_int_divexact(right->row[0][0], sum->row[0][0], right->row[0][0]);
375 isl_seq_clr(sum->row[0]+1, sum->n_col-1);
376 for (i = 1; i < sum->n_row; ++i) {
377 isl_int_mul(sum->row[i][0], left->row[0][0], left->row[i][0]);
378 isl_int_addmul(sum->row[i][0],
379 right->row[0][0], right->row[i][0]);
380 isl_seq_scale(sum->row[i]+1, left->row[i]+1, left->row[0][0],
382 isl_seq_scale(sum->row[i]+left->n_col,
383 right->row[i]+1, right->row[0][0],
387 isl_int_divexact(left->row[0][0], sum->row[0][0], left->row[0][0]);
388 isl_int_divexact(right->row[0][0], sum->row[0][0], right->row[0][0]);
398 static void exchange(struct isl_mat *M, struct isl_mat **U,
399 struct isl_mat **Q, unsigned row, unsigned i, unsigned j)
402 for (r = row; r < M->n_row; ++r)
403 isl_int_swap(M->row[r][i], M->row[r][j]);
405 for (r = 0; r < (*U)->n_row; ++r)
406 isl_int_swap((*U)->row[r][i], (*U)->row[r][j]);
409 isl_mat_swap_rows(*Q, i, j);
412 static void subtract(struct isl_mat *M, struct isl_mat **U,
413 struct isl_mat **Q, unsigned row, unsigned i, unsigned j, isl_int m)
416 for (r = row; r < M->n_row; ++r)
417 isl_int_submul(M->row[r][j], m, M->row[r][i]);
419 for (r = 0; r < (*U)->n_row; ++r)
420 isl_int_submul((*U)->row[r][j], m, (*U)->row[r][i]);
423 for (r = 0; r < (*Q)->n_col; ++r)
424 isl_int_addmul((*Q)->row[i][r], m, (*Q)->row[j][r]);
428 static void oppose(struct isl_mat *M, struct isl_mat **U,
429 struct isl_mat **Q, unsigned row, unsigned col)
432 for (r = row; r < M->n_row; ++r)
433 isl_int_neg(M->row[r][col], M->row[r][col]);
435 for (r = 0; r < (*U)->n_row; ++r)
436 isl_int_neg((*U)->row[r][col], (*U)->row[r][col]);
439 isl_seq_neg((*Q)->row[col], (*Q)->row[col], (*Q)->n_col);
442 /* Given matrix M, compute
447 * with U and Q unimodular matrices and H a matrix in column echelon form
448 * such that on each echelon row the entries in the non-echelon column
449 * are non-negative (if neg == 0) or non-positive (if neg == 1)
450 * and stricly smaller (in absolute value) than the entries in the echelon
452 * If U or Q are NULL, then these matrices are not computed.
454 struct isl_mat *isl_mat_left_hermite(struct isl_mat *M, int neg,
455 struct isl_mat **U, struct isl_mat **Q)
470 *U = isl_mat_identity(M->ctx, M->n_col);
475 *Q = isl_mat_identity(M->ctx, M->n_col);
482 for (row = 0; row < M->n_row; ++row) {
484 first = isl_seq_abs_min_non_zero(M->row[row]+col, M->n_col-col);
489 exchange(M, U, Q, row, first, col);
490 if (isl_int_is_neg(M->row[row][col]))
491 oppose(M, U, Q, row, col);
493 while ((off = isl_seq_first_non_zero(M->row[row]+first,
494 M->n_col-first)) != -1) {
496 isl_int_fdiv_q(c, M->row[row][first], M->row[row][col]);
497 subtract(M, U, Q, row, col, first, c);
498 if (!isl_int_is_zero(M->row[row][first]))
499 exchange(M, U, Q, row, first, col);
503 for (i = 0; i < col; ++i) {
504 if (isl_int_is_zero(M->row[row][i]))
507 isl_int_cdiv_q(c, M->row[row][i], M->row[row][col]);
509 isl_int_fdiv_q(c, M->row[row][i], M->row[row][col]);
510 if (isl_int_is_zero(c))
512 subtract(M, U, Q, row, col, i, c);
531 struct isl_mat *isl_mat_right_kernel(struct isl_mat *mat)
534 struct isl_mat *U = NULL;
537 mat = isl_mat_left_hermite(mat, 0, &U, NULL);
541 for (i = 0, rank = 0; rank < mat->n_col; ++rank) {
542 while (i < mat->n_row && isl_int_is_zero(mat->row[i][rank]))
547 K = isl_mat_alloc(U->ctx, U->n_row, U->n_col - rank);
550 isl_mat_sub_copy(K->ctx, K->row, U->row, U->n_row, 0, rank, U->n_col-rank);
560 struct isl_mat *isl_mat_lin_to_aff(struct isl_mat *mat)
563 struct isl_mat *mat2;
567 mat2 = isl_mat_alloc(mat->ctx, 1+mat->n_row, 1+mat->n_col);
570 isl_int_set_si(mat2->row[0][0], 1);
571 isl_seq_clr(mat2->row[0]+1, mat->n_col);
572 for (i = 0; i < mat->n_row; ++i) {
573 isl_int_set_si(mat2->row[1+i][0], 0);
574 isl_seq_cpy(mat2->row[1+i]+1, mat->row[i], mat->n_col);
583 /* Given two matrices M1 and M2, return the block matrix
588 __isl_give isl_mat *isl_mat_diagonal(__isl_take isl_mat *mat1,
589 __isl_take isl_mat *mat2)
597 mat = isl_mat_alloc(mat1->ctx, mat1->n_row + mat2->n_row,
598 mat1->n_col + mat2->n_col);
601 for (i = 0; i < mat1->n_row; ++i) {
602 isl_seq_cpy(mat->row[i], mat1->row[i], mat1->n_col);
603 isl_seq_clr(mat->row[i] + mat1->n_col, mat2->n_col);
605 for (i = 0; i < mat2->n_row; ++i) {
606 isl_seq_clr(mat->row[mat1->n_row + i], mat1->n_col);
607 isl_seq_cpy(mat->row[mat1->n_row + i] + mat1->n_col,
608 mat2->row[i], mat2->n_col);
619 static int row_first_non_zero(isl_int **row, unsigned n_row, unsigned col)
623 for (i = 0; i < n_row; ++i)
624 if (!isl_int_is_zero(row[i][col]))
629 static int row_abs_min_non_zero(isl_int **row, unsigned n_row, unsigned col)
631 int i, min = row_first_non_zero(row, n_row, col);
634 for (i = min + 1; i < n_row; ++i) {
635 if (isl_int_is_zero(row[i][col]))
637 if (isl_int_abs_lt(row[i][col], row[min][col]))
643 static void inv_exchange(struct isl_mat *left, struct isl_mat *right,
644 unsigned i, unsigned j)
646 left = isl_mat_swap_rows(left, i, j);
647 right = isl_mat_swap_rows(right, i, j);
650 static void inv_oppose(
651 struct isl_mat *left, struct isl_mat *right, unsigned row)
653 isl_seq_neg(left->row[row]+row, left->row[row]+row, left->n_col-row);
654 isl_seq_neg(right->row[row], right->row[row], right->n_col);
657 static void inv_subtract(struct isl_mat *left, struct isl_mat *right,
658 unsigned row, unsigned i, isl_int m)
661 isl_seq_combine(left->row[i]+row,
662 left->ctx->one, left->row[i]+row,
663 m, left->row[row]+row,
665 isl_seq_combine(right->row[i], right->ctx->one, right->row[i],
666 m, right->row[row], right->n_col);
669 /* Compute inv(left)*right
671 struct isl_mat *isl_mat_inverse_product(struct isl_mat *left,
672 struct isl_mat *right)
680 isl_assert(left->ctx, left->n_row == left->n_col, goto error);
681 isl_assert(left->ctx, left->n_row == right->n_row, goto error);
683 if (left->n_row == 0) {
688 left = isl_mat_cow(left);
689 right = isl_mat_cow(right);
695 for (row = 0; row < left->n_row; ++row) {
696 int pivot, first, i, off;
697 pivot = row_abs_min_non_zero(left->row+row, left->n_row-row, row);
701 isl_assert(left->ctx, pivot >= 0, goto error);
705 inv_exchange(left, right, pivot, row);
706 if (isl_int_is_neg(left->row[row][row]))
707 inv_oppose(left, right, row);
709 while ((off = row_first_non_zero(left->row+first,
710 left->n_row-first, row)) != -1) {
712 isl_int_fdiv_q(a, left->row[first][row],
713 left->row[row][row]);
714 inv_subtract(left, right, row, first, a);
715 if (!isl_int_is_zero(left->row[first][row]))
716 inv_exchange(left, right, row, first);
720 for (i = 0; i < row; ++i) {
721 if (isl_int_is_zero(left->row[i][row]))
723 isl_int_gcd(a, left->row[row][row], left->row[i][row]);
724 isl_int_divexact(b, left->row[i][row], a);
725 isl_int_divexact(a, left->row[row][row], a);
727 isl_seq_combine(left->row[i] + i,
729 b, left->row[row] + i,
731 isl_seq_combine(right->row[i], a, right->row[i],
732 b, right->row[row], right->n_col);
737 isl_int_set(a, left->row[0][0]);
738 for (row = 1; row < left->n_row; ++row)
739 isl_int_lcm(a, a, left->row[row][row]);
740 if (isl_int_is_zero(a)){
742 isl_assert(left->ctx, 0, goto error);
744 for (row = 0; row < left->n_row; ++row) {
745 isl_int_divexact(left->row[row][row], a, left->row[row][row]);
746 if (isl_int_is_one(left->row[row][row]))
748 isl_seq_scale(right->row[row], right->row[row],
749 left->row[row][row], right->n_col);
761 void isl_mat_col_scale(struct isl_mat *mat, unsigned col, isl_int m)
765 for (i = 0; i < mat->n_row; ++i)
766 isl_int_mul(mat->row[i][col], mat->row[i][col], m);
769 void isl_mat_col_combine(struct isl_mat *mat, unsigned dst,
770 isl_int m1, unsigned src1, isl_int m2, unsigned src2)
776 for (i = 0; i < mat->n_row; ++i) {
777 isl_int_mul(tmp, m1, mat->row[i][src1]);
778 isl_int_addmul(tmp, m2, mat->row[i][src2]);
779 isl_int_set(mat->row[i][dst], tmp);
784 struct isl_mat *isl_mat_right_inverse(struct isl_mat *mat)
790 mat = isl_mat_cow(mat);
794 inv = isl_mat_identity(mat->ctx, mat->n_col);
795 inv = isl_mat_cow(inv);
801 for (row = 0; row < mat->n_row; ++row) {
802 int pivot, first, i, off;
803 pivot = isl_seq_abs_min_non_zero(mat->row[row]+row, mat->n_col-row);
807 isl_assert(mat->ctx, pivot >= 0, goto error);
811 exchange(mat, &inv, NULL, row, pivot, row);
812 if (isl_int_is_neg(mat->row[row][row]))
813 oppose(mat, &inv, NULL, row, row);
815 while ((off = isl_seq_first_non_zero(mat->row[row]+first,
816 mat->n_col-first)) != -1) {
818 isl_int_fdiv_q(a, mat->row[row][first],
820 subtract(mat, &inv, NULL, row, row, first, a);
821 if (!isl_int_is_zero(mat->row[row][first]))
822 exchange(mat, &inv, NULL, row, row, first);
826 for (i = 0; i < row; ++i) {
827 if (isl_int_is_zero(mat->row[row][i]))
829 isl_int_gcd(a, mat->row[row][row], mat->row[row][i]);
830 isl_int_divexact(b, mat->row[row][i], a);
831 isl_int_divexact(a, mat->row[row][row], a);
833 isl_mat_col_combine(mat, i, a, i, b, row);
834 isl_mat_col_combine(inv, i, a, i, b, row);
839 isl_int_set(a, mat->row[0][0]);
840 for (row = 1; row < mat->n_row; ++row)
841 isl_int_lcm(a, a, mat->row[row][row]);
842 if (isl_int_is_zero(a)){
846 for (row = 0; row < mat->n_row; ++row) {
847 isl_int_divexact(mat->row[row][row], a, mat->row[row][row]);
848 if (isl_int_is_one(mat->row[row][row]))
850 isl_mat_col_scale(inv, row, mat->row[row][row]);
863 struct isl_mat *isl_mat_transpose(struct isl_mat *mat)
865 struct isl_mat *transpose = NULL;
868 if (mat->n_col == mat->n_row) {
869 mat = isl_mat_cow(mat);
872 for (i = 0; i < mat->n_row; ++i)
873 for (j = i + 1; j < mat->n_col; ++j)
874 isl_int_swap(mat->row[i][j], mat->row[j][i]);
877 transpose = isl_mat_alloc(mat->ctx, mat->n_col, mat->n_row);
880 for (i = 0; i < mat->n_row; ++i)
881 for (j = 0; j < mat->n_col; ++j)
882 isl_int_set(transpose->row[j][i], mat->row[i][j]);
890 struct isl_mat *isl_mat_swap_cols(struct isl_mat *mat, unsigned i, unsigned j)
894 mat = isl_mat_cow(mat);
897 isl_assert(mat->ctx, i < mat->n_col, goto error);
898 isl_assert(mat->ctx, j < mat->n_col, goto error);
900 for (r = 0; r < mat->n_row; ++r)
901 isl_int_swap(mat->row[r][i], mat->row[r][j]);
908 struct isl_mat *isl_mat_swap_rows(struct isl_mat *mat, unsigned i, unsigned j)
914 mat = isl_mat_cow(mat);
918 mat->row[i] = mat->row[j];
923 struct isl_mat *isl_mat_product(struct isl_mat *left, struct isl_mat *right)
926 struct isl_mat *prod;
930 isl_assert(left->ctx, left->n_col == right->n_row, goto error);
931 prod = isl_mat_alloc(left->ctx, left->n_row, right->n_col);
934 if (left->n_col == 0) {
935 for (i = 0; i < prod->n_row; ++i)
936 isl_seq_clr(prod->row[i], prod->n_col);
939 for (i = 0; i < prod->n_row; ++i) {
940 for (j = 0; j < prod->n_col; ++j) {
941 isl_int_mul(prod->row[i][j],
942 left->row[i][0], right->row[0][j]);
943 for (k = 1; k < left->n_col; ++k)
944 isl_int_addmul(prod->row[i][j],
945 left->row[i][k], right->row[k][j]);
957 /* Replace the variables x in the rows q by x' given by x = M x',
958 * with M the matrix mat.
960 * If the number of new variables is greater than the original
961 * number of variables, then the rows q have already been
962 * preextended. If the new number is smaller, then the coefficients
963 * of the divs, which are not changed, need to be shifted down.
964 * The row q may be the equalities, the inequalities or the
965 * div expressions. In the latter case, has_div is true and
966 * we need to take into account the extra denominator column.
968 static int preimage(struct isl_ctx *ctx, isl_int **q, unsigned n,
969 unsigned n_div, int has_div, struct isl_mat *mat)
975 if (mat->n_col >= mat->n_row)
978 e = mat->n_row - mat->n_col;
980 for (i = 0; i < n; ++i)
981 isl_int_mul(q[i][0], q[i][0], mat->row[0][0]);
982 t = isl_mat_sub_alloc(mat->ctx, q, 0, n, has_div, mat->n_row);
983 t = isl_mat_product(t, mat);
986 for (i = 0; i < n; ++i) {
987 isl_seq_swp_or_cpy(q[i] + has_div, t->row[i], t->n_col);
988 isl_seq_cpy(q[i] + has_div + t->n_col,
989 q[i] + has_div + t->n_col + e, n_div);
990 isl_seq_clr(q[i] + has_div + t->n_col + n_div, e);
996 /* Replace the variables x in bset by x' given by x = M x', with
999 * If there are fewer variables x' then there are x, then we perform
1000 * the transformation in place, which that, in principle,
1001 * this frees up some extra variables as the number
1002 * of columns remains constant, but we would have to extend
1003 * the div array too as the number of rows in this array is assumed
1004 * to be equal to extra.
1006 struct isl_basic_set *isl_basic_set_preimage(struct isl_basic_set *bset,
1007 struct isl_mat *mat)
1009 struct isl_ctx *ctx;
1015 bset = isl_basic_set_cow(bset);
1019 isl_assert(ctx, bset->dim->nparam == 0, goto error);
1020 isl_assert(ctx, 1+bset->dim->n_out == mat->n_row, goto error);
1021 isl_assert(ctx, mat->n_col > 0, goto error);
1023 if (mat->n_col > mat->n_row) {
1024 bset = isl_basic_set_extend(bset, 0, mat->n_col-1, 0, 0, 0);
1027 } else if (mat->n_col < mat->n_row) {
1028 bset->dim = isl_dim_cow(bset->dim);
1031 bset->dim->n_out -= mat->n_row - mat->n_col;
1034 if (preimage(ctx, bset->eq, bset->n_eq, bset->n_div, 0,
1035 isl_mat_copy(mat)) < 0)
1038 if (preimage(ctx, bset->ineq, bset->n_ineq, bset->n_div, 0,
1039 isl_mat_copy(mat)) < 0)
1042 if (preimage(ctx, bset->div, bset->n_div, bset->n_div, 1, mat) < 0)
1045 ISL_F_CLR(bset, ISL_BASIC_SET_NO_IMPLICIT);
1046 ISL_F_CLR(bset, ISL_BASIC_SET_NO_REDUNDANT);
1047 ISL_F_CLR(bset, ISL_BASIC_SET_NORMALIZED);
1048 ISL_F_CLR(bset, ISL_BASIC_SET_NORMALIZED_DIVS);
1049 ISL_F_CLR(bset, ISL_BASIC_SET_ALL_EQUALITIES);
1051 bset = isl_basic_set_simplify(bset);
1052 bset = isl_basic_set_finalize(bset);
1058 isl_basic_set_free(bset);
1062 struct isl_set *isl_set_preimage(struct isl_set *set, struct isl_mat *mat)
1064 struct isl_ctx *ctx;
1067 set = isl_set_cow(set);
1072 for (i = 0; i < set->n; ++i) {
1073 set->p[i] = isl_basic_set_preimage(set->p[i],
1078 if (mat->n_col != mat->n_row) {
1079 set->dim = isl_dim_cow(set->dim);
1082 set->dim->n_out += mat->n_col;
1083 set->dim->n_out -= mat->n_row;
1086 ISL_F_CLR(set, ISL_SET_NORMALIZED);
1094 void isl_mat_dump(struct isl_mat *mat, FILE *out, int indent)
1099 fprintf(out, "%*snull mat\n", indent, "");
1103 if (mat->n_row == 0)
1104 fprintf(out, "%*s[]\n", indent, "");
1106 for (i = 0; i < mat->n_row; ++i) {
1108 fprintf(out, "%*s[[", indent, "");
1110 fprintf(out, "%*s[", indent+1, "");
1111 for (j = 0; j < mat->n_col; ++j) {
1114 isl_int_print(out, mat->row[i][j], 0);
1116 if (i == mat->n_row-1)
1117 fprintf(out, "]]\n");
1119 fprintf(out, "]\n");
1123 struct isl_mat *isl_mat_drop_cols(struct isl_mat *mat, unsigned col, unsigned n)
1127 mat = isl_mat_cow(mat);
1131 if (col != mat->n_col-n) {
1132 for (r = 0; r < mat->n_row; ++r)
1133 isl_seq_cpy(mat->row[r]+col, mat->row[r]+col+n,
1134 mat->n_col - col - n);
1140 struct isl_mat *isl_mat_drop_rows(struct isl_mat *mat, unsigned row, unsigned n)
1144 mat = isl_mat_cow(mat);
1148 for (r = row; r+n < mat->n_row; ++r)
1149 mat->row[r] = mat->row[r+n];
1155 __isl_give isl_mat *isl_mat_insert_cols(__isl_take isl_mat *mat,
1156 unsigned col, unsigned n)
1165 ext = isl_mat_alloc(mat->ctx, mat->n_row, mat->n_col + n);
1169 isl_mat_sub_copy(mat->ctx, ext->row, mat->row, mat->n_row, 0, 0, col);
1170 isl_mat_sub_copy(mat->ctx, ext->row, mat->row, mat->n_row,
1171 col + n, col, mat->n_col - col);
1180 __isl_give isl_mat *isl_mat_insert_zero_cols(__isl_take isl_mat *mat,
1181 unsigned first, unsigned n)
1187 mat = isl_mat_insert_cols(mat, first, n);
1191 for (i = 0; i < mat->n_row; ++i)
1192 isl_seq_clr(mat->row[i] + first, n);
1197 __isl_give isl_mat *isl_mat_add_zero_cols(__isl_take isl_mat *mat, unsigned n)
1202 return isl_mat_insert_zero_cols(mat, mat->n_col, n);
1205 __isl_give isl_mat *isl_mat_insert_rows(__isl_take isl_mat *mat,
1206 unsigned row, unsigned n)
1215 ext = isl_mat_alloc(mat->ctx, mat->n_row + n, mat->n_col);
1219 isl_mat_sub_copy(mat->ctx, ext->row, mat->row, row, 0, 0, mat->n_col);
1220 isl_mat_sub_copy(mat->ctx, ext->row + row + n, mat->row + row,
1221 mat->n_row - row, 0, 0, mat->n_col);
1230 __isl_give isl_mat *isl_mat_add_rows(__isl_take isl_mat *mat, unsigned n)
1235 return isl_mat_insert_rows(mat, mat->n_row, n);
1238 void isl_mat_col_submul(struct isl_mat *mat,
1239 int dst_col, isl_int f, int src_col)
1243 for (i = 0; i < mat->n_row; ++i)
1244 isl_int_submul(mat->row[i][dst_col], f, mat->row[i][src_col]);
1247 void isl_mat_col_mul(struct isl_mat *mat, int dst_col, isl_int f, int src_col)
1251 for (i = 0; i < mat->n_row; ++i)
1252 isl_int_mul(mat->row[i][dst_col], f, mat->row[i][src_col]);
1255 struct isl_mat *isl_mat_unimodular_complete(struct isl_mat *M, int row)
1258 struct isl_mat *H = NULL, *Q = NULL;
1263 isl_assert(M->ctx, M->n_row == M->n_col, goto error);
1265 H = isl_mat_left_hermite(isl_mat_copy(M), 0, NULL, &Q);
1266 M->n_row = M->n_col;
1269 for (r = 0; r < row; ++r)
1270 isl_assert(M->ctx, isl_int_is_one(H->row[r][r]), goto error);
1271 for (r = row; r < M->n_row; ++r)
1272 isl_seq_cpy(M->row[r], Q->row[r], M->n_col);
1283 __isl_give isl_mat *isl_mat_concat(__isl_take isl_mat *top,
1284 __isl_take isl_mat *bot)
1286 struct isl_mat *mat;
1291 isl_assert(top->ctx, top->n_col == bot->n_col, goto error);
1292 if (top->n_row == 0) {
1296 if (bot->n_row == 0) {
1301 mat = isl_mat_alloc(top->ctx, top->n_row + bot->n_row, top->n_col);
1304 isl_mat_sub_copy(mat->ctx, mat->row, top->row, top->n_row,
1306 isl_mat_sub_copy(mat->ctx, mat->row + top->n_row, bot->row, bot->n_row,
1317 int isl_mat_is_equal(__isl_keep isl_mat *mat1, __isl_keep isl_mat *mat2)
1324 if (mat1->n_row != mat2->n_row)
1327 if (mat1->n_col != mat2->n_col)
1330 for (i = 0; i < mat1->n_row; ++i)
1331 if (!isl_seq_eq(mat1->row[i], mat2->row[i], mat1->n_col))
1337 __isl_give isl_mat *isl_mat_from_row_vec(__isl_take isl_vec *vec)
1339 struct isl_mat *mat;
1343 mat = isl_mat_alloc(vec->ctx, 1, vec->size);
1347 isl_seq_cpy(mat->row[0], vec->el, vec->size);
1356 __isl_give isl_mat *isl_mat_vec_concat(__isl_take isl_mat *top,
1357 __isl_take isl_vec *bot)
1359 return isl_mat_concat(top, isl_mat_from_row_vec(bot));
1362 __isl_give isl_mat *isl_mat_move_cols(__isl_take isl_mat *mat,
1363 unsigned dst_col, unsigned src_col, unsigned n)
1369 if (n == 0 || dst_col == src_col)
1372 res = isl_mat_alloc(mat->ctx, mat->n_row, mat->n_col);
1376 if (dst_col < src_col) {
1377 isl_mat_sub_copy(res->ctx, res->row, mat->row, mat->n_row,
1379 isl_mat_sub_copy(res->ctx, res->row, mat->row, mat->n_row,
1380 dst_col, src_col, n);
1381 isl_mat_sub_copy(res->ctx, res->row, mat->row, mat->n_row,
1382 dst_col + n, dst_col, src_col - dst_col);
1383 isl_mat_sub_copy(res->ctx, res->row, mat->row, mat->n_row,
1384 src_col + n, src_col + n,
1385 res->n_col - src_col - n);
1387 isl_mat_sub_copy(res->ctx, res->row, mat->row, mat->n_row,
1389 isl_mat_sub_copy(res->ctx, res->row, mat->row, mat->n_row,
1390 src_col, src_col + n, dst_col - src_col);
1391 isl_mat_sub_copy(res->ctx, res->row, mat->row, mat->n_row,
1392 dst_col, src_col, n);
1393 isl_mat_sub_copy(res->ctx, res->row, mat->row, mat->n_row,
1394 dst_col + n, dst_col + n,
1395 res->n_col - dst_col - n);
1405 void isl_mat_gcd(__isl_keep isl_mat *mat, isl_int *gcd)
1410 isl_int_set_si(*gcd, 0);
1415 for (i = 0; i < mat->n_row; ++i) {
1416 isl_seq_gcd(mat->row[i], mat->n_col, &g);
1417 isl_int_gcd(*gcd, *gcd, g);
1422 __isl_give isl_mat *isl_mat_scale_down(__isl_take isl_mat *mat, isl_int m)
1429 for (i = 0; i < mat->n_row; ++i)
1430 isl_seq_scale_down(mat->row[i], mat->row[i], m, mat->n_col);
1435 __isl_give isl_mat *isl_mat_normalize(__isl_take isl_mat *mat)
1443 isl_mat_gcd(mat, &gcd);
1444 mat = isl_mat_scale_down(mat, gcd);