1 #include <symengine/matrix.h>
6 #include <symengine/subs.h>
7 #include <symengine/symengine_exception.h>
8 #include <symengine/polys/uexprpoly.h>
10 #include <symengine/test_visitors.h>
11 #if HAVE_SYMENGINE_RTTI
12 #include <symengine/serialize-cereal.h>
19 DenseMatrix::DenseMatrix() {}
21 DenseMatrix::DenseMatrix(
unsigned row,
unsigned col) : row_(row), col_(col)
26 DenseMatrix::DenseMatrix(
unsigned row,
unsigned col,
const vec_basic &l)
27 : m_{l}, row_(row), col_(col){SYMENGINE_ASSERT(m_.size() == row * col)}
29 DenseMatrix::DenseMatrix(
const vec_basic &column_elements)
30 : m_(column_elements), row_(static_cast<unsigned>(column_elements.size())),
36 void DenseMatrix::resize(
unsigned row,
unsigned col)
44 RCP<const Basic> DenseMatrix::get(
unsigned i,
unsigned j)
const
46 SYMENGINE_ASSERT(i < row_ and j < col_);
47 return m_[i * col_ + j];
50 void DenseMatrix::set(
unsigned i,
unsigned j,
const RCP<const Basic> &e)
52 SYMENGINE_ASSERT(i < row_ and j < col_);
57 vec_basic DenseMatrix::as_vec_basic()
const
62 RCP<const Basic> DenseMatrix::trace()
const
64 SYMENGINE_ASSERT(row_ == col_);
67 for (
unsigned i = 0; i < row_; i++) {
68 diag.push_back(m_[offset]);
75 unsigned DenseMatrix::rank()
const
77 throw NotImplementedError(
"Not Implemented");
80 bool DenseMatrix::is_lower()
const
83 unsigned n = A.nrows();
84 for (
unsigned i = 1; i < n; ++i) {
85 for (
unsigned j = 0; j < i; ++j) {
94 bool DenseMatrix::is_upper()
const
97 unsigned n = A.nrows();
98 for (
unsigned i = 0; i < n - 1; ++i) {
99 for (
unsigned j = i + 1; j < n; ++j) {
111 tribool cur = tribool::tritrue;
113 cur = and_tribool(cur, SymEngine::is_zero(*e));
121 tribool DenseMatrix::is_diagonal()
const
124 if (not A.is_square()) {
125 return tribool::trifalse;
127 unsigned ncols = A.ncols();
129 tribool cur = tribool::tritrue;
130 for (
unsigned i = 0; i < ncols; i++) {
132 for (
unsigned j = 0; j < ncols; j++) {
134 auto &e = m_[offset];
135 cur = and_tribool(cur, SymEngine::is_zero(*e));
146 tribool DenseMatrix::is_real(
const Assumptions *assumptions)
const
148 RealVisitor visitor(assumptions);
149 tribool cur = tribool::tritrue;
151 cur = and_tribool(cur, visitor.apply(*e));
159 tribool DenseMatrix::is_symmetric()
const
162 if (not A.is_square()) {
163 return tribool::trifalse;
165 unsigned ncols = A.ncols();
166 tribool cur = tribool::tritrue;
167 for (
unsigned i = 0; i < ncols; i++) {
168 for (
unsigned j = 0; j <= i; j++) {
170 auto &e = m_[i * ncols + j];
171 auto &e2 = m_[j * ncols + i];
172 cur = and_tribool(cur, SymEngine::is_zero(*
sub(e, e2)));
182 tribool DenseMatrix::is_hermitian()
const
185 if (not A.is_square()) {
186 return tribool::trifalse;
188 unsigned ncols = A.ncols();
189 tribool cur = tribool::tritrue;
190 for (
unsigned i = 0; i < ncols; i++) {
191 for (
unsigned j = 0; j <= i; j++) {
192 auto &e = m_[i * ncols + j];
194 auto &e2 = m_[j * ncols + i];
198 cur = and_tribool(cur, SymEngine::is_real(*e));
208 tribool DenseMatrix::is_weakly_diagonally_dominant()
const
211 if (not A.is_square()) {
212 return tribool::trifalse;
215 unsigned ncols = A.ncols();
216 RCP<const Basic> diag;
217 RCP<const Basic> sum;
218 tribool diagdom = tribool::tritrue;
219 for (
unsigned i = 0; i < ncols; i++) {
221 for (
unsigned j = 0; j < ncols; j++) {
222 auto &e = m_[i * ncols + j];
229 diagdom = and_tribool(diagdom, is_nonnegative(*
sub(diag, sum)));
230 if (is_false(diagdom)) {
237 tribool DenseMatrix::is_strictly_diagonally_dominant()
const
240 if (not A.is_square()) {
241 return tribool::trifalse;
244 unsigned ncols = A.ncols();
245 RCP<const Basic> diag;
246 RCP<const Basic> sum;
247 tribool diagdom = tribool::tritrue;
248 for (
unsigned i = 0; i < ncols; i++) {
250 for (
unsigned j = 0; j < ncols; j++) {
251 auto &e = m_[i * ncols + j];
258 diagdom = and_tribool(diagdom, is_positive(*
sub(diag, sum)));
259 if (is_false(diagdom)) {
266 tribool DenseMatrix::shortcut_to_posdef()
const
268 tribool is_diagonal_positive = tribool::tritrue;
270 for (
unsigned i = 0; i < row_; i++) {
272 = and_tribool(is_diagonal_positive, is_positive(*m_[offset]));
273 if (is_false(is_diagonal_positive))
274 return is_diagonal_positive;
277 if (is_true(and_tribool(is_diagonal_positive,
278 this->is_strictly_diagonally_dominant())))
279 return tribool::tritrue;
280 return tribool::indeterminate;
283 tribool DenseMatrix::is_positive_definite_GE()
286 for (
unsigned i = 0; i < size; i++) {
287 auto ispos = is_positive(*m_[i * size + i]);
290 for (
unsigned j = i + 1; j < size; j++) {
291 for (
unsigned k = i + 1; k < size; k++) {
292 m_[j * size + k] =
sub(
mul(m_[i * size + i], m_[j * size + k]),
293 mul(m_[j * size + i], m_[i * size + k]));
297 return tribool::tritrue;
300 tribool DenseMatrix::is_positive_definite()
const
305 const DenseMatrix *H;
306 if (!is_true(A.is_hermitian())) {
308 return tribool::trifalse;
309 DenseMatrix tmp1 = DenseMatrix(A.row_, A.col_);
311 A.conjugate_transpose(tmp1);
312 add_dense_dense(A, tmp1, *B.
get());
318 tribool shortcut = H->shortcut_to_posdef();
319 if (!is_indeterminate(shortcut))
325 return B->is_positive_definite_GE();
328 tribool DenseMatrix::is_negative_definite()
const
330 auto res = DenseMatrix(row_, col_);
331 mul_dense_scalar(*
this,
integer(-1), res);
332 return res.is_positive_definite();
335 RCP<const Basic> DenseMatrix::det()
const
337 return det_bareis(*
this);
340 void DenseMatrix::inv(MatrixBase &result)
const
342 if (is_a<DenseMatrix>(result)) {
343 DenseMatrix &r = down_cast<DenseMatrix &>(result);
344 inverse_pivoted_LU(*
this, r);
348 void DenseMatrix::add_matrix(
const MatrixBase &other, MatrixBase &result)
const
350 SYMENGINE_ASSERT(row_ == result.nrows() and col_ == result.ncols());
352 if (is_a<DenseMatrix>(other) and is_a<DenseMatrix>(result)) {
353 const DenseMatrix &o = down_cast<const DenseMatrix &>(other);
354 DenseMatrix &r = down_cast<DenseMatrix &>(result);
355 add_dense_dense(*
this, o, r);
359 void DenseMatrix::mul_matrix(
const MatrixBase &other, MatrixBase &result)
const
361 SYMENGINE_ASSERT(row_ == result.nrows()
362 and other.ncols() == result.ncols());
364 if (is_a<DenseMatrix>(other) and is_a<DenseMatrix>(result)) {
365 const DenseMatrix &o = down_cast<const DenseMatrix &>(other);
366 DenseMatrix &r = down_cast<DenseMatrix &>(result);
367 mul_dense_dense(*
this, o, r);
371 void DenseMatrix::elementwise_mul_matrix(
const MatrixBase &other,
372 MatrixBase &result)
const
374 SYMENGINE_ASSERT(row_ == result.nrows() and col_ == result.ncols()
375 and row_ == other.nrows() and col_ == other.ncols());
377 if (is_a<DenseMatrix>(other) and is_a<DenseMatrix>(result)) {
378 const DenseMatrix &o = down_cast<const DenseMatrix &>(other);
379 DenseMatrix &r = down_cast<DenseMatrix &>(result);
380 elementwise_mul_dense_dense(*
this, o, r);
385 void DenseMatrix::add_scalar(
const RCP<const Basic> &k,
386 MatrixBase &result)
const
388 if (is_a<DenseMatrix>(result)) {
389 DenseMatrix &r = down_cast<DenseMatrix &>(result);
390 add_dense_scalar(*
this, k, r);
395 void DenseMatrix::mul_scalar(
const RCP<const Basic> &k,
396 MatrixBase &result)
const
398 if (is_a<DenseMatrix>(result)) {
399 DenseMatrix &r = down_cast<DenseMatrix &>(result);
400 mul_dense_scalar(*
this, k, r);
405 void DenseMatrix::conjugate(MatrixBase &result)
const
407 if (is_a<DenseMatrix>(result)) {
408 DenseMatrix &r = down_cast<DenseMatrix &>(result);
409 conjugate_dense(*
this, r);
414 void DenseMatrix::transpose(MatrixBase &result)
const
416 if (is_a<DenseMatrix>(result)) {
417 DenseMatrix &r = down_cast<DenseMatrix &>(result);
418 transpose_dense(*
this, r);
423 void DenseMatrix::conjugate_transpose(MatrixBase &result)
const
425 if (is_a<DenseMatrix>(result)) {
426 DenseMatrix &r = down_cast<DenseMatrix &>(result);
427 conjugate_transpose_dense(*
this, r);
432 void DenseMatrix::submatrix(MatrixBase &result,
unsigned row_start,
433 unsigned col_start,
unsigned row_end,
434 unsigned col_end,
unsigned row_step,
435 unsigned col_step)
const
437 if (is_a<DenseMatrix>(result)) {
438 DenseMatrix &r = down_cast<DenseMatrix &>(result);
439 submatrix_dense(*
this, r, row_start, col_start, row_end, col_end,
445 void DenseMatrix::LU(MatrixBase &L, MatrixBase &U)
const
447 if (is_a<DenseMatrix>(L) and is_a<DenseMatrix>(U)) {
448 DenseMatrix &L_ = down_cast<DenseMatrix &>(L);
449 DenseMatrix &U_ = down_cast<DenseMatrix &>(U);
450 SymEngine::LU(*
this, L_, U_);
455 void DenseMatrix::LDL(MatrixBase &L, MatrixBase &D)
const
457 if (is_a<DenseMatrix>(L) and is_a<DenseMatrix>(D)) {
458 DenseMatrix &L_ = down_cast<DenseMatrix &>(L);
459 DenseMatrix &D_ = down_cast<DenseMatrix &>(D);
460 SymEngine::LDL(*
this, L_, D_);
465 void DenseMatrix::LU_solve(
const MatrixBase &b, MatrixBase &x)
const
467 if (is_a<DenseMatrix>(b) and is_a<DenseMatrix>(x)) {
468 const DenseMatrix &b_ = down_cast<const DenseMatrix &>(b);
469 DenseMatrix &x_ = down_cast<DenseMatrix &>(x);
470 SymEngine::LU_solve(*
this, b_, x_);
475 void DenseMatrix::FFLU(MatrixBase &LU)
const
477 if (is_a<DenseMatrix>(LU)) {
478 DenseMatrix &LU_ = down_cast<DenseMatrix &>(LU);
479 fraction_free_LU(*
this, LU_);
484 void DenseMatrix::FFLDU(MatrixBase &L, MatrixBase &D, MatrixBase &U)
const
486 if (is_a<DenseMatrix>(L) and is_a<DenseMatrix>(D)
487 and is_a<DenseMatrix>(U)) {
488 DenseMatrix &L_ = down_cast<DenseMatrix &>(L);
489 DenseMatrix &D_ = down_cast<DenseMatrix &>(D);
490 DenseMatrix &U_ = down_cast<DenseMatrix &>(U);
491 fraction_free_LDU(*
this, L_, D_, U_);
496 void DenseMatrix::QR(MatrixBase &Q, MatrixBase &R)
const
498 if (is_a<DenseMatrix>(Q) and is_a<DenseMatrix>(R)) {
499 DenseMatrix &Q_ = down_cast<DenseMatrix &>(Q);
500 DenseMatrix &R_ = down_cast<DenseMatrix &>(R);
501 SymEngine::QR(*
this, Q_, R_);
506 void DenseMatrix::cholesky(MatrixBase &L)
const
508 if (is_a<DenseMatrix>(L)) {
509 DenseMatrix &L_ = down_cast<DenseMatrix &>(L);
510 SymEngine::cholesky(*
this, L_);
516 void jacobian(
const DenseMatrix &A,
const DenseMatrix &x, DenseMatrix &result,
519 SYMENGINE_ASSERT(A.col_ == 1);
520 SYMENGINE_ASSERT(x.col_ == 1);
521 SYMENGINE_ASSERT(A.row_ == result.nrows() and x.row_ == result.ncols());
523 #pragma omp parallel for
524 for (
unsigned i = 0; i < result.row_; i++) {
525 for (
unsigned j = 0; j < result.col_; j++) {
526 if (is_a<Symbol>(*(x.m_[j]))) {
527 const RCP<const Symbol> x_
528 = rcp_static_cast<const Symbol>(x.m_[j]);
529 result.m_[i * result.col_ + j] = A.m_[i]->diff(x_, diff_cache);
537 throw SymEngineException(
538 "'x' must contain Symbols only. "
539 "Use sjacobian for SymPy style differentiation");
543 void sjacobian(
const DenseMatrix &A,
const DenseMatrix &x, DenseMatrix &result,
546 SYMENGINE_ASSERT(A.col_ == 1);
547 SYMENGINE_ASSERT(x.col_ == 1);
548 SYMENGINE_ASSERT(A.row_ == result.nrows() and x.row_ == result.ncols());
549 #pragma omp parallel for
550 for (
unsigned i = 0; i < result.row_; i++) {
551 for (
unsigned j = 0; j < result.col_; j++) {
552 if (is_a<Symbol>(*(x.m_[j]))) {
553 const RCP<const Symbol> x_
554 = rcp_static_cast<const Symbol>(x.m_[j]);
555 result.m_[i * result.col_ + j] = A.m_[i]->diff(x_, diff_cache);
558 const RCP<const Symbol> x_ =
symbol(
"x_");
559 result.m_[i * result.col_ + j] =
ssubs(
560 ssubs(A.m_[i], {{x.m_[j], x_}})->diff(x_, diff_cache),
569 void diff(
const DenseMatrix &A,
const RCP<const Symbol> &x, DenseMatrix &result,
572 SYMENGINE_ASSERT(A.row_ == result.nrows() and A.col_ == result.ncols());
573 #pragma omp parallel for
574 for (
unsigned i = 0; i < result.row_; i++) {
575 for (
unsigned j = 0; j < result.col_; j++) {
576 result.m_[i * result.col_ + j]
577 = A.m_[i * result.col_ + j]->diff(x, diff_cache);
582 void sdiff(
const DenseMatrix &A,
const RCP<const Basic> &x, DenseMatrix &result,
585 SYMENGINE_ASSERT(A.row_ == result.nrows() and A.col_ == result.ncols());
586 #pragma omp parallel for
587 for (
unsigned i = 0; i < result.row_; i++) {
588 for (
unsigned j = 0; j < result.col_; j++) {
589 if (is_a<Symbol>(*x)) {
590 const RCP<const Symbol> x_ = rcp_static_cast<const Symbol>(x);
591 result.m_[i * result.col_ + j]
592 = A.m_[i * result.col_ + j]->diff(x_, diff_cache);
595 const RCP<const Symbol> x_ =
symbol(
"_x");
596 result.m_[i * result.col_ + j]
597 =
ssubs(
ssubs(A.m_[i * result.col_ + j], {{x, x_}})
598 ->diff(x_, diff_cache),
606 void conjugate_dense(
const DenseMatrix &A, DenseMatrix &B)
608 SYMENGINE_ASSERT(B.col_ == A.col_ and B.row_ == A.row_);
610 for (
unsigned i = 0; i < A.row_; i++)
611 for (
unsigned j = 0; j < A.col_; j++)
612 B.m_[i * B.col_ + j] =
conjugate(A.m_[i * A.col_ + j]);
616 void transpose_dense(
const DenseMatrix &A, DenseMatrix &B)
618 SYMENGINE_ASSERT(B.row_ == A.col_ and B.col_ == A.row_);
620 for (
unsigned i = 0; i < A.row_; i++)
621 for (
unsigned j = 0; j < A.col_; j++)
622 B.m_[j * B.col_ + i] = A.m_[i * A.col_ + j];
626 void conjugate_transpose_dense(
const DenseMatrix &A, DenseMatrix &B)
628 SYMENGINE_ASSERT(B.row_ == A.col_ and B.col_ == A.row_);
630 for (
unsigned i = 0; i < A.row_; i++)
631 for (
unsigned j = 0; j < A.col_; j++)
632 B.m_[j * B.col_ + i] =
conjugate(A.m_[i * A.col_ + j]);
636 void submatrix_dense(
const DenseMatrix &A, DenseMatrix &B,
unsigned row_start,
637 unsigned col_start,
unsigned row_end,
unsigned col_end,
638 unsigned row_step,
unsigned col_step)
640 SYMENGINE_ASSERT(row_end >= row_start and col_end >= col_start);
641 SYMENGINE_ASSERT(row_end < A.row_);
642 SYMENGINE_ASSERT(col_end < A.col_);
643 SYMENGINE_ASSERT(B.row_ == row_end - row_start + 1
644 and B.col_ == col_end - col_start + 1);
646 unsigned row = B.row_, col = B.col_;
648 for (
unsigned i = 0; i < row; i += row_step)
649 for (
unsigned j = 0; j < col; j += col_step)
650 B.m_[i * col + j] = A.m_[(row_start + i) * A.col_ + col_start + j];
654 void add_dense_dense(
const DenseMatrix &A,
const DenseMatrix &B, DenseMatrix &C)
656 SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_ and A.row_ == C.row_
657 and A.col_ == C.col_);
659 unsigned row = A.row_, col = A.col_;
661 for (
unsigned i = 0; i < row; i++) {
662 for (
unsigned j = 0; j < col; j++) {
663 C.m_[i * col + j] =
add(A.m_[i * col + j], B.m_[i * col + j]);
668 void add_dense_scalar(
const DenseMatrix &A,
const RCP<const Basic> &k,
671 SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
673 unsigned row = A.row_, col = A.col_;
675 for (
unsigned i = 0; i < row; i++) {
676 for (
unsigned j = 0; j < col; j++) {
677 B.m_[i * col + j] =
add(A.m_[i * col + j], k);
683 void mul_dense_dense(
const DenseMatrix &A,
const DenseMatrix &B, DenseMatrix &C)
685 SYMENGINE_ASSERT(A.col_ == B.row_ and C.row_ == A.row_
686 and C.col_ == B.col_);
688 unsigned row = A.row_, col = B.col_;
690 if (&A != &C and &B != &C) {
691 for (
unsigned r = 0; r < row; r++) {
692 for (
unsigned c = 0; c < col; c++) {
693 C.m_[r * col + c] = zero;
694 for (
unsigned k = 0; k < A.col_; k++)
696 =
add(C.m_[r * col + c],
697 mul(A.m_[r * A.col_ + k], B.m_[k * col + c]));
701 DenseMatrix tmp = DenseMatrix(A.row_, B.col_);
702 mul_dense_dense(A, B, tmp);
707 void elementwise_mul_dense_dense(
const DenseMatrix &A,
const DenseMatrix &B,
710 SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_ and A.row_ == C.row_
711 and A.col_ == C.col_);
713 unsigned row = A.row_, col = A.col_;
715 for (
unsigned i = 0; i < row; i++) {
716 for (
unsigned j = 0; j < col; j++) {
717 C.m_[i * col + j] =
mul(A.m_[i * col + j], B.m_[i * col + j]);
722 void mul_dense_scalar(
const DenseMatrix &A,
const RCP<const Basic> &k,
725 SYMENGINE_ASSERT(A.col_ == B.col_ and A.row_ == B.row_);
727 unsigned row = A.row_, col = A.col_;
729 for (
unsigned i = 0; i < row; i++) {
730 for (
unsigned j = 0; j < col; j++) {
731 B.m_[i * col + j] =
mul(A.m_[i * col + j], k);
737 void DenseMatrix::row_join(
const DenseMatrix &B)
739 this->col_insert(B, col_);
742 void DenseMatrix::col_join(
const DenseMatrix &B)
744 this->row_insert(B, row_);
747 void DenseMatrix::row_insert(
const DenseMatrix &B,
unsigned pos)
749 SYMENGINE_ASSERT(col_ == B.col_ and pos <= row_)
751 unsigned row = row_, col = col_;
752 this->resize(row_ + B.row_, col_);
754 for (
unsigned i = row; i-- > pos;) {
755 for (
unsigned j = col; j-- > 0;) {
756 this->m_[(i + B.row_) * col + j] = this->m_[i * col + j];
760 for (
unsigned i = 0; i < B.row_; i++) {
761 for (
unsigned j = 0; j < col; j++) {
762 this->m_[(i + pos) * col + j] = B.m_[i * col + j];
767 void DenseMatrix::col_insert(
const DenseMatrix &B,
unsigned pos)
769 SYMENGINE_ASSERT(row_ == B.row_ and pos <= col_)
771 unsigned row = row_, col = col_;
772 this->resize(row_, col_ + B.col_);
774 for (
unsigned i = row; i-- > 0;) {
775 for (
unsigned j = col; j-- > 0;) {
777 this->m_[i * (col + B.col_) + j + B.col_]
778 = this->m_[i * col + j];
780 this->m_[i * (col + B.col_) + j] = this->m_[i * col + j];
785 for (
unsigned i = 0; i < row; i++) {
786 for (
unsigned j = 0; j < B.col_; j++) {
787 this->m_[i * (col + B.col_) + j + pos] = B.m_[i * B.col_ + j];
792 void DenseMatrix::row_del(
unsigned k)
794 SYMENGINE_ASSERT(k < row_)
799 for (
unsigned i = k; i < row_ - 1; i++) {
800 row_exchange_dense(*
this, i, i + 1);
802 this->resize(row_ - 1, col_);
806 void DenseMatrix::col_del(
unsigned k)
808 SYMENGINE_ASSERT(k < col_)
813 unsigned row = row_, col = col_, m = 0;
815 for (
unsigned i = 0; i < row; i++) {
816 for (
unsigned j = 0; j < col; j++) {
818 this->m_[m] = this->m_[i * col + j];
823 this->resize(row_, col_ - 1);
828 void row_exchange_dense(DenseMatrix &A,
unsigned i,
unsigned j)
830 SYMENGINE_ASSERT(i != j and i < A.row_ and j < A.row_);
832 unsigned col = A.col_;
834 for (
unsigned k = 0; k < A.col_; k++)
835 std::swap(A.m_[i * col + k], A.m_[j * col + k]);
838 void row_mul_scalar_dense(DenseMatrix &A,
unsigned i, RCP<const Basic> &c)
840 SYMENGINE_ASSERT(i < A.row_);
842 unsigned col = A.col_;
844 for (
unsigned j = 0; j < A.col_; j++)
845 A.m_[i * col + j] =
mul(c, A.m_[i * col + j]);
848 void row_add_row_dense(DenseMatrix &A,
unsigned i,
unsigned j,
851 SYMENGINE_ASSERT(i != j and i < A.row_ and j < A.row_);
853 unsigned col = A.col_;
855 for (
unsigned k = 0; k < A.col_; k++)
856 A.m_[i * col + k] =
add(A.m_[i * col + k],
mul(c, A.m_[j * col + k]));
859 void permuteFwd(DenseMatrix &A, permutelist &pl)
862 row_exchange_dense(A, p.first, p.second);
867 void column_exchange_dense(DenseMatrix &A,
unsigned i,
unsigned j)
869 SYMENGINE_ASSERT(i != j and i < A.col_ and j < A.col_);
871 unsigned col = A.col_;
873 for (
unsigned k = 0; k < A.row_; k++)
874 std::swap(A.m_[k * col + i], A.m_[k * col + j]);
878 void pivoted_gaussian_elimination(
const DenseMatrix &A, DenseMatrix &B,
881 SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
883 unsigned row = A.row_, col = A.col_;
884 unsigned index = 0, i, j, k;
887 RCP<const Basic> scale;
889 for (i = 0; i < col - 1; i++) {
893 k = pivot(B, index, i);
897 row_exchange_dense(B, k, index);
898 pl.push_back({k, index});
901 scale =
div(one, B.m_[index * col + i]);
902 row_mul_scalar_dense(B, index, scale);
904 for (j = i + 1; j < row; j++) {
905 for (k = i + 1; k < col; k++)
907 =
sub(B.m_[j * col + k],
908 mul(B.m_[j * col + i], B.m_[i * col + k]));
909 B.m_[j * col + i] = zero;
919 void fraction_free_gaussian_elimination(
const DenseMatrix &A, DenseMatrix &B)
921 SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
923 unsigned col = A.col_;
926 for (
unsigned i = 0; i < col - 1; i++)
927 for (
unsigned j = i + 1; j < A.row_; j++) {
928 for (
unsigned k = i + 1; k < col; k++) {
930 =
sub(
mul(B.m_[i * col + i], B.m_[j * col + k]),
931 mul(B.m_[j * col + i], B.m_[i * col + k]));
934 =
div(B.m_[j * col + k], B.m_[i * col - col + i - 1]);
936 B.m_[j * col + i] = zero;
941 void pivoted_fraction_free_gaussian_elimination(
const DenseMatrix &A,
942 DenseMatrix &B, permutelist &pl)
944 SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
946 unsigned col = A.col_, row = A.row_;
947 unsigned index = 0, i, k, j;
950 for (i = 0; i < col - 1; i++) {
954 k = pivot(B, index, i);
958 row_exchange_dense(B, k, index);
959 pl.push_back({k, index});
962 for (j = i + 1; j < row; j++) {
963 for (k = i + 1; k < col; k++) {
965 =
sub(
mul(B.m_[i * col + i], B.m_[j * col + k]),
966 mul(B.m_[j * col + i], B.m_[i * col + k]));
969 =
div(B.m_[j * col + k], B.m_[i * col - col + i - 1]);
971 B.m_[j * col + i] = zero;
978 void pivoted_gauss_jordan_elimination(
const DenseMatrix &A, DenseMatrix &B,
981 SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
983 unsigned row = A.row_, col = A.col_;
984 unsigned index = 0, i, j, k;
985 RCP<const Basic> scale;
988 for (i = 0; i < col; i++) {
992 k = pivot(B, index, i);
996 row_exchange_dense(B, k, index);
997 pl.push_back({k, index});
1000 scale =
div(one, B.m_[index * col + i]);
1001 row_mul_scalar_dense(B, index, scale);
1003 for (j = 0; j < row; j++) {
1007 scale =
mul(minus_one, B.m_[j * col + i]);
1008 row_add_row_dense(B, j, index, scale);
1018 void fraction_free_gauss_jordan_elimination(
const DenseMatrix &A,
1021 SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
1023 unsigned row = A.row_, col = A.col_;
1029 for (i = 0; i < col; i++) {
1031 d = B.m_[i * col - col + i - 1];
1032 for (j = 0; j < row; j++)
1034 for (k = 0; k < col; k++) {
1037 =
sub(
mul(B.m_[i * col + i], B.m_[j * col + k]),
1038 mul(B.m_[j * col + i], B.m_[i * col + k]));
1040 B.m_[j * col + k] =
div(B.m_[j * col + k], d);
1043 for (j = 0; j < row; j++)
1045 B.m_[j * col + i] = zero;
1049 void pivoted_fraction_free_gauss_jordan_elimination(
const DenseMatrix &A,
1053 SYMENGINE_ASSERT(A.row_ == B.row_ and A.col_ == B.col_);
1055 unsigned row = A.row_, col = A.col_;
1056 unsigned index = 0, i, k, j;
1061 for (i = 0; i < col; i++) {
1065 k = pivot(B, index, i);
1069 row_exchange_dense(B, k, index);
1070 pl.push_back({k, index});
1073 for (j = 0; j < row; j++) {
1076 for (k = 0; k < col; k++) {
1080 =
sub(
mul(B.m_[index * col + i], B.m_[j * col + k]),
1081 mul(B.m_[j * col + i], B.m_[index * col + k]));
1084 B.m_[j * col + k] =
div(B.m_[j * col + k], d);
1088 d = B.m_[index * col + i];
1090 for (j = 0; j < row; j++) {
1093 B.m_[j * col + i] = zero;
1099 unsigned pivot(DenseMatrix &B,
unsigned r,
unsigned c)
1101 for (
unsigned k = r; k < B.row_; k++) {
1102 if (!is_true(
is_zero(*(B.m_[k * B.col_ + c])))) {
1109 void reduced_row_echelon_form(
const DenseMatrix &A, DenseMatrix &b,
1110 vec_uint &pivot_cols,
bool normalize_last)
1113 if (normalize_last) {
1114 pivoted_fraction_free_gauss_jordan_elimination(A, b, pl);
1116 pivoted_gauss_jordan_elimination(A, b, pl);
1119 for (
unsigned col = 0; col < b.col_ && row < b.row_; col++) {
1120 if (is_true(
is_zero(*b.get(row, col))))
1122 pivot_cols.push_back(col);
1123 if (row == 0 and normalize_last) {
1124 RCP<const Basic> m =
div(one, b.get(row, col));
1133 void diagonal_solve(
const DenseMatrix &A,
const DenseMatrix &b, DenseMatrix &x)
1135 SYMENGINE_ASSERT(A.row_ == A.col_);
1136 SYMENGINE_ASSERT(x.row_ == A.col_ and x.col_ == b.col_);
1138 const unsigned sys = b.col_;
1141 for (
unsigned k = 0; k < sys; k++) {
1142 for (
unsigned i = 0; i < A.col_; i++) {
1143 x.m_[i * sys + k] =
div(b.m_[i * sys + k], A.m_[i * A.col_ + i]);
1148 void back_substitution(
const DenseMatrix &U,
const DenseMatrix &b,
1151 SYMENGINE_ASSERT(U.row_ == U.col_);
1152 SYMENGINE_ASSERT(b.row_ == U.row_);
1153 SYMENGINE_ASSERT(x.row_ == U.col_ and x.col_ == b.col_);
1155 const unsigned col = U.col_;
1156 const unsigned sys = b.col_;
1159 for (
unsigned k = 0; k < sys; k++) {
1160 for (
int i = col - 1; i >= 0; i--) {
1161 for (
unsigned j = i + 1; j < col; j++)
1163 =
sub(x.m_[i * sys + k],
1164 mul(U.m_[i * col + j], x.m_[j * sys + k]));
1165 x.m_[i * sys + k] =
div(x.m_[i * sys + k], U.m_[i * col + i]);
1170 void forward_substitution(
const DenseMatrix &A,
const DenseMatrix &b,
1173 SYMENGINE_ASSERT(A.row_ == A.col_);
1174 SYMENGINE_ASSERT(b.row_ == A.row_);
1175 SYMENGINE_ASSERT(x.row_ == A.col_ and x.col_ == b.col_);
1177 unsigned col = A.col_;
1178 const unsigned sys = b.col_;
1181 for (
unsigned k = 0; k < b.col_; k++) {
1182 for (
unsigned i = 0; i < col - 1; i++) {
1183 for (
unsigned j = i + 1; j < col; j++) {
1185 =
sub(
mul(A.m_[i * col + i], x.m_[j * sys + k]),
1186 mul(A.m_[j * col + i], x.m_[i * sys + k]));
1189 =
div(x.m_[j * sys + k], A.m_[(i - 1) * col + i - 1]);
1195 void fraction_free_gaussian_elimination_solve(
const DenseMatrix &A,
1196 const DenseMatrix &b,
1199 SYMENGINE_ASSERT(A.row_ == A.col_);
1200 SYMENGINE_ASSERT(b.row_ == A.row_ and x.row_ == A.row_);
1201 SYMENGINE_ASSERT(x.col_ == b.col_);
1203 int i, j, k, col = A.col_, bcol = b.col_;
1204 DenseMatrix A_ = DenseMatrix(A.row_, A.col_, A.m_);
1205 DenseMatrix b_ = DenseMatrix(b.row_, b.col_, b.m_);
1207 for (i = 0; i < col - 1; i++)
1208 for (j = i + 1; j < col; j++) {
1209 for (k = 0; k < bcol; k++) {
1211 =
sub(
mul(A_.m_[i * col + i], b_.m_[j * bcol + k]),
1212 mul(A_.m_[j * col + i], b_.m_[i * bcol + k]));
1214 b_.m_[j * bcol + k] =
div(b_.m_[j * bcol + k],
1215 A_.m_[i * col - col + i - 1]);
1218 for (k = i + 1; k < col; k++) {
1220 =
sub(
mul(A_.m_[i * col + i], A_.m_[j * col + k]),
1221 mul(A_.m_[j * col + i], A_.m_[i * col + k]));
1224 =
div(A_.m_[j * col + k], A_.m_[i * col - col + i - 1]);
1226 A_.m_[j * col + i] = zero;
1229 for (i = 0; i < col * bcol; i++)
1232 for (k = 0; k < bcol; k++) {
1233 for (i = col - 1; i >= 0; i--) {
1234 for (j = i + 1; j < col; j++)
1236 =
sub(b_.m_[i * bcol + k],
1237 mul(A_.m_[i * col + j], x.m_[j * bcol + k]));
1238 x.m_[i * bcol + k] =
div(b_.m_[i * bcol + k], A_.m_[i * col + i]);
1243 void fraction_free_gauss_jordan_solve(
const DenseMatrix &A,
1244 const DenseMatrix &b, DenseMatrix &x,
1247 SYMENGINE_ASSERT(A.row_ == A.col_);
1248 SYMENGINE_ASSERT(b.row_ == A.row_ and x.row_ == A.row_);
1249 SYMENGINE_ASSERT(x.col_ == b.col_);
1251 unsigned i, j, k, p, col = A.col_, bcol = b.col_;
1253 DenseMatrix A_ = DenseMatrix(A.row_, A.col_, A.m_);
1254 DenseMatrix b_ = DenseMatrix(b.row_, b.col_, b.m_);
1256 for (i = 0; i < col; i++) {
1258 d = A_.m_[i * col - col + i - 1];
1261 while (p < col and
eq(*A_.m_[p * col + i], *zero)) {
1264 SYMENGINE_ASSERT(p != col);
1267 for (k = i; k < col; k++) {
1268 std::swap(A_.m_[p * col + k], A_.m_[i * col + k]);
1270 for (k = 0; k < bcol; k++) {
1271 std::swap(b_.m_[p * bcol + k], b_.m_[i * bcol + k]);
1275 for (j = 0; j < col; j++)
1277 for (k = 0; k < bcol; k++) {
1279 =
sub(
mul(A_.m_[i * col + i], b_.m_[j * bcol + k]),
1280 mul(A_.m_[j * col + i], b_.m_[i * bcol + k]));
1282 b_.m_[j * bcol + k] =
div(b_.m_[j * bcol + k], d);
1285 for (k = 0; k < col; k++) {
1288 =
sub(
mul(A_.m_[i * col + i], A_.m_[j * col + k]),
1289 mul(A_.m_[j * col + i], A_.m_[i * col + k]));
1291 A_.m_[j * col + k] =
div(A_.m_[j * col + k], d);
1296 for (j = 0; j < col; j++)
1298 A_.m_[j * col + i] = zero;
1302 for (k = 0; k < bcol; k++)
1303 for (i = 0; i < col; i++)
1304 x.m_[i * bcol + k] =
div(b_.m_[i * bcol + k], A_.m_[i * col + i]);
1307 void fraction_free_LU_solve(
const DenseMatrix &A,
const DenseMatrix &b,
1310 DenseMatrix LU = DenseMatrix(A.nrows(), A.ncols());
1311 DenseMatrix x_ = DenseMatrix(b.nrows(), b.ncols());
1313 fraction_free_LU(A, LU);
1314 forward_substitution(LU, b, x_);
1315 back_substitution(LU, x_, x);
1318 void LU_solve(
const DenseMatrix &A,
const DenseMatrix &b, DenseMatrix &x)
1320 DenseMatrix L = DenseMatrix(A.nrows(), A.ncols());
1321 DenseMatrix U = DenseMatrix(A.nrows(), A.ncols());
1322 DenseMatrix x_ = DenseMatrix(b.nrows(), b.ncols());
1325 forward_substitution(L, b, x_);
1326 back_substitution(U, x_, x);
1329 void pivoted_LU_solve(
const DenseMatrix &A,
const DenseMatrix &b,
1332 DenseMatrix L = DenseMatrix(A.nrows(), A.ncols());
1333 DenseMatrix U = DenseMatrix(A.nrows(), A.ncols());
1334 DenseMatrix x_ = DenseMatrix(b);
1337 pivoted_LU(A, L, U, pl);
1339 forward_substitution(L, x_, x_);
1340 back_substitution(U, x_, x);
1343 void LDL_solve(
const DenseMatrix &A,
const DenseMatrix &b, DenseMatrix &x)
1345 DenseMatrix L = DenseMatrix(A.nrows(), A.ncols());
1346 DenseMatrix D = DenseMatrix(A.nrows(), A.ncols());
1347 DenseMatrix x_ = DenseMatrix(b.nrows(), b.ncols());
1349 if (not is_symmetric_dense(A))
1350 throw SymEngineException(
"Matrix must be symmetric");
1353 forward_substitution(L, b, x);
1354 diagonal_solve(D, x, x_);
1355 transpose_dense(L, D);
1356 back_substitution(D, x_, x);
1367 void fraction_free_LU(
const DenseMatrix &A, DenseMatrix &LU)
1369 SYMENGINE_ASSERT(A.row_ == A.col_ and LU.row_ == LU.col_
1370 and A.row_ == LU.row_);
1372 unsigned n = A.row_;
1377 for (i = 0; i < n - 1; i++)
1378 for (j = i + 1; j < n; j++)
1379 for (k = i + 1; k < n; k++) {
1380 LU.m_[j * n + k] =
sub(
mul(LU.m_[i * n + i], LU.m_[j * n + k]),
1381 mul(LU.m_[j * n + i], LU.m_[i * n + k]));
1384 =
div(LU.m_[j * n + k], LU.m_[i * n - n + i - 1]);
1391 void LU(
const DenseMatrix &A, DenseMatrix &L, DenseMatrix &U)
1393 SYMENGINE_ASSERT(A.row_ == A.col_ and L.row_ == L.col_
1394 and U.row_ == U.col_);
1395 SYMENGINE_ASSERT(A.row_ == L.row_ and A.row_ == U.row_);
1397 unsigned n = A.row_;
1399 RCP<const Basic> scale;
1403 for (j = 0; j < n; j++) {
1404 for (i = 0; i < j; i++)
1405 for (k = 0; k < i; k++)
1406 U.m_[i * n + j] =
sub(U.m_[i * n + j],
1407 mul(U.m_[i * n + k], U.m_[k * n + j]));
1409 for (i = j; i < n; i++) {
1410 for (k = 0; k < j; k++)
1411 U.m_[i * n + j] =
sub(U.m_[i * n + j],
1412 mul(U.m_[i * n + k], U.m_[k * n + j]));
1415 scale =
div(one, U.m_[j * n + j]);
1417 for (i = j + 1; i < n; i++)
1418 U.m_[i * n + j] =
mul(U.m_[i * n + j], scale);
1421 for (i = 0; i < n; i++) {
1422 for (j = 0; j < i; j++) {
1423 L.m_[i * n + j] = U.m_[i * n + j];
1424 U.m_[i * n + j] = zero;
1426 L.m_[i * n + i] = one;
1427 for (j = i + 1; j < n; j++)
1428 L.m_[i * n + j] = zero;
1435 void pivoted_LU(
const DenseMatrix &A, DenseMatrix &LU, permutelist &pl)
1437 SYMENGINE_ASSERT(A.row_ == A.col_ and LU.row_ == LU.col_);
1438 SYMENGINE_ASSERT(A.row_ == LU.row_);
1440 unsigned n = A.row_;
1442 RCP<const Basic> scale;
1447 for (j = 0; j < n; j++) {
1448 for (i = 0; i < j; i++)
1449 for (k = 0; k < i; k++)
1450 LU.m_[i * n + j] =
sub(LU.m_[i * n + j],
1451 mul(LU.m_[i * n + k], LU.m_[k * n + j]));
1453 for (i = j; i < n; i++) {
1454 for (k = 0; k < j; k++) {
1455 LU.m_[i * n + j] =
sub(LU.m_[i * n + j],
1456 mul(LU.m_[i * n + k], LU.m_[k * n + j]));
1458 if (pivot == -1 and !is_true(
is_zero(*LU.m_[i * n + j])))
1462 throw SymEngineException(
"Matrix is rank deficient");
1463 if (pivot - j != 0) {
1464 row_exchange_dense(LU, pivot, j);
1465 pl.push_back({pivot, j});
1467 scale =
div(one, LU.m_[j * n + j]);
1469 for (i = j + 1; i < n; i++)
1470 LU.m_[i * n + j] =
mul(LU.m_[i * n + j], scale);
1477 void pivoted_LU(
const DenseMatrix &A, DenseMatrix &L, DenseMatrix &U,
1480 SYMENGINE_ASSERT(A.row_ == A.col_ and L.row_ == L.col_
1481 and U.row_ == U.col_);
1482 SYMENGINE_ASSERT(A.row_ == L.row_ and A.row_ == U.row_);
1484 pivoted_LU(A, U, pl);
1485 unsigned n = A.col_;
1486 for (
unsigned i = 0; i < n; i++) {
1487 for (
unsigned j = 0; j < i; j++) {
1488 L.m_[i * n + j] = U.m_[i * n + j];
1489 U.m_[i * n + j] = zero;
1491 L.m_[i * n + i] = one;
1492 for (
unsigned j = i + 1; j < n; j++)
1493 L.m_[i * n + j] = zero;
1502 void fraction_free_LDU(
const DenseMatrix &A, DenseMatrix &L, DenseMatrix &D,
1505 SYMENGINE_ASSERT(A.row_ == L.row_ and A.row_ == U.row_);
1506 SYMENGINE_ASSERT(A.col_ == L.col_ and A.col_ == U.col_);
1508 unsigned row = A.row_, col = A.col_;
1510 RCP<const Basic> old =
integer(1);
1515 for (i = 0; i < row; i++)
1516 for (j = 0; j < row; j++)
1518 L.m_[i * col + j] = zero;
1520 L.m_[i * col + i] = one;
1523 for (i = 0; i < row * row; i++)
1526 for (k = 0; k < row - 1; k++) {
1527 L.m_[k * col + k] = U.m_[k * col + k];
1528 D.m_[k * col + k] =
mul(old, U.m_[k * col + k]);
1530 for (i = k + 1; i < row; i++) {
1531 L.m_[i * col + k] = U.m_[i * col + k];
1532 for (j = k + 1; j < col; j++)
1534 =
div(
sub(
mul(U.m_[k * col + k], U.m_[i * col + j]),
1535 mul(U.m_[k * col + j], U.m_[i * col + k])),
1537 U.m_[i * col + k] = zero;
1540 old = U.m_[k * col + k];
1543 D.m_[row * col - col + row - 1] = old;
1548 void QR(
const DenseMatrix &A, DenseMatrix &Q, DenseMatrix &R)
1550 unsigned row = A.row_;
1551 unsigned col = A.col_;
1553 SYMENGINE_ASSERT(Q.row_ == row and Q.col_ == col and R.row_ == col
1561 for (i = 0; i < row * col; i++)
1565 for (i = 0; i < col * col; i++)
1568 for (j = 0; j < col; j++) {
1570 for (k = 0; k < row; k++)
1571 tmp[k] = A.m_[k * col + j];
1573 for (i = 0; i < j; i++) {
1575 for (k = 0; k < row; k++)
1576 t =
add(t,
mul(A.m_[k * col + j], Q.m_[k * col + i]));
1577 for (k = 0; k < row; k++)
1578 tmp[k] =
expand(
sub(tmp[k],
mul(Q.m_[k * col + i], t)));
1583 for (k = 0; k < row; k++)
1588 R.m_[j * col + j] = t;
1589 for (k = 0; k < row; k++)
1590 Q.m_[k * col + j] =
div(tmp[k], t);
1592 for (i = 0; i < j; i++) {
1594 for (k = 0; k < row; k++)
1595 t =
add(t,
mul(Q.m_[k * col + i], A.m_[k * col + j]));
1596 R.m_[i * col + j] = t;
1603 void LDL(
const DenseMatrix &A, DenseMatrix &L, DenseMatrix &D)
1605 SYMENGINE_ASSERT(A.row_ == A.col_);
1606 SYMENGINE_ASSERT(L.row_ == A.row_ and L.col_ == A.row_);
1607 SYMENGINE_ASSERT(D.row_ == A.row_ and D.col_ == A.row_);
1609 unsigned col = A.col_;
1611 RCP<const Basic> sum;
1612 RCP<const Basic> i2 =
integer(2);
1615 for (i = 0; i < col; i++)
1616 for (j = 0; j < col; j++)
1617 D.m_[i * col + j] = zero;
1620 for (i = 0; i < col; i++)
1621 for (j = 0; j < col; j++)
1622 L.m_[i * col + j] = (i != j) ? zero : one;
1624 for (i = 0; i < col; i++) {
1625 for (j = 0; j < i; j++) {
1627 for (k = 0; k < j; k++)
1628 sum =
add(sum,
mul(
mul(L.m_[i * col + k], L.m_[j * col + k]),
1629 D.m_[k * col + k]));
1631 =
mul(
div(one, D.m_[j * col + j]),
sub(A.m_[i * col + j], sum));
1634 for (k = 0; k < i; k++)
1635 sum =
add(sum,
mul(
pow(L.m_[i * col + k], i2), D.m_[k * col + k]));
1636 D.m_[i * col + i] =
sub(A.m_[i * col + i], sum);
1641 void cholesky(
const DenseMatrix &A, DenseMatrix &L)
1643 SYMENGINE_ASSERT(A.row_ == A.col_);
1644 SYMENGINE_ASSERT(L.row_ == A.row_ and L.col_ == A.row_);
1646 unsigned col = A.col_;
1648 RCP<const Basic> sum;
1649 RCP<const Basic> i2 =
integer(2);
1650 RCP<const Basic> half =
div(one, i2);
1653 for (i = 0; i < col; i++)
1654 for (j = 0; j < col; j++)
1655 L.m_[i * col + j] = zero;
1657 for (i = 0; i < col; i++) {
1658 for (j = 0; j < i; j++) {
1660 for (k = 0; k < j; k++)
1661 sum =
add(sum,
mul(L.m_[i * col + k], L.m_[j * col + k]));
1664 =
mul(
div(one, L.m_[j * col + j]),
sub(A.m_[i * col + j], sum));
1667 for (k = 0; k < i; k++)
1668 sum =
add(sum,
pow(L.m_[i * col + k], i2));
1669 L.m_[i * col + i] =
pow(
sub(A.m_[i * col + i], sum), half);
1674 bool is_symmetric_dense(
const DenseMatrix &A)
1676 if (A.col_ != A.row_)
1679 unsigned col = A.col_;
1682 for (
unsigned i = 0; i < col; i++)
1683 for (
unsigned j = i + 1; j < col; j++)
1684 if (not
eq(*(A.m_[j * col + i]), *(A.m_[i * col + j]))) {
1693 RCP<const Basic> det_bareis(
const DenseMatrix &A)
1695 SYMENGINE_ASSERT(A.row_ == A.col_);
1697 unsigned n = A.row_;
1701 }
else if (n == 2) {
1703 return sub(
mul(A.m_[0], A.m_[3]),
mul(A.m_[1], A.m_[2]));
1704 }
else if (n == 3) {
1708 mul(
mul(A.m_[1], A.m_[5]), A.m_[6])),
1709 mul(
mul(A.m_[2], A.m_[3]), A.m_[7])),
1711 mul(
mul(A.m_[1], A.m_[3]), A.m_[8])),
1712 mul(
mul(A.m_[0], A.m_[5]), A.m_[7])));
1715 if (A.is_lower() or A.is_upper()) {
1716 RCP<const Basic> det = A.m_[0];
1717 for (
unsigned i = 1; i < n; ++i) {
1718 det =
mul(det, A.m_[i * n + i]);
1723 DenseMatrix B = DenseMatrix(n, n, A.m_);
1724 unsigned i,
sign = 1;
1727 for (
unsigned k = 0; k < n - 1; k++) {
1728 if (is_true(
is_zero(*B.m_[k * n + k]))) {
1729 for (i = k + 1; i < n; i++)
1730 if (!is_true(
is_zero(*B.m_[i * n + k]))) {
1731 row_exchange_dense(B, i, k);
1739 for (i = k + 1; i < n; i++) {
1740 for (
unsigned j = k + 1; j < n; j++) {
1741 d =
sub(
mul(B.m_[k * n + k], B.m_[i * n + j]),
1742 mul(B.m_[i * n + k], B.m_[k * n + j]));
1744 d =
div(d, B.m_[(k - 1) * n + k - 1]);
1745 B.m_[i * n + j] = d;
1750 return (sign == 1) ? B.m_[n * n - 1] :
mul(minus_one, B.m_[n * n - 1]);
1760 SYMENGINE_ASSERT(A.row_ == A.col_);
1762 unsigned col = A.col_;
1763 unsigned i, k, l, m;
1769 for (
unsigned n = col; n > 1; n--) {
1772 DenseMatrix T = DenseMatrix(n + 1, n);
1773 DenseMatrix C = DenseMatrix(k, 1);
1776 for (i = 0; i < n * (n + 1); i++)
1778 for (i = 0; i < k; i++)
1779 C.m_[i] = A.m_[i * col + k];
1782 for (i = 0; i < n - 2; i++) {
1783 DenseMatrix B = DenseMatrix(k, 1);
1784 for (l = 0; l < k; l++) {
1786 for (m = 0; m < k; m++)
1788 =
add(B.m_[l],
mul(A.m_[l * col + m], items[i].m_[m]));
1794 for (i = 0; i < n - 1; i++) {
1795 RCP<const Basic> element = zero;
1796 for (l = 0; l < k; l++)
1797 element =
add(element,
mul(A.m_[k * col + l], items[i].m_[l]));
1800 items_.
insert(items_.
begin(),
mul(minus_one, A.m_[k * col + k]));
1803 for (i = 0; i < n; i++) {
1804 for (l = 0; l < n - i + 1; l++)
1805 T.m_[(i + l) * n + i] = items_[l];
1811 polys.
push_back(DenseMatrix(2, 1, {one,
mul(A.m_[0], minus_one)}));
1813 for (i = 0; i < col - 1; i++) {
1814 unsigned t_row = transforms[col - 2 - i].nrows();
1815 unsigned t_col = transforms[col - 2 - i].ncols();
1816 DenseMatrix B = DenseMatrix(t_row, 1);
1818 for (l = 0; l < t_row; l++) {
1820 for (m = 0; m < t_col; m++) {
1821 B.m_[l] =
add(B.m_[l],
1822 mul(transforms[col - 2 - i].m_[l * t_col + m],
1824 B.m_[l] =
expand(B.m_[l]);
1831 RCP<const Basic> det_berkowitz(
const DenseMatrix &A)
1835 berkowitz(A, polys);
1836 DenseMatrix poly = polys[polys.
size() - 1];
1838 if (polys.
size() % 2 == 1)
1839 return mul(minus_one, poly.get(poly.nrows() - 1, 0));
1841 return poly.get(poly.nrows() - 1, 0);
1844 void char_poly(
const DenseMatrix &A, DenseMatrix &B)
1846 SYMENGINE_ASSERT(B.ncols() == 1 and B.nrows() == A.nrows() + 1);
1847 SYMENGINE_ASSERT(A.nrows() == A.ncols());
1851 berkowitz(A, polys);
1852 B = polys[polys.
size() - 1];
1855 void inverse_fraction_free_LU(
const DenseMatrix &A, DenseMatrix &B)
1857 SYMENGINE_ASSERT(A.row_ == A.col_ and B.row_ == B.col_
1858 and B.row_ == A.row_);
1860 unsigned n = A.row_, i;
1861 DenseMatrix LU = DenseMatrix(n, n);
1862 DenseMatrix e = DenseMatrix(n, 1);
1863 DenseMatrix x = DenseMatrix(n, 1);
1864 DenseMatrix x_ = DenseMatrix(n, 1);
1867 for (i = 0; i < n * n; i++) {
1871 for (i = 0; i < n; i++) {
1877 fraction_free_LU(A, LU);
1882 for (
unsigned j = 0; j < n; j++) {
1885 forward_substitution(LU, e, x_);
1886 back_substitution(LU, x_, x);
1888 for (i = 0; i < n; i++)
1889 B.m_[i * n + j] = x.m_[i];
1895 void inverse_LU(
const DenseMatrix &A, DenseMatrix &B)
1897 SYMENGINE_ASSERT(A.row_ == A.col_ and B.row_ == B.col_
1898 and B.row_ == A.row_);
1899 DenseMatrix e = DenseMatrix(A.row_, A.col_);
1904 void inverse_pivoted_LU(
const DenseMatrix &A, DenseMatrix &B)
1906 SYMENGINE_ASSERT(A.row_ == A.col_ and B.row_ == B.col_
1907 and B.row_ == A.row_);
1908 DenseMatrix e = DenseMatrix(A.row_, A.col_);
1910 pivoted_LU_solve(A, e, B);
1913 void inverse_gauss_jordan(
const DenseMatrix &A, DenseMatrix &B)
1915 SYMENGINE_ASSERT(A.row_ == A.col_ and B.row_ == B.col_
1916 and B.row_ == A.row_);
1918 unsigned n = A.row_;
1919 DenseMatrix e = DenseMatrix(n, n);
1922 for (
unsigned i = 0; i < n; i++)
1923 for (
unsigned j = 0; j < n; j++) {
1925 e.m_[i * n + j] = zero;
1927 e.m_[i * n + i] = one;
1929 B.m_[i * n + j] = zero;
1932 fraction_free_gauss_jordan_solve(A, e, B);
1937 void dot(
const DenseMatrix &A,
const DenseMatrix &B, DenseMatrix &C)
1939 if (A.col_ == B.row_) {
1941 DenseMatrix tmp1 = DenseMatrix(A.col_, A.row_);
1943 DenseMatrix tmp2 = DenseMatrix(B.col_, B.row_);
1945 C.resize(tmp1.row_, tmp2.col_);
1946 mul_dense_dense(tmp1, tmp2, C);
1948 C.resize(A.row_, B.col_);
1949 mul_dense_dense(A, B, C);
1951 C.resize(1, C.row_ * C.col_);
1952 }
else if (A.col_ == B.col_) {
1953 DenseMatrix tmp2 = DenseMatrix(B.col_, B.row_);
1956 }
else if (A.row_ == B.row_) {
1957 DenseMatrix tmp1 = DenseMatrix(A.col_, A.row_);
1961 throw SymEngineException(
"Dimensions incorrect for dot product");
1965 void cross(
const DenseMatrix &A,
const DenseMatrix &B, DenseMatrix &C)
1967 SYMENGINE_ASSERT((A.row_ * A.col_ == 3 and B.row_ * B.col_ == 3)
1968 and (A.row_ == C.row_ and A.col_ == C.col_));
1970 C.m_[0] =
sub(
mul(A.m_[1], B.m_[2]),
mul(A.m_[2], B.m_[1]));
1971 C.m_[1] =
sub(
mul(A.m_[2], B.m_[0]),
mul(A.m_[0], B.m_[2]));
1972 C.m_[2] =
sub(
mul(A.m_[0], B.m_[1]),
mul(A.m_[1], B.m_[0]));
1975 RCP<const Set> eigen_values(
const DenseMatrix &A)
1977 unsigned n = A.nrows();
1978 if (A.is_lower() or A.is_upper()) {
1979 RCP<const Set> eigenvals =
emptyset();
1981 for (
unsigned i = 0; i < n; ++i) {
1988 DenseMatrix B = DenseMatrix(A.nrows() + 1, 1);
1990 map_int_Expr coeffs;
1991 auto nr = A.nrows();
1992 for (
unsigned i = 0; i <= nr; i++)
1993 insert(coeffs, nr - i, B.get(i, 0));
1994 auto lambda =
symbol(
"lambda");
1995 return solve_poly(uexpr_poly(lambda, coeffs), lambda);
2001 void eye(DenseMatrix &A,
int k)
2003 if ((k >= 0 and (
unsigned) k >= A.col_) or k + A.row_ <= 0) {
2007 vec_basic v = vec_basic(k > 0 ? A.col_ - k : A.row_ + k, one);
2013 void diag(DenseMatrix &A, vec_basic &v,
int k)
2015 SYMENGINE_ASSERT(v.size() > 0);
2017 unsigned k_ = std::abs(k);
2020 for (
unsigned i = 0; i < A.row_; i++) {
2021 for (
unsigned j = 0; j < A.col_; j++) {
2022 if (j != (
unsigned)k) {
2023 A.m_[i * A.col_ + j] = zero;
2025 A.m_[i * A.col_ + j] = v[k - k_];
2033 for (
unsigned j = 0; j < A.col_; j++) {
2034 for (
unsigned i = 0; i < A.row_; i++) {
2035 if (i != (
unsigned)k) {
2036 A.m_[i * A.col_ + j] = zero;
2038 A.m_[i * A.col_ + j] = v[k - k_];
2047 void ones(DenseMatrix &A)
2049 for (
unsigned i = 0; i < A.row_ * A.col_; i++) {
2055 void zeros(DenseMatrix &A)
2057 for (
unsigned i = 0; i < A.row_ * A.col_; i++) {
2064 #if HAVE_SYMENGINE_RTTI
2066 unsigned short major = SYMENGINE_MAJOR_VERSION;
2067 unsigned short minor = SYMENGINE_MINOR_VERSION;
2069 major, minor, row_, col_, m_);
2072 throw NotImplementedError(
"Serialization not implemented in no-rtti mode");
2078 #if HAVE_SYMENGINE_RTTI
2079 unsigned short major, minor;
2084 iarchive(major, minor);
2085 if (major != SYMENGINE_MAJOR_VERSION or minor != SYMENGINE_MINOR_VERSION) {
2086 throw SerializationError(StreamFmt()
2087 <<
"SymEngine-" << SYMENGINE_MAJOR_VERSION
2088 <<
"." << SYMENGINE_MINOR_VERSION
2089 <<
" was asked to deserialize an object "
2090 <<
"created using SymEngine-" << major <<
"."
2093 iarchive(row, col, obj);
2096 throw NotImplementedError(
"Serialization not implemented in no-rtti mode");
Classes and functions relating to the binary operation of addition.
RCP< const Basic > sub(const RCP< const Basic > &a, const RCP< const Basic > &b)
Substracts b from a.
RCP< const Basic > add(const RCP< const Basic > &a, const RCP< const Basic > &b)
Adds two objects (safely).
Main namespace for SymEngine package.
std::enable_if< std::is_integral< T >::value, RCP< const Integer > >::type integer(T i)
bool is_number_and_zero(const Basic &b)
RCP< const Symbol > symbol(const std::string &name)
inline version to return Symbol
bool eq(const Basic &a, const Basic &b)
Checks equality for a and b
RCP< const EmptySet > emptyset()
RCP< const Basic > sign(const RCP< const Basic > &arg)
Canonicalize Sign.
RCP< const Basic > abs(const RCP< const Basic > &arg)
Canonicalize Abs:
RCP< const Basic > mul(const RCP< const Basic > &a, const RCP< const Basic > &b)
Multiplication.
void insert(T1 &m, const T2 &first, const T3 &second)
tribool is_zero(const Basic &b, const Assumptions *assumptions=nullptr)
Check if a number is zero.
RCP< const Set > finiteset(const set_basic &container)
RCP< const Basic > expand(const RCP< const Basic > &self, bool deep=true)
Expands self
RCP< const Basic > ssubs(const RCP< const Basic > &x, const map_basic_basic &subs_dict, bool cache=true)
SymPy compatible subs.
RCP< const Basic > conjugate(const RCP< const Basic > &arg)
Canonicalize Conjugate.