// SymPy's cholesky decomposition void cholesky(const DenseMatrix &A, DenseMatrix &L) { SYMENGINE_ASSERT(A.row_ == A.col_); SYMENGINE_ASSERT(L.row_ == A.row_ and L.col_ == A.row_); unsigned col = A.col_; unsigned i, j, k; RCP<const Basic> sum; RCP<const Basic> i2 = integer(2); RCP<const Basic> half = div(one, i2); // Initialize L for (i = 0; i < col; i++) for (j = 0; j < col; j++) L.m_[i * col + j] = zero; for (i = 0; i < col; i++) { for (j = 0; j < i; j++) { sum = zero; for (k = 0; k < j; k++) sum = add(sum, mul(L.m_[i * col + k], L.m_[j * col + k])); L.m_[i * col + j] = mul(div(one, L.m_[j * col + j]), sub(A.m_[i * col + j], sum)); } sum = zero; for (k = 0; k < i; k++) sum = add(sum, pow(L.m_[i * col + k], i2)); L.m_[i * col + i] = pow(sub(A.m_[i * col + i], sum), half); } }
void jacobian(const DenseMatrix &A, const DenseMatrix &x, DenseMatrix &result) { SYMENGINE_ASSERT(A.col_ == 1); SYMENGINE_ASSERT(x.col_ == 1); SYMENGINE_ASSERT(A.row_ == result.nrows() and x.row_ == result.ncols()); bool error = false; #pragma omp parallel for for (unsigned i = 0; i < result.row_; i++) { for (unsigned j = 0; j < result.col_; j++) { if (is_a<Symbol>(*(x.m_[j]))) { const RCP<const Symbol> x_ = rcp_static_cast<const Symbol>(x.m_[j]); result.m_[i * result.col_ + j] = A.m_[i]->diff(x_); } else { error = true; break; } } } if (error) { throw SymEngineException( "'x' must contain Symbols only. " "Use sjacobian for SymPy style differentiation"); } }
void char_poly(const DenseMatrix &A, DenseMatrix &B) { SYMENGINE_ASSERT(B.ncols() == 1 and B.nrows() == A.nrows() + 1); SYMENGINE_ASSERT(A.nrows() == A.ncols()); std::vector<DenseMatrix> polys; berkowitz(A, polys); B = polys[polys.size() - 1]; }
// --------------------------- Solve Ax = b ---------------------------------// // Assuming A is a diagonal square matrix void diagonal_solve(const DenseMatrix &A, const DenseMatrix &b, DenseMatrix &x) { SYMENGINE_ASSERT(A.row_ == A.col_); SYMENGINE_ASSERT(x.row_ == A.col_ and x.col_ == b.col_); const unsigned sys = b.col_; // No checks are done to see if the diagonal entries are zero for (unsigned k = 0; k < sys; k++) { for (unsigned i = 0; i < A.col_; i++) { x.m_[i * sys + k] = div(b.m_[i * sys + k], A.m_[i * A.col_ + i]); } } }
// ------------------------------- Submatrix ---------------------------------// void submatrix_dense(const DenseMatrix &A, DenseMatrix &B, unsigned row_start, unsigned col_start, unsigned row_end, unsigned col_end, unsigned row_step, unsigned col_step) { SYMENGINE_ASSERT(row_end >= row_start and col_end >= col_start); SYMENGINE_ASSERT(row_end < A.row_); SYMENGINE_ASSERT(col_end < A.col_); SYMENGINE_ASSERT(B.row_ == row_end - row_start + 1 and B.col_ == col_end - col_start + 1); unsigned row = B.row_, col = B.col_; for (unsigned i = 0; i < row; i += row_step) for (unsigned j = 0; j < col; j += col_step) B.m_[i * col + j] = A.m_[(row_start + i) * A.col_ + col_start + j]; }
// Extract main diagonal of CSR matrix A void csr_diagonal(const CSRMatrix &A, DenseMatrix &D) { unsigned N = std::min(A.row_, A.col_); SYMENGINE_ASSERT(D.nrows() == N and D.ncols() == 1); unsigned row_start; unsigned row_end; RCP<const Basic> diag; for (unsigned i = 0; i < N; i++) { row_start = A.p_[i]; row_end = A.p_[i + 1]; diag = zero; unsigned jj; while (row_start <= row_end) { jj = (row_start + row_end) / 2; if (A.j_[jj] == i) { diag = A.x_[jj]; break; } else if (A.j_[jj] < i) { row_start = jj + 1; } else { row_end = jj - 1; } } D.set(i, 0, diag); } }
// Create diagonal matrices directly void diag(DenseMatrix &A, vec_basic &v, int k) { SYMENGINE_ASSERT(v.size() > 0); unsigned k_ = std::abs(k); if (k >= 0) { for (unsigned i = 0; i < A.row_; i++) { for (unsigned j = 0; j < A.col_; j++) { if (j != (unsigned)k) { A.m_[i * A.col_ + j] = zero; } else { A.m_[i * A.col_ + j] = v[k - k_]; } } k++; } } else { k = -k; for (unsigned j = 0; j < A.col_; j++) { for (unsigned i = 0; i < A.row_; i++) { if (i != (unsigned)k) { A.m_[i * A.col_ + j] = zero; } else { A.m_[i * A.col_ + j] = v[k - k_]; } } k++; } } }
// Get and set elements RCP<const Basic> CSRMatrix::get(unsigned i, unsigned j) const { SYMENGINE_ASSERT(i < row_ and j < col_); unsigned row_start = p_[i]; unsigned row_end = p_[i + 1]; unsigned k; if (row_start == row_end) { return zero; } while (row_start <= row_end) { k = (row_start + row_end) / 2; if (j_[k] == j) { return x_[k]; } else if (j_[k] < j) { row_start = k + 1; } else { row_end = k - 1; } } return zero; }
// Algorithm 2, page 13, Nakos, G. C., Turner, P. R., Williams, R. M. (1997). // Fraction-free algorithms for linear and polynomial equations. // ACM SIGSAM Bulletin, 31(3), 11–19. doi:10.1145/271130.271133. void fraction_free_gauss_jordan_elimination(const DenseMatrix &A, DenseMatrix &B) { SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_); unsigned row = A.row_, col = A.col_; unsigned i, j, k; RCP<const Basic> d; B.m_ = A.m_; for (i = 0; i < col; i++) { if (i > 0) d = B.m_[i * col - col + i - 1]; for (j = 0; j < row; j++) if (j != i) for (k = 0; k < col; k++) { if (k != i) { B.m_[j * col + k] = sub(mul(B.m_[i * col + i], B.m_[j * col + k]), mul(B.m_[j * col + i], B.m_[i * col + k])); if (i > 0) B.m_[j * col + k] = div(B.m_[j * col + k], d); } } for (j = 0; j < row; j++) if (j != i) B.m_[j * col + i] = zero; } }
void Mul::as_base_exp(const RCP<const Basic> &self, const Ptr<RCP<const Basic>> &exp, const Ptr<RCP<const Basic>> &base) { if (is_a_Number(*self)) { // Always ensure it is of form |num| > |den| // in case of Integers den = 1 if (is_a<Rational>(*self)) { RCP<const Rational> self_new = rcp_static_cast<const Rational>(self); if (abs(self_new->i.get_num()) < abs(self_new->i.get_den())) { *exp = minus_one; *base = self_new->rdiv(*rcp_static_cast<const Number>(one)); } else { *exp = one; *base = self; } } else { *exp = one; *base = self; } } else if (is_a<Pow>(*self)) { *exp = rcp_static_cast<const Pow>(self)->exp_; *base = rcp_static_cast<const Pow>(self)->base_; } else { SYMENGINE_ASSERT(!is_a<Mul>(*self)); *exp = one; *base = self; } }
void fraction_free_gaussian_elimination_solve(const DenseMatrix &A, const DenseMatrix &b, DenseMatrix &x) { SYMENGINE_ASSERT(A.row_ == A.col_); SYMENGINE_ASSERT(b.row_ == A.row_ and x.row_ == A.row_); SYMENGINE_ASSERT(x.col_ == b.col_); int i, j, k, col = A.col_, bcol = b.col_; DenseMatrix A_ = DenseMatrix(A.row_, A.col_, A.m_); DenseMatrix b_ = DenseMatrix(b.row_, b.col_, b.m_); for (i = 0; i < col - 1; i++) for (j = i + 1; j < col; j++) { for (k = 0; k < bcol; k++) { b_.m_[j * bcol + k] = sub(mul(A_.m_[i * col + i], b_.m_[j * bcol + k]), mul(A_.m_[j * col + i], b_.m_[i * bcol + k])); if (i > 0) b_.m_[j * bcol + k] = div(b_.m_[j * bcol + k], A_.m_[i * col - col + i - 1]); } for (k = i + 1; k < col; k++) { A_.m_[j * col + k] = sub(mul(A_.m_[i * col + i], A_.m_[j * col + k]), mul(A_.m_[j * col + i], A_.m_[i * col + k])); if (i > 0) A_.m_[j * col + k] = div(A_.m_[j * col + k], A_.m_[i * col - col + i - 1]); } A_.m_[j * col + i] = zero; } for (i = 0; i < col * bcol; i++) x.m_[i] = zero; // Integer zero; for (k = 0; k < bcol; k++) { for (i = col - 1; i >= 0; i--) { for (j = i + 1; j < col; j++) b_.m_[i * bcol + k] = sub(b_.m_[i * bcol + k], mul(A_.m_[i * col + j], x.m_[j * bcol + k])); x.m_[i * bcol + k] = div(b_.m_[i * bcol + k], A_.m_[i * col + i]); } } }
void fraction_free_gauss_jordan_solve(const DenseMatrix &A, const DenseMatrix &b, DenseMatrix &x) { SYMENGINE_ASSERT(A.row_ == A.col_); SYMENGINE_ASSERT(b.row_ == A.row_ and x.row_ == A.row_); SYMENGINE_ASSERT(x.col_ == b.col_); unsigned i, j, k, col = A.col_, bcol = b.col_; RCP<const Basic> d; DenseMatrix A_ = DenseMatrix(A.row_, A.col_, A.m_); DenseMatrix b_ = DenseMatrix(b.row_, b.col_, b.m_); for (i = 0; i < col; i++) { if (i > 0) d = A_.m_[i * col - col + i - 1]; for (j = 0; j < col; j++) if (j != i) { for (k = 0; k < bcol; k++) { b_.m_[j * bcol + k] = sub(mul(A_.m_[i * col + i], b_.m_[j * bcol + k]), mul(A_.m_[j * col + i], b_.m_[i * bcol + k])); if (i > 0) b_.m_[j * bcol + k] = div(b_.m_[j * bcol + k], d); } for (k = 0; k < col; k++) { if (k != i) { A_.m_[j * col + k] = sub(mul(A_.m_[i * col + i], A_.m_[j * col + k]), mul(A_.m_[j * col + i], A_.m_[i * col + k])); if (i > 0) A_.m_[j * col + k] = div(A_.m_[j * col + k], d); } } } for (j = 0; j < col; j++) if (j != i) A_.m_[j * col + i] = zero; } // No checks are done to see if the diagonal entries are zero for (k = 0; k < bcol; k++) for (i = 0; i < col; i++) x.m_[i * bcol + k] = div(b_.m_[i * bcol + k], A_.m_[i * col + i]); }
// ----------------------------- Matrix Transpose ----------------------------// void transpose_dense(const DenseMatrix &A, DenseMatrix &B) { SYMENGINE_ASSERT(B.row_ == A.col_ and B.col_ == A.row_); for (unsigned i = 0; i < A.row_; i++) for (unsigned j = 0; j < A.col_; j++) B.m_[j * B.col_ + i] = A.m_[i * A.col_ + j]; }
void row_mul_scalar_dense(DenseMatrix &A, unsigned i, RCP<const Basic> &c) { SYMENGINE_ASSERT(i < A.row_); unsigned col = A.col_; for (unsigned j = 0; j < A.col_; j++) A.m_[i * col + j] = mul(c, A.m_[i * col + j]); }
// -------------------------------- Row Operations ---------------------------// void row_exchange_dense(DenseMatrix &A, unsigned i, unsigned j) { SYMENGINE_ASSERT(i != j and i < A.row_ and j < A.row_); unsigned col = A.col_; for (unsigned k = 0; k < A.col_; k++) std::swap(A.m_[i * col + k], A.m_[j * col + k]); }
void diff(const DenseMatrix &A, const RCP<const Symbol> &x, DenseMatrix &result) { SYMENGINE_ASSERT(A.row_ == result.nrows() and A.col_ == result.ncols()); #pragma omp parallel for for (unsigned i = 0; i < result.row_; i++) { for (unsigned j = 0; j < result.col_; j++) { result.m_[i * result.col_ + j] = A.m_[i * result.col_ + j]->diff(x); } } }
void row_add_row_dense(DenseMatrix &A, unsigned i, unsigned j, RCP<const Basic> &c) { SYMENGINE_ASSERT(i != j and i < A.row_ and j < A.row_); unsigned col = A.col_; for (unsigned k = 0; k < A.col_; k++) A.m_[i * col + k] = add(A.m_[i * col + k], mul(c, A.m_[j * col + k])); }
// SymPy LUDecomposition algorithm, in // sympy.matrices.matrices.Matrix.LUdecomposition_Simple with pivoting. // P must be an initialized matrix and will be permuted. void pivoted_LU(const DenseMatrix &A, DenseMatrix &L, DenseMatrix &U, permutelist &pl) { SYMENGINE_ASSERT(A.row_ == A.col_ and L.row_ == L.col_ and U.row_ == U.col_); SYMENGINE_ASSERT(A.row_ == L.row_ and A.row_ == U.row_); pivoted_LU(A, U, pl); unsigned n = A.col_; for (unsigned i = 0; i < n; i++) { for (unsigned j = 0; j < i; j++) { L.m_[i * n + j] = U.m_[i * n + j]; U.m_[i * n + j] = zero; // Integer zero } L.m_[i * n + i] = one; // Integer one for (unsigned j = i + 1; j < n; j++) L.m_[i * n + j] = zero; // Integer zero } }
// SymPy's fraction free LU decomposition, without pivoting // sympy.matrices.matrices.MatrixBase.LUdecompositionFF // W. Zhou & D.J. Jeffrey, "Fraction-free matrix factors: new forms for LU and // QR factors". // Frontiers in Computer Science in China, Vol 2, no. 1, pp. 67-80, 2008. void fraction_free_LDU(const DenseMatrix &A, DenseMatrix &L, DenseMatrix &D, DenseMatrix &U) { SYMENGINE_ASSERT(A.row_ == L.row_ and A.row_ == U.row_); SYMENGINE_ASSERT(A.col_ == L.col_ and A.col_ == U.col_); unsigned row = A.row_, col = A.col_; unsigned i, j, k; RCP<const Basic> old = integer(1); U.m_ = A.m_; // Initialize L for (i = 0; i < row; i++) for (j = 0; j < row; j++) if (i != j) L.m_[i * col + j] = zero; else L.m_[i * col + i] = one; // Initialize D for (i = 0; i < row * row; i++) D.m_[i] = zero; // Integer zero for (k = 0; k < row - 1; k++) { L.m_[k * col + k] = U.m_[k * col + k]; D.m_[k * col + k] = mul(old, U.m_[k * col + k]); for (i = k + 1; i < row; i++) { L.m_[i * col + k] = U.m_[i * col + k]; for (j = k + 1; j < col; j++) U.m_[i * col + j] = div(sub(mul(U.m_[k * col + k], U.m_[i * col + j]), mul(U.m_[k * col + j], U.m_[i * col + k])), old); U.m_[i * col + k] = zero; // Integer zero } old = U.m_[k * col + k]; } D.m_[row * col - col + row - 1] = old; }
void DenseMatrix::add_matrix(const MatrixBase &other, MatrixBase &result) const { SYMENGINE_ASSERT(row_ == result.nrows() and col_ == result.ncols()); if (is_a<DenseMatrix>(other) and is_a<DenseMatrix>(result)) { const DenseMatrix &o = static_cast<const DenseMatrix &>(other); DenseMatrix &r = static_cast<DenseMatrix &>(result); add_dense_dense(*this, o, r); } }
// Scale the rows of a CSR matrix *in place* // A[i, :] *= X[i] void csr_scale_rows(CSRMatrix &A, const DenseMatrix &X) { SYMENGINE_ASSERT(A.row_ == X.nrows() and X.ncols() == 1); for (unsigned i = 0; i < A.row_; i++) { if (eq(*(X.get(i, 0)), *zero)) throw SymEngineException("Scaling factor can't be zero"); for (unsigned jj = A.p_[i]; jj < A.p_[i + 1]; jj++) A.x_[jj] = mul(A.x_[jj], X.get(i, 0)); } }
void sjacobian(const DenseMatrix &A, const DenseMatrix &x, DenseMatrix &result) { SYMENGINE_ASSERT(A.col_ == 1); SYMENGINE_ASSERT(x.col_ == 1); SYMENGINE_ASSERT(A.row_ == result.nrows() and x.row_ == result.ncols()); #pragma omp parallel for for (unsigned i = 0; i < result.row_; i++) { for (unsigned j = 0; j < result.col_; j++) { if (is_a<Symbol>(*(x.m_[j]))) { const RCP<const Symbol> x_ = rcp_static_cast<const Symbol>(x.m_[j]); result.m_[i * result.col_ + j] = A.m_[i]->diff(x_); } else { // TODO: Use a dummy symbol const RCP<const Symbol> x_ = symbol("x_"); result.m_[i * result.col_ + j] = ssubs( ssubs(A.m_[i], {{x.m_[j], x_}})->diff(x_), {{x_, x.m_[j]}}); } } } }
void back_substitution(const DenseMatrix &U, const DenseMatrix &b, DenseMatrix &x) { SYMENGINE_ASSERT(U.row_ == U.col_); SYMENGINE_ASSERT(b.row_ == U.row_); SYMENGINE_ASSERT(x.row_ == U.col_ and x.col_ == b.col_); const unsigned col = U.col_; const unsigned sys = b.col_; x.m_ = b.m_; for (unsigned k = 0; k < sys; k++) { for (int i = col - 1; i >= 0; i--) { for (unsigned j = i + 1; j < col; j++) x.m_[i * sys + k] = sub(x.m_[i * sys + k], mul(U.m_[i * col + j], x.m_[j * sys + k])); x.m_[i * sys + k] = div(x.m_[i * sys + k], U.m_[i * col + i]); } } }
// SymPy LUDecomposition algorithm, in // sympy.matrices.matrices.Matrix.LUdecomposition // with no pivoting void LU(const DenseMatrix &A, DenseMatrix &L, DenseMatrix &U) { SYMENGINE_ASSERT(A.row_ == A.col_ and L.row_ == L.col_ and U.row_ == U.col_); SYMENGINE_ASSERT(A.row_ == L.row_ and A.row_ == U.row_); unsigned n = A.row_; unsigned i, j, k; RCP<const Basic> scale; U.m_ = A.m_; for (j = 0; j < n; j++) { for (i = 0; i < j; i++) for (k = 0; k < i; k++) U.m_[i * n + j] = sub(U.m_[i * n + j], mul(U.m_[i * n + k], U.m_[k * n + j])); for (i = j; i < n; i++) { for (k = 0; k < j; k++) U.m_[i * n + j] = sub(U.m_[i * n + j], mul(U.m_[i * n + k], U.m_[k * n + j])); } scale = div(one, U.m_[j * n + j]); for (i = j + 1; i < n; i++) U.m_[i * n + j] = mul(U.m_[i * n + j], scale); } for (i = 0; i < n; i++) { for (j = 0; j < i; j++) { L.m_[i * n + j] = U.m_[i * n + j]; U.m_[i * n + j] = zero; // Integer zero } L.m_[i * n + i] = one; // Integer one for (j = i + 1; j < n; j++) L.m_[i * n + j] = zero; // Integer zero } }
// SymPy's QRecomposition in sympy.matrices.matrices.MatrixBase.QRdecomposition // Rank check is not performed void QR(const DenseMatrix &A, DenseMatrix &Q, DenseMatrix &R) { unsigned row = A.row_; unsigned col = A.col_; SYMENGINE_ASSERT(Q.row_ == row and Q.col_ == col and R.row_ == col and R.col_ == col); unsigned i, j, k; RCP<const Basic> t; std::vector<RCP<const Basic>> tmp(row); // Initialize Q for (i = 0; i < row * col; i++) Q.m_[i] = zero; // Initialize R for (i = 0; i < col * col; i++) R.m_[i] = zero; for (j = 0; j < col; j++) { // Use submatrix for this for (k = 0; k < row; k++) tmp[k] = A.m_[k * col + j]; for (i = 0; i < j; i++) { t = zero; for (k = 0; k < row; k++) t = add(t, mul(A.m_[k * col + j], Q.m_[k * col + i])); for (k = 0; k < row; k++) tmp[k] = expand(sub(tmp[k], mul(Q.m_[k * col + i], t))); } // calculate norm t = zero; for (k = 0; k < row; k++) t = add(t, pow(tmp[k], integer(2))); t = pow(t, div(one, integer(2))); R.m_[j * col + j] = t; for (k = 0; k < row; k++) Q.m_[k * col + j] = div(tmp[k], t); for (i = 0; i < j; i++) { t = zero; for (k = 0; k < row; k++) t = add(t, mul(Q.m_[k * col + i], A.m_[k * col + j])); R.m_[i * col + j] = t; } } }
// ------------------------------- Matrix Addition ---------------------------// void add_dense_dense(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &C) { SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_ and A.row_ == C.row_ and A.col_ == C.col_); unsigned row = A.row_, col = A.col_; for (unsigned i = 0; i < row; i++) { for (unsigned j = 0; j < col; j++) { C.m_[i * col + j] = add(A.m_[i * col + j], B.m_[i * col + j]); } } }
void mul_dense_scalar(const DenseMatrix &A, const RCP<const Basic> &k, DenseMatrix &B) { SYMENGINE_ASSERT(A.col_ == B.col_ and A.row_ == B.row_); unsigned row = A.row_, col = A.col_; for (unsigned i = 0; i < row; i++) { for (unsigned j = 0; j < col; j++) { B.m_[i * col + j] = mul(A.m_[i * col + j], k); } } }
// ----------------------------- Determinant ---------------------------------// RCP<const Basic> det_bareis(const DenseMatrix &A) { SYMENGINE_ASSERT(A.row_ == A.col_); unsigned n = A.row_; if (n == 1) { return A.m_[0]; } else if (n == 2) { // If A = [[a, b], [c, d]] then det(A) = ad - bc return sub(mul(A.m_[0], A.m_[3]), mul(A.m_[1], A.m_[2])); } else if (n == 3) { // if A = [[a, b, c], [d, e, f], [g, h, i]] then // det(A) = (aei + bfg + cdh) - (ceg + bdi + afh) return sub(add(add(mul(mul(A.m_[0], A.m_[4]), A.m_[8]), mul(mul(A.m_[1], A.m_[5]), A.m_[6])), mul(mul(A.m_[2], A.m_[3]), A.m_[7])), add(add(mul(mul(A.m_[2], A.m_[4]), A.m_[6]), mul(mul(A.m_[1], A.m_[3]), A.m_[8])), mul(mul(A.m_[0], A.m_[5]), A.m_[7]))); } else { DenseMatrix B = DenseMatrix(n, n, A.m_); unsigned i, sign = 1; RCP<const Basic> d; for (unsigned k = 0; k < n - 1; k++) { if (eq(*(B.m_[k * n + k]), *zero)) { for (i = k + 1; i < n; i++) if (neq(*(B.m_[i * n + k]), *zero)) { row_exchange_dense(B, i, k); sign *= -1; break; } if (i == n) return zero; } for (i = k + 1; i < n; i++) { for (unsigned j = k + 1; j < n; j++) { d = sub(mul(B.m_[k * n + k], B.m_[i * n + j]), mul(B.m_[i * n + k], B.m_[k * n + j])); if (k > 0) d = div(d, B.m_[(k - 1) * n + k - 1]); B.m_[i * n + j] = d; } } } return (sign == 1) ? B.m_[n * n - 1] : mul(minus_one, B.m_[n * n - 1]); } }
// SymPy LUDecomposition_Simple algorithm, in // sympy.matrices.matrices.Matrix.LUdecomposition_Simple with pivoting. // P must be an initialized matrix and will be permuted. void pivoted_LU(const DenseMatrix &A, DenseMatrix &LU, permutelist &pl) { SYMENGINE_ASSERT(A.row_ == A.col_ and LU.row_ == LU.col_); SYMENGINE_ASSERT(A.row_ == LU.row_); unsigned n = A.row_; unsigned i, j, k; RCP<const Basic> scale; int pivot; LU.m_ = A.m_; for (j = 0; j < n; j++) { for (i = 0; i < j; i++) for (k = 0; k < i; k++) LU.m_[i * n + j] = sub(LU.m_[i * n + j], mul(LU.m_[i * n + k], LU.m_[k * n + j])); pivot = -1; for (i = j; i < n; i++) { for (k = 0; k < j; k++) { LU.m_[i * n + j] = sub(LU.m_[i * n + j], mul(LU.m_[i * n + k], LU.m_[k * n + j])); } if (pivot == -1 and neq(*LU.m_[i * n + j], *zero)) pivot = i; } if (pivot == -1) throw SymEngineException("Matrix is rank deficient"); if (pivot - j != 0) { // row must be swapped row_exchange_dense(LU, pivot, j); pl.push_back({pivot, j}); } scale = div(one, LU.m_[j * n + j]); for (i = j + 1; i < n; i++) LU.m_[i * n + j] = mul(LU.m_[i * n + j], scale); } }
bool order(const DenseMatrix &t, const std::vector<DenseMatrix> &basis, unsigned k) { bool eq = true; for (unsigned j = 0; j < t.ncols(); j++) { SYMENGINE_ASSERT(is_a<Integer>(*t.get(0, j))); integer_class t_ = rcp_static_cast<const Integer>(t.get(0, j))->as_integer_class(); SYMENGINE_ASSERT(is_a<Integer>(*basis[k].get(0, j))); integer_class b_ = rcp_static_cast<const Integer>(basis[k].get(0, j)) ->as_integer_class(); if (t_ < b_) { return false; } if (t_ > b_) { eq = false; } } return not eq; }